1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
12 // Lesser General Public  License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
17 
18 
19 
20 // Local includes
21 #include "libmesh/dof_map.h"
22 
23 // libMesh includes
24 #include "libmesh/coupling_matrix.h"
25 #include "libmesh/default_coupling.h"
26 #include "libmesh/dense_matrix.h"
27 #include "libmesh/dense_vector_base.h"
28 #include "libmesh/dirichlet_boundaries.h"
29 #include "libmesh/elem.h"
30 #include "libmesh/enum_to_string.h"
31 #include "libmesh/fe_interface.h"
32 #include "libmesh/fe_type.h"
33 #include "libmesh/fe_base.h" // FEBase::build() for continuity test
34 #include "libmesh/ghosting_functor.h"
35 #include "libmesh/int_range.h"
36 #include "libmesh/libmesh_logging.h"
37 #include "libmesh/mesh_base.h"
38 #include "libmesh/mesh_subdivision_support.h"
39 #include "libmesh/mesh_tools.h"
40 #include "libmesh/numeric_vector.h"
41 #include "libmesh/periodic_boundaries.h"
42 #include "libmesh/sparse_matrix.h"
43 #include "libmesh/sparsity_pattern.h"
44 #include "libmesh/threads.h"
45 #include "libmesh/auto_ptr.h" // libmesh_make_unique
46 
47 // TIMPI includes
48 #include "timpi/parallel_implementation.h"
49 #include "timpi/parallel_sync.h"
50 
51 // C++ Includes
52 #include <set>
53 #include <algorithm> // for std::fill, std::equal_range, std::max, std::lower_bound, etc.
54 #include <sstream>
55 #include <unordered_map>
56 
57 namespace libMesh
58 {
59 
60 // ------------------------------------------------------------
61 // DofMap member functions
62 std::unique_ptr<SparsityPattern::Build>
build_sparsity(const MeshBase & mesh)63 DofMap::build_sparsity (const MeshBase & mesh) const
64 {
65   libmesh_assert (mesh.is_prepared());
66 
67   LOG_SCOPE("build_sparsity()", "DofMap");
68 
69   // Compute the sparsity structure of the global matrix.  This can be
70   // fed into a PetscMatrix to allocate exactly the number of nonzeros
71   // necessary to store the matrix.  This algorithm should be linear
72   // in the (# of elements)*(# nodes per element)
73 
74   // We can be more efficient in the threaded sparsity pattern assembly
75   // if we don't need the exact pattern.  For some sparse matrix formats
76   // a good upper bound will suffice.
77 
78   // See if we need to include sparsity pattern entries for coupling
79   // between neighbor dofs
80   bool implicit_neighbor_dofs = this->use_coupled_neighbor_dofs(mesh);
81 
82   // We can compute the sparsity pattern in parallel on multiple
83   // threads.  The goal is for each thread to compute the full sparsity
84   // pattern for a subset of elements.  These sparsity patterns can
85   // be efficiently merged in the SparsityPattern::Build::join()
86   // method, especially if there is not too much overlap between them.
87   // Even better, if the full sparsity pattern is not needed then
88   // the number of nonzeros per row can be estimated from the
89   // sparsity patterns created on each thread.
90   auto sp = libmesh_make_unique<SparsityPattern::Build>
91     (mesh,
92      *this,
93      this->_dof_coupling,
94      this->_coupling_functors,
95      implicit_neighbor_dofs,
96      need_full_sparsity_pattern);
97 
98   Threads::parallel_reduce (ConstElemRange (mesh.active_local_elements_begin(),
99                                             mesh.active_local_elements_end()), *sp);
100 
101   sp->parallel_sync();
102 
103 #ifndef NDEBUG
104   // Avoid declaring these variables unless asserts are enabled.
105   const processor_id_type proc_id        = mesh.processor_id();
106   const dof_id_type n_dofs_on_proc = this->n_dofs_on_processor(proc_id);
107 #endif
108   libmesh_assert_equal_to (sp->sparsity_pattern.size(), n_dofs_on_proc);
109 
110   // Check to see if we have any extra stuff to add to the sparsity_pattern
111   if (_extra_sparsity_function)
112     {
113       if (_augment_sparsity_pattern)
114         {
115           libmesh_here();
116           libMesh::out << "WARNING:  You have specified both an extra sparsity function and object.\n"
117                        << "          Are you sure this is what you meant to do??"
118                        << std::endl;
119         }
120 
121       _extra_sparsity_function
122         (sp->sparsity_pattern, sp->n_nz,
123          sp->n_oz, _extra_sparsity_context);
124     }
125 
126   if (_augment_sparsity_pattern)
127     _augment_sparsity_pattern->augment_sparsity_pattern
128       (sp->sparsity_pattern, sp->n_nz, sp->n_oz);
129 
130   return std::unique_ptr<SparsityPattern::Build>(sp.release());
131 }
132 
133 
134 
DofMap(const unsigned int number,MeshBase & mesh)135 DofMap::DofMap(const unsigned int number,
136                MeshBase & mesh) :
137   ParallelObject (mesh.comm()),
138   _dof_coupling(nullptr),
139   _error_on_constraint_loop(false),
140   _variables(),
141   _variable_groups(),
142   _variable_group_numbers(),
143   _sys_number(number),
144   _mesh(mesh),
145   _matrices(),
146   _first_df(),
147   _end_df(),
148   _first_scalar_df(),
149   _send_list(),
150   _augment_sparsity_pattern(nullptr),
151   _extra_sparsity_function(nullptr),
152   _extra_sparsity_context(nullptr),
153   _augment_send_list(nullptr),
154   _extra_send_list_function(nullptr),
155   _extra_send_list_context(nullptr),
156   _default_coupling(libmesh_make_unique<DefaultCoupling>()),
157   _default_evaluating(libmesh_make_unique<DefaultCoupling>()),
158   need_full_sparsity_pattern(false),
159   _n_nz(nullptr),
160   _n_oz(nullptr),
161   _n_dfs(0),
162   _n_SCALAR_dofs(0)
163 #ifdef LIBMESH_ENABLE_AMR
164   , _n_old_dfs(0),
165   _first_old_df(),
166   _end_old_df(),
167   _first_old_scalar_df()
168 #endif
169 #ifdef LIBMESH_ENABLE_CONSTRAINTS
170   , _dof_constraints()
171   , _stashed_dof_constraints()
172   , _primal_constraint_values()
173   , _adjoint_constraint_values()
174 #endif
175 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
176   , _node_constraints()
177 #endif
178 #ifdef LIBMESH_ENABLE_PERIODIC
179   , _periodic_boundaries(libmesh_make_unique<PeriodicBoundaries>())
180 #endif
181 #ifdef LIBMESH_ENABLE_DIRICHLET
182   , _dirichlet_boundaries(libmesh_make_unique<DirichletBoundaries>())
183   , _adjoint_dirichlet_boundaries()
184 #endif
185   , _implicit_neighbor_dofs_initialized(false),
186   _implicit_neighbor_dofs(false)
187 {
188   _matrices.clear();
189 
190   _default_coupling->set_mesh(&_mesh);
191   _default_evaluating->set_mesh(&_mesh);
192   _default_evaluating->set_n_levels(1);
193 
194 #ifdef LIBMESH_ENABLE_PERIODIC
195   _default_coupling->set_periodic_boundaries(_periodic_boundaries.get());
196   _default_evaluating->set_periodic_boundaries(_periodic_boundaries.get());
197 #endif
198 
199   this->add_coupling_functor(*_default_coupling);
200   this->add_algebraic_ghosting_functor(*_default_evaluating);
201 }
202 
203 
204 
205 // Destructor
~DofMap()206 DofMap::~DofMap()
207 {
208   this->clear();
209 
210   // clear() resets all but the default DofMap-based functors.  We
211   // need to remove those from the mesh too before we die.
212   _mesh.remove_ghosting_functor(*_default_coupling);
213   _mesh.remove_ghosting_functor(*_default_evaluating);
214 
215 #ifdef LIBMESH_ENABLE_DIRICHLET
216   for (auto & bnd : _adjoint_dirichlet_boundaries)
217     delete bnd;
218 #endif
219 }
220 
221 
222 #ifdef LIBMESH_ENABLE_PERIODIC
223 
is_periodic_boundary(const boundary_id_type boundaryid)224 bool DofMap::is_periodic_boundary (const boundary_id_type boundaryid) const
225 {
226   if (_periodic_boundaries->count(boundaryid) != 0)
227     return true;
228 
229   return false;
230 }
231 
232 #endif
233 
234 
235 
236 // void DofMap::add_variable (const Variable & var)
237 // {
238 //   libmesh_not_implemented();
239 //   _variables.push_back (var);
240 // }
241 
242 
set_error_on_cyclic_constraint(bool error_on_cyclic_constraint)243 void DofMap::set_error_on_cyclic_constraint(bool error_on_cyclic_constraint)
244 {
245   // This function will eventually be officially libmesh_deprecated();
246   // Call DofMap::set_error_on_constraint_loop() instead.
247   set_error_on_constraint_loop(error_on_cyclic_constraint);
248 }
249 
set_error_on_constraint_loop(bool error_on_constraint_loop)250 void DofMap::set_error_on_constraint_loop(bool error_on_constraint_loop)
251 {
252   _error_on_constraint_loop = error_on_constraint_loop;
253 }
254 
255 
256 
add_variable_group(const VariableGroup & var_group)257 void DofMap::add_variable_group (const VariableGroup & var_group)
258 {
259   const unsigned int vg = cast_int<unsigned int>(_variable_groups.size());
260 
261   _variable_groups.push_back(var_group);
262 
263   VariableGroup & new_var_group = _variable_groups.back();
264 
265   for (auto var : make_range(new_var_group.n_variables()))
266     {
267       _variables.push_back (new_var_group(var));
268       _variable_group_numbers.push_back (vg);
269     }
270 }
271 
272 
273 
attach_matrix(SparseMatrix<Number> & matrix)274 void DofMap::attach_matrix (SparseMatrix<Number> & matrix)
275 {
276   parallel_object_only();
277 
278   // We shouldn't be trying to re-attach the same matrices repeatedly
279   libmesh_assert (std::find(_matrices.begin(), _matrices.end(),
280                             &matrix) == _matrices.end());
281 
282   _matrices.push_back(&matrix);
283 
284   matrix.attach_dof_map (*this);
285 
286   // If we've already computed sparsity, then it's too late
287   // to wait for "compute_sparsity" to help with sparse matrix
288   // initialization, and we need to handle this matrix individually
289   bool computed_sparsity_already =
290     ((_n_nz && !_n_nz->empty()) ||
291      (_n_oz && !_n_oz->empty()));
292   this->comm().max(computed_sparsity_already);
293   if (computed_sparsity_already &&
294       matrix.need_full_sparsity_pattern())
295     {
296       // We'd better have already computed the full sparsity pattern
297       // if we need it here
298       libmesh_assert(need_full_sparsity_pattern);
299       libmesh_assert(_sp.get());
300 
301       matrix.update_sparsity_pattern (_sp->sparsity_pattern);
302     }
303 
304   if (matrix.need_full_sparsity_pattern())
305     need_full_sparsity_pattern = true;
306 }
307 
308 
309 
is_attached(SparseMatrix<Number> & matrix)310 bool DofMap::is_attached (SparseMatrix<Number> & matrix)
311 {
312   return (std::find(_matrices.begin(), _matrices.end(),
313                     &matrix) != _matrices.end());
314 }
315 
316 
317 
node_ptr(MeshBase & mesh,dof_id_type i)318 DofObject * DofMap::node_ptr(MeshBase & mesh, dof_id_type i) const
319 {
320   return mesh.node_ptr(i);
321 }
322 
323 
324 
elem_ptr(MeshBase & mesh,dof_id_type i)325 DofObject * DofMap::elem_ptr(MeshBase & mesh, dof_id_type i) const
326 {
327   return mesh.elem_ptr(i);
328 }
329 
330 
331 
332 template <typename iterator_type>
set_nonlocal_dof_objects(iterator_type objects_begin,iterator_type objects_end,MeshBase & mesh,dofobject_accessor objects)333 void DofMap::set_nonlocal_dof_objects(iterator_type objects_begin,
334                                       iterator_type objects_end,
335                                       MeshBase & mesh,
336                                       dofobject_accessor objects)
337 {
338   // This function must be run on all processors at once
339   parallel_object_only();
340 
341   // First, iterate over local objects to find out how many
342   // are on each processor
343   std::unordered_map<processor_id_type, dof_id_type> ghost_objects_from_proc;
344 
345   iterator_type it  = objects_begin;
346 
347   for (; it != objects_end; ++it)
348     {
349       DofObject * obj = *it;
350 
351       if (obj)
352         {
353           processor_id_type obj_procid = obj->processor_id();
354           // We'd better be completely partitioned by now
355           libmesh_assert_not_equal_to (obj_procid, DofObject::invalid_processor_id);
356           ghost_objects_from_proc[obj_procid]++;
357         }
358     }
359 
360   // Request sets to send to each processor
361   std::map<processor_id_type, std::vector<dof_id_type>>
362     requested_ids;
363 
364   // We know how many of our objects live on each processor, so
365   // reserve() space for requests from each.
366   for (auto pair : ghost_objects_from_proc)
367     {
368       const processor_id_type p = pair.first;
369       if (p != this->processor_id())
370         requested_ids[p].reserve(pair.second);
371     }
372 
373   for (it = objects_begin; it != objects_end; ++it)
374     {
375       DofObject * obj = *it;
376       if (obj->processor_id() != DofObject::invalid_processor_id)
377         requested_ids[obj->processor_id()].push_back(obj->id());
378     }
379 #ifdef DEBUG
380   for (auto p : make_range(this->n_processors()))
381     {
382       if (ghost_objects_from_proc.count(p))
383         libmesh_assert_equal_to (requested_ids[p].size(), ghost_objects_from_proc[p]);
384       else
385         libmesh_assert(!requested_ids.count(p));
386     }
387 #endif
388 
389   typedef std::vector<dof_id_type> datum;
390 
391   auto gather_functor =
392     [this, &mesh, &objects]
393     (processor_id_type,
394      const std::vector<dof_id_type> & ids,
395      std::vector<datum> & data)
396     {
397       // Fill those requests
398       const unsigned int
399         sys_num      = this->sys_number(),
400         n_var_groups = this->n_variable_groups();
401 
402       const std::size_t query_size = ids.size();
403 
404       data.resize(query_size);
405       for (auto & d : data)
406         d.resize(2 * n_var_groups);
407 
408       for (std::size_t i=0; i != query_size; ++i)
409         {
410           DofObject * requested = (this->*objects)(mesh, ids[i]);
411           libmesh_assert(requested);
412           libmesh_assert_equal_to (requested->processor_id(), this->processor_id());
413           libmesh_assert_equal_to (requested->n_var_groups(sys_num), n_var_groups);
414           for (unsigned int vg=0; vg != n_var_groups; ++vg)
415             {
416               unsigned int n_comp_g =
417                 requested->n_comp_group(sys_num, vg);
418               data[i][vg] = n_comp_g;
419               dof_id_type my_first_dof = n_comp_g ?
420                 requested->vg_dof_base(sys_num, vg) : 0;
421               libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
422               data[i][n_var_groups+vg] = my_first_dof;
423             }
424         }
425     };
426 
427   auto action_functor =
428     [this, &mesh, &objects]
429     (processor_id_type libmesh_dbg_var(pid),
430      const std::vector<dof_id_type> & ids,
431      const std::vector<datum> & data)
432     {
433       const unsigned int
434         sys_num      = this->sys_number(),
435         n_var_groups = this->n_variable_groups();
436 
437       // Copy the id changes we've now been informed of
438       for (auto i : index_range(ids))
439         {
440           DofObject * requested = (this->*objects)(mesh, ids[i]);
441           libmesh_assert(requested);
442           libmesh_assert_equal_to (requested->processor_id(), pid);
443           for (unsigned int vg=0; vg != n_var_groups; ++vg)
444             {
445               unsigned int n_comp_g =
446                 cast_int<unsigned int>(data[i][vg]);
447               requested->set_n_comp_group(sys_num, vg, n_comp_g);
448               if (n_comp_g)
449                 {
450                   dof_id_type my_first_dof = data[i][n_var_groups+vg];
451                   libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
452                   requested->set_vg_dof_base
453                     (sys_num, vg, my_first_dof);
454                 }
455             }
456         }
457     };
458 
459   datum * ex = nullptr;
460   Parallel::pull_parallel_vector_data
461     (this->comm(), requested_ids, gather_functor, action_functor, ex);
462 
463 #ifdef DEBUG
464   // Double check for invalid dofs
465   for (it = objects_begin; it != objects_end; ++it)
466     {
467       DofObject * obj = *it;
468       libmesh_assert (obj);
469       unsigned int num_variables = obj->n_vars(this->sys_number());
470       for (unsigned int v=0; v != num_variables; ++v)
471         {
472           unsigned int n_comp =
473             obj->n_comp(this->sys_number(), v);
474           dof_id_type my_first_dof = n_comp ?
475             obj->dof_number(this->sys_number(), v, 0) : 0;
476           libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
477         }
478     }
479 #endif
480 }
481 
482 
483 
reinit(MeshBase & mesh)484 void DofMap::reinit(MeshBase & mesh)
485 {
486   libmesh_assert (mesh.is_prepared());
487 
488   LOG_SCOPE("reinit()", "DofMap");
489 
490   // We ought to reconfigure our default coupling functor.
491   //
492   // The user might have removed it from our coupling functors set,
493   // but if so, who cares, this reconfiguration is cheap.
494 
495   // Avoid calling set_dof_coupling() with an empty/non-nullptr
496   // _dof_coupling matrix which may happen when there are actually no
497   // variables on the system.
498   if (this->_dof_coupling && this->_dof_coupling->empty() && !this->n_variables())
499     this->_dof_coupling = nullptr;
500   _default_coupling->set_dof_coupling(this->_dof_coupling);
501 
502   // By default we may want 0 or 1 levels of coupling
503   unsigned int standard_n_levels =
504     this->use_coupled_neighbor_dofs(mesh);
505   _default_coupling->set_n_levels
506     (std::max(_default_coupling->n_levels(), standard_n_levels));
507 
508   // But we *don't* want to restrict to a CouplingMatrix unless the
509   // user does so manually; the original libMesh behavior was to put
510   // ghost indices on the send_list regardless of variable.
511   //_default_evaluating->set_dof_coupling(this->_dof_coupling);
512 
513   const unsigned int
514     sys_num      = this->sys_number(),
515     n_var_groups = this->n_variable_groups();
516 
517   // The DofObjects need to know how many variable groups we have, and
518   // how many variables there are in each group.
519   std::vector<unsigned int> n_vars_per_group; /**/ n_vars_per_group.reserve (n_var_groups);
520 
521   for (unsigned int vg=0; vg<n_var_groups; vg++)
522     n_vars_per_group.push_back (this->variable_group(vg).n_variables());
523 
524 #ifdef LIBMESH_ENABLE_AMR
525 
526   //------------------------------------------------------------
527   // Clear the old_dof_objects for all the nodes
528   // and elements so that we can overwrite them
529   for (auto & node : mesh.node_ptr_range())
530     {
531       node->clear_old_dof_object();
532       libmesh_assert (!node->old_dof_object);
533     }
534 
535   for (auto & elem : mesh.element_ptr_range())
536     {
537       elem->clear_old_dof_object();
538       libmesh_assert (!elem->old_dof_object);
539     }
540 
541 
542   //------------------------------------------------------------
543   // Set the old_dof_objects for the elements that
544   // weren't just created, if these old dof objects
545   // had variables
546   for (auto & elem : mesh.element_ptr_range())
547     {
548       // Skip the elements that were just refined
549       if (elem->refinement_flag() == Elem::JUST_REFINED)
550         continue;
551 
552       for (Node & node : elem->node_ref_range())
553         if (node.old_dof_object == nullptr)
554           if (node.has_dofs(sys_num))
555             node.set_old_dof_object();
556 
557       libmesh_assert (!elem->old_dof_object);
558 
559       if (elem->has_dofs(sys_num))
560         elem->set_old_dof_object();
561     }
562 
563 #endif // #ifdef LIBMESH_ENABLE_AMR
564 
565 
566   //------------------------------------------------------------
567   // Then set the number of variables for each \p DofObject
568   // equal to n_variables() for this system.  This will
569   // handle new \p DofObjects that may have just been created
570 
571   // All the nodes
572   for (auto & node : mesh.node_ptr_range())
573     node->set_n_vars_per_group(sys_num, n_vars_per_group);
574 
575   // All the elements
576   for (auto & elem : mesh.element_ptr_range())
577     elem->set_n_vars_per_group(sys_num, n_vars_per_group);
578 
579   // Zero _n_SCALAR_dofs, it will be updated below.
580   this->_n_SCALAR_dofs = 0;
581 
582   //------------------------------------------------------------
583   // Next allocate space for the DOF indices
584   for (unsigned int vg=0; vg<n_var_groups; vg++)
585     {
586       const VariableGroup & vg_description = this->variable_group(vg);
587 
588       const unsigned int n_var_in_group = vg_description.n_variables();
589       const FEType & base_fe_type        = vg_description.type();
590 
591       // Don't need to loop over elements for a SCALAR variable
592       // Just increment _n_SCALAR_dofs
593       if (base_fe_type.family == SCALAR)
594         {
595           this->_n_SCALAR_dofs += base_fe_type.order.get_order()*n_var_in_group;
596           continue;
597         }
598 
599       // This should be constant even on p-refined elements
600       const bool extra_hanging_dofs =
601         FEInterface::extra_hanging_dofs(base_fe_type);
602 
603       // For all the active elements, count vertex degrees of freedom.
604       for (auto & elem : mesh.active_element_ptr_range())
605         {
606           libmesh_assert(elem);
607 
608           // Skip the numbering if this variable is
609           // not active on this element's subdomain
610           if (!vg_description.active_on_subdomain(elem->subdomain_id()))
611             continue;
612 
613           FEType fe_type = base_fe_type;
614 
615 #ifdef LIBMESH_ENABLE_AMR
616           const ElemType type = elem->type();
617 
618           // Make sure we haven't done more p refinement than we can
619           // handle
620           if (elem->p_level() + base_fe_type.order >
621               FEInterface::max_order(base_fe_type, type))
622             {
623               libmesh_assert_less_msg(base_fe_type.order.get_order(),
624                                       FEInterface::max_order(base_fe_type,type),
625                                       "ERROR: Finite element "
626                                       << Utility::enum_to_string(base_fe_type.family)
627                                       << " on geometric element "
628                                       << Utility::enum_to_string(type)
629                                       << "\nonly supports FEInterface::max_order = "
630                                       << FEInterface::max_order(base_fe_type,type)
631                                       << ", not fe_type.order = "
632                                       << base_fe_type.order);
633 
634 #  ifdef DEBUG
635               libMesh::err << "WARNING: Finite element "
636                            << Utility::enum_to_string(base_fe_type.family)
637                            << " on geometric element "
638                            << Utility::enum_to_string(type) << std::endl
639                            << "could not be p refined past FEInterface::max_order = "
640                            << FEInterface::max_order(base_fe_type,type)
641                            << std::endl;
642 #  endif
643               elem->set_p_level(FEInterface::max_order(base_fe_type,type)
644                                 - base_fe_type.order);
645             }
646 #endif
647 
648           // Allocate the vertex DOFs
649           for (auto n : elem->node_index_range())
650             {
651               Node & node = elem->node_ref(n);
652 
653               if (elem->is_vertex(n))
654                 {
655                   const unsigned int old_node_dofs =
656                     node.n_comp_group(sys_num, vg);
657 
658                   const unsigned int vertex_dofs =
659                     std::max(FEInterface::n_dofs_at_node(fe_type, elem, n),
660                              old_node_dofs);
661 
662                   // Some discontinuous FEs have no vertex dofs
663                   if (vertex_dofs > old_node_dofs)
664                     {
665                       node.set_n_comp_group(sys_num, vg,
666                                             vertex_dofs);
667 
668                       // Abusing dof_number to set a "this is a
669                       // vertex" flag
670                       node.set_vg_dof_base(sys_num, vg,
671                                            vertex_dofs);
672 
673                       // libMesh::out << "sys_num,vg,old_node_dofs,vertex_dofs="
674                       //       << sys_num << ","
675                       //       << vg << ","
676                       //       << old_node_dofs << ","
677                       //       << vertex_dofs << '\n',
678                       // node.debug_buffer();
679 
680                       // libmesh_assert_equal_to (vertex_dofs, node.n_comp(sys_num, vg));
681                       // libmesh_assert_equal_to (vertex_dofs, node.vg_dof_base(sys_num, vg));
682                     }
683                 }
684             }
685         } // done counting vertex dofs
686 
687       // count edge & face dofs next
688       for (auto & elem : mesh.active_element_ptr_range())
689         {
690           libmesh_assert(elem);
691 
692           // Skip the numbering if this variable is
693           // not active on this element's subdomain
694           if (!vg_description.active_on_subdomain(elem->subdomain_id()))
695             continue;
696 
697           // Allocate the edge and face DOFs
698           for (auto n : elem->node_index_range())
699             {
700               Node & node = elem->node_ref(n);
701 
702               const unsigned int old_node_dofs =
703                 node.n_comp_group(sys_num, vg);
704 
705               const unsigned int vertex_dofs = old_node_dofs?
706                 cast_int<unsigned int>(node.vg_dof_base (sys_num,vg)):0;
707 
708               const unsigned int new_node_dofs =
709                 FEInterface::n_dofs_at_node(base_fe_type, elem, n);
710 
711               // We've already allocated vertex DOFs
712               if (elem->is_vertex(n))
713                 {
714                   libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
715                   // //if (vertex_dofs < new_node_dofs)
716                   //   libMesh::out << "sys_num,vg,old_node_dofs,vertex_dofs,new_node_dofs="
717                   //                << sys_num << ","
718                   //                << vg << ","
719                   //                << old_node_dofs << ","
720                   //                << vertex_dofs << ","
721                   //                << new_node_dofs << '\n',
722                   //     node.debug_buffer();
723 
724                   libmesh_assert_greater_equal (vertex_dofs,   new_node_dofs);
725                 }
726               // We need to allocate the rest
727               else
728                 {
729                   // If this has no dofs yet, it needs no vertex
730                   // dofs, so we just give it edge or face dofs
731                   if (!old_node_dofs)
732                     {
733                       node.set_n_comp_group(sys_num, vg,
734                                             new_node_dofs);
735                       // Abusing dof_number to set a "this has no
736                       // vertex dofs" flag
737                       if (new_node_dofs)
738                         node.set_vg_dof_base(sys_num, vg, 0);
739                     }
740 
741                   // If this has dofs, but has no vertex dofs,
742                   // it may still need more edge or face dofs if
743                   // we're p-refined.
744                   else if (vertex_dofs == 0)
745                     {
746                       if (new_node_dofs > old_node_dofs)
747                         {
748                           node.set_n_comp_group(sys_num, vg,
749                                                 new_node_dofs);
750 
751                           node.set_vg_dof_base(sys_num, vg,
752                                                vertex_dofs);
753                         }
754                     }
755                   // If this is another element's vertex,
756                   // add more (non-overlapping) edge/face dofs if
757                   // necessary
758                   else if (extra_hanging_dofs)
759                     {
760                       if (new_node_dofs > old_node_dofs - vertex_dofs)
761                         {
762                           node.set_n_comp_group(sys_num, vg,
763                                                 vertex_dofs + new_node_dofs);
764 
765                           node.set_vg_dof_base(sys_num, vg,
766                                                vertex_dofs);
767                         }
768                     }
769                   // If this is another element's vertex, add any
770                   // (overlapping) edge/face dofs if necessary
771                   else
772                     {
773                       libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
774                       if (new_node_dofs > old_node_dofs)
775                         {
776                           node.set_n_comp_group(sys_num, vg,
777                                                 new_node_dofs);
778 
779                           node.set_vg_dof_base (sys_num, vg,
780                                                 vertex_dofs);
781                         }
782                     }
783                 }
784             }
785           // Allocate the element DOFs
786           const unsigned int dofs_per_elem =
787             FEInterface::n_dofs_per_elem(base_fe_type, elem);
788 
789           elem->set_n_comp_group(sys_num, vg, dofs_per_elem);
790 
791         }
792     } // end loop over variable groups
793 
794   // Calling DofMap::reinit() by itself makes little sense,
795   // so we won't bother with nonlocal DofObjects.
796   // Those will be fixed by distribute_dofs
797 
798   //------------------------------------------------------------
799   // Finally, clear all the current DOF indices
800   // (distribute_dofs expects them cleared!)
801   this->invalidate_dofs(mesh);
802 }
803 
804 
805 
invalidate_dofs(MeshBase & mesh)806 void DofMap::invalidate_dofs(MeshBase & mesh) const
807 {
808   const unsigned int sys_num = this->sys_number();
809 
810   // All the nodes
811   for (auto & node : mesh.node_ptr_range())
812     node->invalidate_dofs(sys_num);
813 
814   // All the active elements.
815   for (auto & elem : mesh.active_element_ptr_range())
816     elem->invalidate_dofs(sys_num);
817 }
818 
819 
820 
clear()821 void DofMap::clear()
822 {
823   // we don't want to clear
824   // the coupling matrix!
825   // It should not change...
826   //_dof_coupling->clear();
827   //
828   // But it would be inconsistent to leave our coupling settings
829   // through a clear()...
830   _dof_coupling = nullptr;
831 
832   // Reset ghosting functor statuses
833   {
834     for (const auto & gf : _coupling_functors)
835       {
836         libmesh_assert(gf);
837         _mesh.remove_ghosting_functor(*gf);
838       }
839     this->_coupling_functors.clear();
840 
841     // Go back to default coupling
842 
843     _default_coupling->set_dof_coupling(this->_dof_coupling);
844     _default_coupling->set_n_levels(this->use_coupled_neighbor_dofs(this->_mesh));
845 
846     this->add_coupling_functor(*_default_coupling);
847   }
848 
849 
850   {
851     for (const auto & gf : _algebraic_ghosting_functors)
852       {
853         libmesh_assert(gf);
854         _mesh.remove_ghosting_functor(*gf);
855       }
856     this->_algebraic_ghosting_functors.clear();
857 
858     // Go back to default send_list generation
859 
860     // _default_evaluating->set_dof_coupling(this->_dof_coupling);
861     _default_evaluating->set_n_levels(1);
862     this->add_algebraic_ghosting_functor(*_default_evaluating);
863   }
864 
865   this->_shared_functors.clear();
866 
867   _variables.clear();
868   _variable_groups.clear();
869   _variable_group_numbers.clear();
870   _first_df.clear();
871   _end_df.clear();
872   _first_scalar_df.clear();
873   this->clear_send_list();
874   this->clear_sparsity();
875   need_full_sparsity_pattern = false;
876 
877 #ifdef LIBMESH_ENABLE_AMR
878 
879   _dof_constraints.clear();
880   _stashed_dof_constraints.clear();
881   _primal_constraint_values.clear();
882   _adjoint_constraint_values.clear();
883   _n_old_dfs = 0;
884   _first_old_df.clear();
885   _end_old_df.clear();
886   _first_old_scalar_df.clear();
887 
888 #endif
889 
890   _matrices.clear();
891 
892   _n_dfs = 0;
893 }
894 
895 
896 
distribute_dofs(MeshBase & mesh)897 void DofMap::distribute_dofs (MeshBase & mesh)
898 {
899   // This function must be run on all processors at once
900   parallel_object_only();
901 
902   // Log how long it takes to distribute the degrees of freedom
903   LOG_SCOPE("distribute_dofs()", "DofMap");
904 
905   libmesh_assert (mesh.is_prepared());
906 
907   const processor_id_type proc_id = this->processor_id();
908   const processor_id_type n_proc  = this->n_processors();
909 
910   //  libmesh_assert_greater (this->n_variables(), 0);
911   libmesh_assert_less (proc_id, n_proc);
912 
913   // re-init in case the mesh has changed
914   this->reinit(mesh);
915 
916   // By default distribute variables in a
917   // var-major fashion, but allow run-time
918   // specification
919   bool node_major_dofs = libMesh::on_command_line ("--node-major-dofs");
920 
921   // The DOF counter, will be incremented as we encounter
922   // new degrees of freedom
923   dof_id_type next_free_dof = 0;
924 
925   // Clear the send list before we rebuild it
926   this->clear_send_list();
927 
928   // Set temporary DOF indices on this processor
929   if (node_major_dofs)
930     this->distribute_local_dofs_node_major (next_free_dof, mesh);
931   else
932     this->distribute_local_dofs_var_major (next_free_dof, mesh);
933 
934   // Get DOF counts on all processors
935   std::vector<dof_id_type> dofs_on_proc(n_proc, 0);
936   this->comm().allgather(next_free_dof, dofs_on_proc);
937 
938   // Resize and fill the _first_df and _end_df arrays
939 #ifdef LIBMESH_ENABLE_AMR
940   _first_old_df = _first_df;
941   _end_old_df = _end_df;
942 #endif
943 
944   _first_df.resize(n_proc);
945   _end_df.resize (n_proc);
946 
947   // Get DOF offsets
948   _first_df[0] = 0;
949   for (processor_id_type i=1; i < n_proc; ++i)
950     _first_df[i] = _end_df[i-1] = _first_df[i-1] + dofs_on_proc[i-1];
951   _end_df[n_proc-1] = _first_df[n_proc-1] + dofs_on_proc[n_proc-1];
952 
953   // Clear all the current DOF indices
954   // (distribute_dofs expects them cleared!)
955   this->invalidate_dofs(mesh);
956 
957   next_free_dof = _first_df[proc_id];
958 
959   // Set permanent DOF indices on this processor
960   if (node_major_dofs)
961     this->distribute_local_dofs_node_major (next_free_dof, mesh);
962   else
963     this->distribute_local_dofs_var_major (next_free_dof, mesh);
964 
965   libmesh_assert_equal_to (next_free_dof, _end_df[proc_id]);
966 
967   //------------------------------------------------------------
968   // At this point, all n_comp and dof_number values on local
969   // DofObjects should be correct, but a DistributedMesh might have
970   // incorrect values on non-local DofObjects.  Let's request the
971   // correct values from each other processor.
972 
973   if (this->n_processors() > 1)
974     {
975       this->set_nonlocal_dof_objects(mesh.nodes_begin(),
976                                      mesh.nodes_end(),
977                                      mesh, &DofMap::node_ptr);
978 
979       this->set_nonlocal_dof_objects(mesh.elements_begin(),
980                                      mesh.elements_end(),
981                                      mesh, &DofMap::elem_ptr);
982     }
983 
984 #ifdef DEBUG
985   {
986     const unsigned int
987       sys_num = this->sys_number();
988 
989     // Processors should all agree on DoF ids for the newly numbered
990     // system.
991     MeshTools::libmesh_assert_valid_dof_ids(mesh, sys_num);
992 
993     // DoF processor ids should match DofObject processor ids
994     for (auto & node : mesh.node_ptr_range())
995       {
996         DofObject const * const dofobj = node;
997         const processor_id_type obj_proc_id = dofobj->processor_id();
998 
999         for (auto v : make_range(dofobj->n_vars(sys_num)))
1000           for (auto c : make_range(dofobj->n_comp(sys_num,v)))
1001             {
1002               const dof_id_type dofid = dofobj->dof_number(sys_num,v,c);
1003               libmesh_assert_greater_equal (dofid, this->first_dof(obj_proc_id));
1004               libmesh_assert_less (dofid, this->end_dof(obj_proc_id));
1005             }
1006       }
1007 
1008     for (auto & elem : mesh.element_ptr_range())
1009       {
1010         DofObject const * const dofobj = elem;
1011         const processor_id_type obj_proc_id = dofobj->processor_id();
1012 
1013         for (auto v : make_range(dofobj->n_vars(sys_num)))
1014           for (auto c : make_range(dofobj->n_comp(sys_num,v)))
1015             {
1016               const dof_id_type dofid = dofobj->dof_number(sys_num,v,c);
1017               libmesh_assert_greater_equal (dofid, this->first_dof(obj_proc_id));
1018               libmesh_assert_less (dofid, this->end_dof(obj_proc_id));
1019             }
1020       }
1021   }
1022 #endif
1023 
1024   // Set the total number of degrees of freedom, then start finding
1025   // SCALAR degrees of freedom
1026 #ifdef LIBMESH_ENABLE_AMR
1027   _n_old_dfs = _n_dfs;
1028   _first_old_scalar_df = _first_scalar_df;
1029 #endif
1030   _n_dfs = _end_df[n_proc-1];
1031   _first_scalar_df.clear();
1032   _first_scalar_df.resize(this->n_variables(), DofObject::invalid_id);
1033   dof_id_type current_SCALAR_dof_index = n_dofs() - n_SCALAR_dofs();
1034 
1035   // Calculate and cache the initial DoF indices for SCALAR variables.
1036   // This is an O(N_vars) calculation so we want to do it once per
1037   // renumbering rather than once per SCALAR_dof_indices() call
1038 
1039   for (auto v : make_range(this->n_variables()))
1040     if (this->variable(v).type().family == SCALAR)
1041       {
1042         _first_scalar_df[v] = current_SCALAR_dof_index;
1043         current_SCALAR_dof_index += this->variable(v).type().order.get_order();
1044       }
1045 
1046   // Allow our GhostingFunctor objects to reinit if necessary
1047   for (const auto & gf : _algebraic_ghosting_functors)
1048     {
1049       libmesh_assert(gf);
1050       gf->dofmap_reinit();
1051     }
1052 
1053   for (const auto & gf : _coupling_functors)
1054     {
1055       libmesh_assert(gf);
1056       gf->dofmap_reinit();
1057     }
1058 
1059   // Note that in the add_neighbors_to_send_list nodes on processor
1060   // boundaries that are shared by multiple elements are added for
1061   // each element.
1062   this->add_neighbors_to_send_list(mesh);
1063 
1064   // Here we used to clean up that data structure; now System and
1065   // EquationSystems call that for us, after we've added constraint
1066   // dependencies to the send_list too.
1067   // this->sort_send_list ();
1068 }
1069 
1070 
local_variable_indices(std::vector<dof_id_type> & idx,const MeshBase & mesh,unsigned int var_num)1071 void DofMap::local_variable_indices(std::vector<dof_id_type> & idx,
1072                                     const MeshBase & mesh,
1073                                     unsigned int var_num) const
1074 {
1075   // Count dofs in the *exact* order that distribute_dofs numbered
1076   // them, so that we can assume ascending indices and use push_back
1077   // instead of find+insert.
1078 
1079   const unsigned int sys_num       = this->sys_number();
1080 
1081   // If this isn't a SCALAR variable, we need to find all its field
1082   // dofs on the mesh
1083   if (this->variable_type(var_num).family != SCALAR)
1084     {
1085       const Variable & var(this->variable(var_num));
1086 
1087       for (auto & elem : mesh.active_local_element_ptr_range())
1088         {
1089           if (!var.active_on_subdomain(elem->subdomain_id()))
1090             continue;
1091 
1092           // Only count dofs connected to active
1093           // elements on this processor.
1094           const unsigned int n_nodes = elem->n_nodes();
1095 
1096           // First get any new nodal DOFS
1097           for (unsigned int n=0; n<n_nodes; n++)
1098             {
1099               Node & node = elem->node_ref(n);
1100 
1101               if (node.processor_id() != this->processor_id())
1102                 continue;
1103 
1104               const unsigned int n_comp = node.n_comp(sys_num, var_num);
1105               for(unsigned int i=0; i<n_comp; i++)
1106                 {
1107                   const dof_id_type index = node.dof_number(sys_num,var_num,i);
1108                   libmesh_assert (this->local_index(index));
1109 
1110                   if (idx.empty() || index > idx.back())
1111                     idx.push_back(index);
1112                 }
1113             }
1114 
1115           // Next get any new element DOFS
1116           const unsigned int n_comp = elem->n_comp(sys_num, var_num);
1117           for (unsigned int i=0; i<n_comp; i++)
1118             {
1119               const dof_id_type index = elem->dof_number(sys_num,var_num,i);
1120               if (idx.empty() || index > idx.back())
1121                 idx.push_back(index);
1122             }
1123         } // done looping over elements
1124 
1125 
1126       // we may have missed assigning DOFs to nodes that we own
1127       // but to which we have no connected elements matching our
1128       // variable restriction criterion.  this will happen, for example,
1129       // if variable V is restricted to subdomain S.  We may not own
1130       // any elements which live in S, but we may own nodes which are
1131       // *connected* to elements which do.  in this scenario these nodes
1132       // will presently have unnumbered DOFs. we need to take care of
1133       // them here since we own them and no other processor will touch them.
1134       for (const auto & node : mesh.local_node_ptr_range())
1135         {
1136           libmesh_assert(node);
1137 
1138           const unsigned int n_comp = node->n_comp(sys_num, var_num);
1139           for (unsigned int i=0; i<n_comp; i++)
1140             {
1141               const dof_id_type index = node->dof_number(sys_num,var_num,i);
1142               if (idx.empty() || index > idx.back())
1143                 idx.push_back(index);
1144             }
1145         }
1146     }
1147   // Otherwise, count up the SCALAR dofs, if we're on the processor
1148   // that holds this SCALAR variable
1149   else if (this->processor_id() == (this->n_processors()-1))
1150     {
1151       std::vector<dof_id_type> di_scalar;
1152       this->SCALAR_dof_indices(di_scalar,var_num);
1153       idx.insert( idx.end(), di_scalar.begin(), di_scalar.end());
1154     }
1155 }
1156 
1157 
distribute_local_dofs_node_major(dof_id_type & next_free_dof,MeshBase & mesh)1158 void DofMap::distribute_local_dofs_node_major(dof_id_type & next_free_dof,
1159                                               MeshBase & mesh)
1160 {
1161   const unsigned int sys_num       = this->sys_number();
1162   const unsigned int n_var_groups  = this->n_variable_groups();
1163 
1164   // Our numbering here must be kept consistent with the numbering
1165   // scheme assumed by DofMap::local_variable_indices!
1166 
1167   //-------------------------------------------------------------------------
1168   // First count and assign temporary numbers to local dofs
1169   for (auto & elem : mesh.active_local_element_ptr_range())
1170     {
1171       // Only number dofs connected to active
1172       // elements on this processor.
1173       const unsigned int n_nodes = elem->n_nodes();
1174 
1175       // First number the nodal DOFS
1176       for (unsigned int n=0; n<n_nodes; n++)
1177         {
1178           Node & node = elem->node_ref(n);
1179 
1180           for (unsigned vg=0; vg<n_var_groups; vg++)
1181             {
1182               const VariableGroup & vg_description(this->variable_group(vg));
1183 
1184               if ((vg_description.type().family != SCALAR) &&
1185                   (vg_description.active_on_subdomain(elem->subdomain_id())))
1186                 {
1187                   // assign dof numbers (all at once) if this is
1188                   // our node and if they aren't already there
1189                   if ((node.n_comp_group(sys_num,vg) > 0) &&
1190                       (node.processor_id() == this->processor_id()) &&
1191                       (node.vg_dof_base(sys_num,vg) ==
1192                        DofObject::invalid_id))
1193                     {
1194                       node.set_vg_dof_base(sys_num, vg,
1195                                            next_free_dof);
1196                       next_free_dof += (vg_description.n_variables()*
1197                                         node.n_comp_group(sys_num,vg));
1198                       //node.debug_buffer();
1199                     }
1200                 }
1201             }
1202         }
1203 
1204       // Now number the element DOFS
1205       for (unsigned vg=0; vg<n_var_groups; vg++)
1206         {
1207           const VariableGroup & vg_description(this->variable_group(vg));
1208 
1209           if ((vg_description.type().family != SCALAR) &&
1210               (vg_description.active_on_subdomain(elem->subdomain_id())))
1211             if (elem->n_comp_group(sys_num,vg) > 0)
1212               {
1213                 libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
1214                                          DofObject::invalid_id);
1215 
1216                 elem->set_vg_dof_base(sys_num,
1217                                       vg,
1218                                       next_free_dof);
1219 
1220                 next_free_dof += (vg_description.n_variables()*
1221                                   elem->n_comp(sys_num,vg));
1222               }
1223         }
1224     } // done looping over elements
1225 
1226 
1227   // we may have missed assigning DOFs to nodes that we own
1228   // but to which we have no connected elements matching our
1229   // variable restriction criterion.  this will happen, for example,
1230   // if variable V is restricted to subdomain S.  We may not own
1231   // any elements which live in S, but we may own nodes which are
1232   // *connected* to elements which do.  in this scenario these nodes
1233   // will presently have unnumbered DOFs. we need to take care of
1234   // them here since we own them and no other processor will touch them.
1235   for (auto & node : mesh.local_node_ptr_range())
1236     for (unsigned vg=0; vg<n_var_groups; vg++)
1237       {
1238         const VariableGroup & vg_description(this->variable_group(vg));
1239 
1240         if (node->n_comp_group(sys_num,vg))
1241           if (node->vg_dof_base(sys_num,vg) == DofObject::invalid_id)
1242             {
1243               node->set_vg_dof_base (sys_num,
1244                                      vg,
1245                                      next_free_dof);
1246 
1247               next_free_dof += (vg_description.n_variables()*
1248                                 node->n_comp(sys_num,vg));
1249             }
1250       }
1251 
1252   // Finally, count up the SCALAR dofs
1253   this->_n_SCALAR_dofs = 0;
1254   for (unsigned vg=0; vg<n_var_groups; vg++)
1255     {
1256       const VariableGroup & vg_description(this->variable_group(vg));
1257 
1258       if (vg_description.type().family == SCALAR)
1259         {
1260           this->_n_SCALAR_dofs += (vg_description.n_variables()*
1261                                    vg_description.type().order.get_order());
1262           continue;
1263         }
1264     }
1265 
1266   // Only increment next_free_dof if we're on the processor
1267   // that holds this SCALAR variable
1268   if (this->processor_id() == (this->n_processors()-1))
1269     next_free_dof += _n_SCALAR_dofs;
1270 
1271 #ifdef DEBUG
1272   {
1273     // libMesh::out << "next_free_dof=" << next_free_dof << std::endl
1274     //       << "_n_SCALAR_dofs=" << _n_SCALAR_dofs << std::endl;
1275 
1276     // Make sure we didn't miss any nodes
1277     MeshTools::libmesh_assert_valid_procids<Node>(mesh);
1278 
1279     for (auto & node : mesh.local_node_ptr_range())
1280       {
1281         unsigned int n_var_g = node->n_var_groups(this->sys_number());
1282         for (unsigned int vg=0; vg != n_var_g; ++vg)
1283           {
1284             unsigned int n_comp_g =
1285               node->n_comp_group(this->sys_number(), vg);
1286             dof_id_type my_first_dof = n_comp_g ?
1287               node->vg_dof_base(this->sys_number(), vg) : 0;
1288             libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
1289           }
1290       }
1291   }
1292 #endif // DEBUG
1293 }
1294 
1295 
1296 
distribute_local_dofs_var_major(dof_id_type & next_free_dof,MeshBase & mesh)1297 void DofMap::distribute_local_dofs_var_major(dof_id_type & next_free_dof,
1298                                              MeshBase & mesh)
1299 {
1300   const unsigned int sys_num      = this->sys_number();
1301   const unsigned int n_var_groups = this->n_variable_groups();
1302 
1303   // Our numbering here must be kept consistent with the numbering
1304   // scheme assumed by DofMap::local_variable_indices!
1305 
1306   //-------------------------------------------------------------------------
1307   // First count and assign temporary numbers to local dofs
1308   for (unsigned vg=0; vg<n_var_groups; vg++)
1309     {
1310       const VariableGroup & vg_description(this->variable_group(vg));
1311 
1312       const unsigned int n_vars_in_group = vg_description.n_variables();
1313 
1314       // Skip the SCALAR dofs
1315       if (vg_description.type().family == SCALAR)
1316         continue;
1317 
1318       for (auto & elem : mesh.active_local_element_ptr_range())
1319         {
1320           // Only number dofs connected to active elements on this
1321           // processor and only variables which are active on on this
1322           // element's subdomain.
1323           if (!vg_description.active_on_subdomain(elem->subdomain_id()))
1324             continue;
1325 
1326           const unsigned int n_nodes = elem->n_nodes();
1327 
1328           // First number the nodal DOFS
1329           for (unsigned int n=0; n<n_nodes; n++)
1330             {
1331               Node & node = elem->node_ref(n);
1332 
1333               // assign dof numbers (all at once) if this is
1334               // our node and if they aren't already there
1335               if ((node.n_comp_group(sys_num,vg) > 0) &&
1336                   (node.processor_id() == this->processor_id()) &&
1337                   (node.vg_dof_base(sys_num,vg) ==
1338                    DofObject::invalid_id))
1339                 {
1340                   node.set_vg_dof_base(sys_num, vg, next_free_dof);
1341 
1342                   next_free_dof += (n_vars_in_group*
1343                                     node.n_comp_group(sys_num,vg));
1344                 }
1345             }
1346 
1347           // Now number the element DOFS
1348           if (elem->n_comp_group(sys_num,vg) > 0)
1349             {
1350               libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
1351                                        DofObject::invalid_id);
1352 
1353               elem->set_vg_dof_base(sys_num,
1354                                     vg,
1355                                     next_free_dof);
1356 
1357               next_free_dof += (n_vars_in_group*
1358                                 elem->n_comp_group(sys_num,vg));
1359             }
1360         } // end loop on elements
1361 
1362       // we may have missed assigning DOFs to nodes that we own
1363       // but to which we have no connected elements matching our
1364       // variable restriction criterion.  this will happen, for example,
1365       // if variable V is restricted to subdomain S.  We may not own
1366       // any elements which live in S, but we may own nodes which are
1367       // *connected* to elements which do.  in this scenario these nodes
1368       // will presently have unnumbered DOFs. we need to take care of
1369       // them here since we own them and no other processor will touch them.
1370       for (auto & node : mesh.local_node_ptr_range())
1371         if (node->n_comp_group(sys_num,vg))
1372           if (node->vg_dof_base(sys_num,vg) == DofObject::invalid_id)
1373             {
1374               node->set_vg_dof_base (sys_num,
1375                                      vg,
1376                                      next_free_dof);
1377 
1378               next_free_dof += (n_vars_in_group*
1379                                 node->n_comp_group(sys_num,vg));
1380             }
1381     } // end loop on variable groups
1382 
1383   // Finally, count up the SCALAR dofs
1384   this->_n_SCALAR_dofs = 0;
1385   for (unsigned vg=0; vg<n_var_groups; vg++)
1386     {
1387       const VariableGroup & vg_description(this->variable_group(vg));
1388 
1389       if (vg_description.type().family == SCALAR)
1390         {
1391           this->_n_SCALAR_dofs += (vg_description.n_variables()*
1392                                    vg_description.type().order.get_order());
1393           continue;
1394         }
1395     }
1396 
1397   // Only increment next_free_dof if we're on the processor
1398   // that holds this SCALAR variable
1399   if (this->processor_id() == (this->n_processors()-1))
1400     next_free_dof += _n_SCALAR_dofs;
1401 
1402 #ifdef DEBUG
1403   {
1404     // Make sure we didn't miss any nodes
1405     MeshTools::libmesh_assert_valid_procids<Node>(mesh);
1406 
1407     for (auto & node : mesh.local_node_ptr_range())
1408       {
1409         unsigned int n_var_g = node->n_var_groups(this->sys_number());
1410         for (unsigned int vg=0; vg != n_var_g; ++vg)
1411           {
1412             unsigned int n_comp_g =
1413               node->n_comp_group(this->sys_number(), vg);
1414             dof_id_type my_first_dof = n_comp_g ?
1415               node->vg_dof_base(this->sys_number(), vg) : 0;
1416             libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
1417           }
1418       }
1419   }
1420 #endif // DEBUG
1421 }
1422 
1423 
1424 
1425 void
1426 DofMap::
merge_ghost_functor_outputs(GhostingFunctor::map_type & elements_to_ghost,std::set<CouplingMatrix * > & temporary_coupling_matrices,const std::set<GhostingFunctor * >::iterator & gf_begin,const std::set<GhostingFunctor * >::iterator & gf_end,const MeshBase::const_element_iterator & elems_begin,const MeshBase::const_element_iterator & elems_end,processor_id_type p)1427 merge_ghost_functor_outputs(GhostingFunctor::map_type & elements_to_ghost,
1428                             std::set<CouplingMatrix *> & temporary_coupling_matrices,
1429                             const std::set<GhostingFunctor *>::iterator & gf_begin,
1430                             const std::set<GhostingFunctor *>::iterator & gf_end,
1431                             const MeshBase::const_element_iterator & elems_begin,
1432                             const MeshBase::const_element_iterator & elems_end,
1433                             processor_id_type p)
1434 {
1435   for (const auto & gf : as_range(gf_begin, gf_end))
1436     {
1437       GhostingFunctor::map_type more_elements_to_ghost;
1438 
1439       libmesh_assert(gf);
1440       (*gf)(elems_begin, elems_end, p, more_elements_to_ghost);
1441 
1442       for (const auto & pr : more_elements_to_ghost)
1443         {
1444           GhostingFunctor::map_type::iterator existing_it =
1445             elements_to_ghost.find (pr.first);
1446           if (existing_it == elements_to_ghost.end())
1447             elements_to_ghost.insert(pr);
1448           else
1449             {
1450               if (existing_it->second)
1451                 {
1452                   if (pr.second)
1453                     {
1454                       // If this isn't already a temporary
1455                       // then we need to make one so we'll
1456                       // have a non-const matrix to merge
1457                       if (temporary_coupling_matrices.empty() ||
1458                           temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second)) == temporary_coupling_matrices.end())
1459                         {
1460                           CouplingMatrix * cm = new CouplingMatrix(*existing_it->second);
1461                           temporary_coupling_matrices.insert(cm);
1462                           existing_it->second = cm;
1463                         }
1464                       const_cast<CouplingMatrix &>(*existing_it->second) &= *pr.second;
1465                     }
1466                   else
1467                     {
1468                       // Any existing_it matrix merged with a full
1469                       // matrix (symbolized as nullptr) gives another
1470                       // full matrix (symbolizable as nullptr).
1471 
1472                       // So if existing_it->second is a temporary then
1473                       // we don't need it anymore; we might as well
1474                       // remove it to keep the set of temporaries
1475                       // small.
1476                       std::set<CouplingMatrix *>::iterator temp_it =
1477                         temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second));
1478                       if (temp_it != temporary_coupling_matrices.end())
1479                       {
1480                         delete *temp_it;
1481                         temporary_coupling_matrices.erase(temp_it);
1482                       }
1483 
1484                       existing_it->second = nullptr;
1485                     }
1486                 }
1487               // else we have a nullptr already, then we have a full
1488               // coupling matrix, already, and merging with anything
1489               // else won't change that, so we're done.
1490             }
1491         }
1492     }
1493 }
1494 
1495 
1496 
add_neighbors_to_send_list(MeshBase & mesh)1497 void DofMap::add_neighbors_to_send_list(MeshBase & mesh)
1498 {
1499   LOG_SCOPE("add_neighbors_to_send_list()", "DofMap");
1500 
1501   // Return immediately if there's no ghost data
1502   if (this->n_processors() == 1)
1503     return;
1504 
1505   const unsigned int n_var  = this->n_variables();
1506 
1507   MeshBase::const_element_iterator       local_elem_it
1508     = mesh.active_local_elements_begin();
1509   const MeshBase::const_element_iterator local_elem_end
1510     = mesh.active_local_elements_end();
1511 
1512   GhostingFunctor::map_type elements_to_send;
1513 
1514   // Man, I wish we had guaranteed unique_ptr availability...
1515   std::set<CouplingMatrix *> temporary_coupling_matrices;
1516 
1517   // We need to add dofs to the send list if they've been directly
1518   // requested by an algebraic ghosting functor or they've been
1519   // indirectly requested by a coupling functor.
1520   this->merge_ghost_functor_outputs(elements_to_send,
1521                                     temporary_coupling_matrices,
1522                                     this->algebraic_ghosting_functors_begin(),
1523                                     this->algebraic_ghosting_functors_end(),
1524                                     local_elem_it, local_elem_end, mesh.processor_id());
1525 
1526   this->merge_ghost_functor_outputs(elements_to_send,
1527                                     temporary_coupling_matrices,
1528                                     this->coupling_functors_begin(),
1529                                     this->coupling_functors_end(),
1530                                     local_elem_it, local_elem_end, mesh.processor_id());
1531 
1532   // Making a list of non-zero coupling matrix columns is an
1533   // O(N_var^2) operation.  We cache it so we only have to do it once
1534   // per CouplingMatrix and not once per element.
1535   std::map<const CouplingMatrix *, std::vector<unsigned int>>
1536     column_variable_lists;
1537 
1538   for (auto & pr : elements_to_send)
1539     {
1540       const Elem * const partner = pr.first;
1541 
1542       // We asked ghosting functors not to give us local elements
1543       libmesh_assert_not_equal_to
1544         (partner->processor_id(), this->processor_id());
1545 
1546       const CouplingMatrix * ghost_coupling = pr.second;
1547 
1548       // Loop over any present coupling matrix column variables if we
1549       // have a coupling matrix, or just add all variables to
1550       // send_list if not.
1551       if (ghost_coupling)
1552         {
1553           libmesh_assert_equal_to (ghost_coupling->size(), n_var);
1554 
1555           // Try to find a cached list of column variables.
1556           std::map<const CouplingMatrix *, std::vector<unsigned int>>::const_iterator
1557             column_variable_list = column_variable_lists.find(ghost_coupling);
1558 
1559           // If we didn't find it, then we need to create it.
1560           if (column_variable_list == column_variable_lists.end())
1561             {
1562               auto inserted_variable_list_pair =
1563                 column_variable_lists.emplace(ghost_coupling, std::vector<unsigned int>());
1564               column_variable_list = inserted_variable_list_pair.first;
1565 
1566               std::vector<unsigned int> & new_variable_list =
1567                 inserted_variable_list_pair.first->second;
1568 
1569               std::vector<unsigned char> has_variable(n_var, false);
1570 
1571               for (unsigned int vi = 0; vi != n_var; ++vi)
1572                 {
1573                   ConstCouplingRow ccr(vi, *ghost_coupling);
1574 
1575                   for (const auto & vj : ccr)
1576                     has_variable[vj] = true;
1577                 }
1578               for (unsigned int vj = 0; vj != n_var; ++vj)
1579                 {
1580                   if (has_variable[vj])
1581                     new_variable_list.push_back(vj);
1582                 }
1583             }
1584 
1585           const std::vector<unsigned int> & variable_list =
1586             column_variable_list->second;
1587 
1588           for (const auto & vj : variable_list)
1589             {
1590               std::vector<dof_id_type> di;
1591               this->dof_indices (partner, di, vj);
1592 
1593               // Insert the remote DOF indices into the send list
1594               for (auto d : di)
1595                 if (d != DofObject::invalid_id &&
1596                     !this->local_index(d))
1597                   {
1598                     libmesh_assert_less(d, this->n_dofs());
1599                     _send_list.push_back(d);
1600                   }
1601             }
1602         }
1603       else
1604         {
1605           std::vector<dof_id_type> di;
1606           this->dof_indices (partner, di);
1607 
1608           // Insert the remote DOF indices into the send list
1609           for (const auto & dof : di)
1610             if (dof != DofObject::invalid_id &&
1611                 !this->local_index(dof))
1612               {
1613                 libmesh_assert_less(dof, this->n_dofs());
1614                 _send_list.push_back(dof);
1615               }
1616         }
1617 
1618     }
1619 
1620   // We're now done with any merged coupling matrices we had to create.
1621   for (auto & mat : temporary_coupling_matrices)
1622     delete mat;
1623 
1624   //-------------------------------------------------------------------------
1625   // Our coupling functors added dofs from neighboring elements to the
1626   // send list, but we may still need to add non-local dofs from local
1627   // elements.
1628   //-------------------------------------------------------------------------
1629 
1630   // Loop over the active local elements, adding all active elements
1631   // that neighbor an active local element to the send list.
1632   for ( ; local_elem_it != local_elem_end; ++local_elem_it)
1633     {
1634       const Elem * elem = *local_elem_it;
1635 
1636       std::vector<dof_id_type> di;
1637       this->dof_indices (elem, di);
1638 
1639       // Insert the remote DOF indices into the send list
1640       for (const auto & dof : di)
1641         if (dof != DofObject::invalid_id &&
1642             !this->local_index(dof))
1643           {
1644             libmesh_assert_less(dof, this->n_dofs());
1645             _send_list.push_back(dof);
1646           }
1647     }
1648 }
1649 
1650 
1651 
prepare_send_list()1652 void DofMap::prepare_send_list ()
1653 {
1654   LOG_SCOPE("prepare_send_list()", "DofMap");
1655 
1656   // Return immediately if there's no ghost data
1657   if (this->n_processors() == 1)
1658     return;
1659 
1660   // Check to see if we have any extra stuff to add to the send_list
1661   if (_extra_send_list_function)
1662     {
1663       if (_augment_send_list)
1664         {
1665           libmesh_here();
1666           libMesh::out << "WARNING:  You have specified both an extra send list function and object.\n"
1667                        << "          Are you sure this is what you meant to do??"
1668                        << std::endl;
1669         }
1670 
1671       _extra_send_list_function(_send_list, _extra_send_list_context);
1672     }
1673 
1674   if (_augment_send_list)
1675     _augment_send_list->augment_send_list (_send_list);
1676 
1677   // First sort the send list.  After this
1678   // duplicated elements will be adjacent in the
1679   // vector
1680   std::sort(_send_list.begin(), _send_list.end());
1681 
1682   // Now use std::unique to remove duplicate entries
1683   std::vector<dof_id_type>::iterator new_end =
1684     std::unique (_send_list.begin(), _send_list.end());
1685 
1686   // Remove the end of the send_list.  Use the "swap trick"
1687   // from Effective STL
1688   std::vector<dof_id_type> (_send_list.begin(), new_end).swap (_send_list);
1689 
1690   // Make sure the send list has nothing invalid in it.
1691   libmesh_assert(_send_list.empty() || _send_list.back() < this->n_dofs());
1692 }
1693 
reinit_send_list(MeshBase & mesh)1694 void DofMap::reinit_send_list (MeshBase & mesh)
1695 {
1696   this->clear_send_list();
1697   this->add_neighbors_to_send_list(mesh);
1698 
1699 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1700   // This is assuming that we only need to recommunicate
1701   // the constraints and no new ones have been added since
1702   // a previous call to reinit_constraints.
1703   this->process_constraints(mesh);
1704 #endif
1705   this->prepare_send_list();
1706 }
1707 
set_implicit_neighbor_dofs(bool implicit_neighbor_dofs)1708 void DofMap::set_implicit_neighbor_dofs(bool implicit_neighbor_dofs)
1709 {
1710   _implicit_neighbor_dofs_initialized = true;
1711   _implicit_neighbor_dofs = implicit_neighbor_dofs;
1712 }
1713 
1714 
use_coupled_neighbor_dofs(const MeshBase & mesh)1715 bool DofMap::use_coupled_neighbor_dofs(const MeshBase & mesh) const
1716 {
1717   // If we were asked on the command line, then we need to
1718   // include sensitivities between neighbor degrees of freedom
1719   bool implicit_neighbor_dofs =
1720     libMesh::on_command_line ("--implicit-neighbor-dofs");
1721 
1722   // If the user specifies --implicit-neighbor-dofs 0, then
1723   // presumably he knows what he is doing and we won't try to
1724   // automatically turn it on even when all the variables are
1725   // discontinuous.
1726   if (implicit_neighbor_dofs)
1727     {
1728       // No flag provided defaults to 'true'
1729       int flag = 1;
1730       flag = libMesh::command_line_next ("--implicit-neighbor-dofs", flag);
1731 
1732       if (!flag)
1733         {
1734           // The user said --implicit-neighbor-dofs 0, so he knows
1735           // what he is doing and really doesn't want it.
1736           return false;
1737         }
1738     }
1739 
1740   // Possibly override the commandline option, if set_implicit_neighbor_dofs
1741   // has been called.
1742   if (_implicit_neighbor_dofs_initialized)
1743     {
1744       implicit_neighbor_dofs = _implicit_neighbor_dofs;
1745 
1746       // Again, if the user explicitly says implicit_neighbor_dofs = false,
1747       // then we return here.
1748       if (!implicit_neighbor_dofs)
1749         return false;
1750     }
1751 
1752   // Look at all the variables in this system.  If every one is
1753   // discontinuous then the user must be doing DG/FVM, so be nice
1754   // and force implicit_neighbor_dofs=true.
1755   {
1756     bool all_discontinuous_dofs = true;
1757 
1758     for (auto var : make_range(this->n_variables()))
1759       if (FEAbstract::build (mesh.mesh_dimension(),
1760                              this->variable_type(var))->get_continuity() !=  DISCONTINUOUS)
1761         all_discontinuous_dofs = false;
1762 
1763     if (all_discontinuous_dofs)
1764       implicit_neighbor_dofs = true;
1765   }
1766 
1767   return implicit_neighbor_dofs;
1768 }
1769 
1770 
1771 
compute_sparsity(const MeshBase & mesh)1772 void DofMap::compute_sparsity(const MeshBase & mesh)
1773 {
1774   _sp = this->build_sparsity(mesh);
1775 
1776   // It is possible that some \p SparseMatrix implementations want to
1777   // see the sparsity pattern before we throw it away.  If so, we
1778   // share a view of its arrays, and we pass it in to the matrices.
1779   if (need_full_sparsity_pattern)
1780     {
1781       _n_nz = &_sp->n_nz;
1782       _n_oz = &_sp->n_oz;
1783 
1784       for (const auto & mat : _matrices)
1785         mat->update_sparsity_pattern (_sp->sparsity_pattern);
1786     }
1787   // If we don't need the full sparsity pattern anymore, steal the
1788   // arrays we do need and free the rest of the memory
1789   else
1790     {
1791       if (!_n_nz)
1792         _n_nz = new std::vector<dof_id_type>();
1793       _n_nz->swap(_sp->n_nz);
1794       if (!_n_oz)
1795         _n_oz = new std::vector<dof_id_type>();
1796       _n_oz->swap(_sp->n_oz);
1797 
1798       _sp.reset();
1799     }
1800 }
1801 
1802 
1803 
clear_sparsity()1804 void DofMap::clear_sparsity()
1805 {
1806   if (need_full_sparsity_pattern)
1807     {
1808       libmesh_assert(_sp.get());
1809       libmesh_assert(!_n_nz || _n_nz == &_sp->n_nz);
1810       libmesh_assert(!_n_oz || _n_oz == &_sp->n_oz);
1811       _sp.reset();
1812     }
1813   else
1814     {
1815       libmesh_assert(!_sp.get());
1816       delete _n_nz;
1817       delete _n_oz;
1818     }
1819   _n_nz = nullptr;
1820   _n_oz = nullptr;
1821 }
1822 
1823 
1824 
remove_default_ghosting()1825 void DofMap::remove_default_ghosting()
1826 {
1827   this->remove_coupling_functor(this->default_coupling());
1828   this->remove_algebraic_ghosting_functor(this->default_algebraic_ghosting());
1829 }
1830 
1831 
1832 
add_default_ghosting()1833 void DofMap::add_default_ghosting()
1834 {
1835   this->add_coupling_functor(this->default_coupling());
1836   this->add_algebraic_ghosting_functor(this->default_algebraic_ghosting());
1837 }
1838 
1839 
1840 
1841 void
add_coupling_functor(GhostingFunctor & coupling_functor,bool to_mesh)1842 DofMap::add_coupling_functor(GhostingFunctor & coupling_functor,
1843                              bool to_mesh)
1844 {
1845   _coupling_functors.insert(&coupling_functor);
1846   coupling_functor.set_mesh(&_mesh);
1847   if (to_mesh)
1848     _mesh.add_ghosting_functor(coupling_functor);
1849 }
1850 
1851 
1852 
1853 void
remove_coupling_functor(GhostingFunctor & coupling_functor)1854 DofMap::remove_coupling_functor(GhostingFunctor & coupling_functor)
1855 {
1856   _coupling_functors.erase(&coupling_functor);
1857   _mesh.remove_ghosting_functor(coupling_functor);
1858 
1859   auto it = _shared_functors.find(&coupling_functor);
1860   if (it != _shared_functors.end())
1861     _shared_functors.erase(it);
1862 }
1863 
1864 
1865 
1866 void
add_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor,bool to_mesh)1867 DofMap::add_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor,
1868                                        bool to_mesh)
1869 {
1870   _algebraic_ghosting_functors.insert(&evaluable_functor);
1871   evaluable_functor.set_mesh(&_mesh);
1872   if (to_mesh)
1873     _mesh.add_ghosting_functor(evaluable_functor);
1874 }
1875 
1876 
1877 
1878 void
remove_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor)1879 DofMap::remove_algebraic_ghosting_functor(GhostingFunctor & evaluable_functor)
1880 {
1881   _algebraic_ghosting_functors.erase(&evaluable_functor);
1882   _mesh.remove_ghosting_functor(evaluable_functor);
1883 
1884   auto it = _shared_functors.find(&evaluable_functor);
1885   if (it != _shared_functors.end())
1886     _shared_functors.erase(it);
1887 }
1888 
1889 
1890 
extract_local_vector(const NumericVector<Number> & Ug,const std::vector<dof_id_type> & dof_indices_in,DenseVectorBase<Number> & Ue)1891 void DofMap::extract_local_vector (const NumericVector<Number> & Ug,
1892                                    const std::vector<dof_id_type> & dof_indices_in,
1893                                    DenseVectorBase<Number> & Ue) const
1894 {
1895   const unsigned int n_original_dofs = dof_indices_in.size();
1896 
1897 #ifdef LIBMESH_ENABLE_AMR
1898 
1899   // Trivial mapping
1900   libmesh_assert_equal_to (dof_indices_in.size(), Ue.size());
1901   bool has_constrained_dofs = false;
1902 
1903   for (unsigned int il=0; il != n_original_dofs; ++il)
1904     {
1905       const dof_id_type ig = dof_indices_in[il];
1906 
1907       if (this->is_constrained_dof (ig)) has_constrained_dofs = true;
1908 
1909       libmesh_assert_less (ig, Ug.size());
1910 
1911       Ue.el(il) = Ug(ig);
1912     }
1913 
1914   // If the element has any constrained DOFs then we need
1915   // to account for them in the mapping.  This will handle
1916   // the case that the input vector is not constrained.
1917   if (has_constrained_dofs)
1918     {
1919       // Copy the input DOF indices.
1920       std::vector<dof_id_type> constrained_dof_indices(dof_indices_in);
1921 
1922       DenseMatrix<Number> C;
1923       DenseVector<Number> H;
1924 
1925       this->build_constraint_matrix_and_vector (C, H, constrained_dof_indices);
1926 
1927       libmesh_assert_equal_to (dof_indices_in.size(), C.m());
1928       libmesh_assert_equal_to (constrained_dof_indices.size(), C.n());
1929 
1930       // zero-out Ue
1931       Ue.zero();
1932 
1933       // compute Ue = C Ug, with proper mapping.
1934       for (unsigned int i=0; i != n_original_dofs; i++)
1935         {
1936           Ue.el(i) = H(i);
1937 
1938           const unsigned int n_constrained =
1939             cast_int<unsigned int>(constrained_dof_indices.size());
1940           for (unsigned int j=0; j<n_constrained; j++)
1941             {
1942               const dof_id_type jg = constrained_dof_indices[j];
1943 
1944               //          If Ug is a serial or ghosted vector, then this assert is
1945               //          overzealous.  If Ug is a parallel vector, then this assert
1946               //          is redundant.
1947               //    libmesh_assert ((jg >= Ug.first_local_index()) &&
1948               //    (jg <  Ug.last_local_index()));
1949 
1950               Ue.el(i) += C(i,j)*Ug(jg);
1951             }
1952         }
1953     }
1954 
1955 #else
1956 
1957   // Trivial mapping
1958 
1959   libmesh_assert_equal_to (n_original_dofs, Ue.size());
1960 
1961   for (unsigned int il=0; il<n_original_dofs; il++)
1962     {
1963       const dof_id_type ig = dof_indices_in[il];
1964 
1965       libmesh_assert ((ig >= Ug.first_local_index()) && (ig <  Ug.last_local_index()));
1966 
1967       Ue.el(il) = Ug(ig);
1968     }
1969 
1970 #endif
1971 }
1972 
dof_indices(const Elem * const elem,std::vector<dof_id_type> & di)1973 void DofMap::dof_indices (const Elem * const elem,
1974                           std::vector<dof_id_type> & di) const
1975 {
1976   // We now allow elem==nullptr to request just SCALAR dofs
1977   // libmesh_assert(elem);
1978 
1979   // If we are asking for current indices on an element, it ought to
1980   // be an active element (or a Side proxy, which also thinks it's
1981   // active)
1982   libmesh_assert(!elem || elem->active());
1983 
1984   LOG_SCOPE("dof_indices()", "DofMap");
1985 
1986   // Clear the DOF indices vector
1987   di.clear();
1988 
1989   const unsigned int n_var_groups  = this->n_variable_groups();
1990 
1991 #ifdef DEBUG
1992   // Check that sizes match in DEBUG mode
1993   std::size_t tot_size = 0;
1994 #endif
1995 
1996   if (elem && elem->type() == TRI3SUBDIVISION)
1997     {
1998       // Subdivision surface FE require the 1-ring around elem
1999       const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2000 
2001       // Ghost subdivision elements have no real dofs
2002       if (!sd_elem->is_ghost())
2003         {
2004           // Determine the nodes contributing to element elem
2005           std::vector<const Node *> elem_nodes;
2006           MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2007 
2008           // Get the dof numbers
2009           for (unsigned int vg=0; vg<n_var_groups; vg++)
2010             {
2011               const VariableGroup & var = this->variable_group(vg);
2012               const unsigned int vars_in_group = var.n_variables();
2013 
2014               if (var.type().family == SCALAR &&
2015                   var.active_on_subdomain(elem->subdomain_id()))
2016                 {
2017                   for (unsigned int vig=0; vig != vars_in_group; ++vig)
2018                     {
2019 #ifdef DEBUG
2020                       tot_size += var.type().order;
2021 #endif
2022                       std::vector<dof_id_type> di_new;
2023                       this->SCALAR_dof_indices(di_new,var.number(vig));
2024                       di.insert( di.end(), di_new.begin(), di_new.end());
2025                     }
2026                 }
2027               else
2028                 for (unsigned int vig=0; vig != vars_in_group; ++vig)
2029                   {
2030                     _dof_indices(*elem, elem->p_level(), di, vg, vig,
2031                                  elem_nodes.data(),
2032                                  cast_int<unsigned int>(elem_nodes.size())
2033 #ifdef DEBUG
2034                                  , var.number(vig), tot_size
2035 #endif
2036                                  );
2037                   }
2038             }
2039         }
2040 
2041       return;
2042     }
2043 
2044   // Get the dof numbers for each variable
2045   const unsigned int n_nodes = elem ? elem->n_nodes() : 0;
2046   for (unsigned int vg=0; vg<n_var_groups; vg++)
2047     {
2048       const VariableGroup & var = this->variable_group(vg);
2049       const unsigned int vars_in_group = var.n_variables();
2050 
2051       if (var.type().family == SCALAR &&
2052           (!elem ||
2053            var.active_on_subdomain(elem->subdomain_id())))
2054         {
2055           for (unsigned int vig=0; vig != vars_in_group; ++vig)
2056             {
2057 #ifdef DEBUG
2058               tot_size += var.type().order;
2059 #endif
2060               std::vector<dof_id_type> di_new;
2061               this->SCALAR_dof_indices(di_new,var.number(vig));
2062               di.insert( di.end(), di_new.begin(), di_new.end());
2063             }
2064         }
2065       else if (elem)
2066         for (unsigned int vig=0; vig != vars_in_group; ++vig)
2067           {
2068             _dof_indices(*elem, elem->p_level(), di, vg, vig,
2069                          elem->get_nodes(), n_nodes
2070 #ifdef DEBUG
2071                          , var.number(vig), tot_size
2072 #endif
2073                      );
2074           }
2075     }
2076 
2077 #ifdef DEBUG
2078   libmesh_assert_equal_to (tot_size, di.size());
2079 #endif
2080 }
2081 
2082 
dof_indices(const Elem * const elem,std::vector<dof_id_type> & di,const unsigned int vn,int p_level)2083 void DofMap::dof_indices (const Elem * const elem,
2084                           std::vector<dof_id_type> & di,
2085                           const unsigned int vn,
2086                           int p_level) const
2087 {
2088   // We now allow elem==nullptr to request just SCALAR dofs
2089   // libmesh_assert(elem);
2090 
2091   LOG_SCOPE("dof_indices()", "DofMap");
2092 
2093   // Clear the DOF indices vector
2094   di.clear();
2095 
2096   // Use the default p refinement level?
2097   if (p_level == -12345)
2098     p_level = elem ? elem->p_level() : 0;
2099 
2100   const unsigned int vg = this->_variable_group_numbers[vn];
2101   const VariableGroup & var = this->variable_group(vg);
2102   const unsigned int vig = vn - var.number();
2103 
2104 #ifdef DEBUG
2105   // Check that sizes match in DEBUG mode
2106   std::size_t tot_size = 0;
2107 #endif
2108 
2109   if (elem && elem->type() == TRI3SUBDIVISION)
2110     {
2111       // Subdivision surface FE require the 1-ring around elem
2112       const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2113 
2114       // Ghost subdivision elements have no real dofs
2115       if (!sd_elem->is_ghost())
2116         {
2117           // Determine the nodes contributing to element elem
2118           std::vector<const Node *> elem_nodes;
2119           MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2120 
2121           _dof_indices(*elem, p_level, di, vg, vig, elem_nodes.data(),
2122                        cast_int<unsigned int>(elem_nodes.size())
2123 #ifdef DEBUG
2124                        , vn, tot_size
2125 #endif
2126                        );
2127         }
2128 
2129       return;
2130     }
2131 
2132   // Get the dof numbers
2133   if (var.type().family == SCALAR &&
2134       (!elem ||
2135        var.active_on_subdomain(elem->subdomain_id())))
2136     {
2137 #ifdef DEBUG
2138       tot_size += var.type().order;
2139 #endif
2140       std::vector<dof_id_type> di_new;
2141       this->SCALAR_dof_indices(di_new,vn);
2142       di.insert( di.end(), di_new.begin(), di_new.end());
2143     }
2144   else if (elem)
2145     _dof_indices(*elem, p_level, di, vg, vig, elem->get_nodes(),
2146                  elem->n_nodes()
2147 #ifdef DEBUG
2148                  , vn, tot_size
2149 #endif
2150                  );
2151 
2152 #ifdef DEBUG
2153   libmesh_assert_equal_to (tot_size, di.size());
2154 #endif
2155 }
2156 
2157 
dof_indices(const Node * const node,std::vector<dof_id_type> & di)2158 void DofMap::dof_indices (const Node * const node,
2159                           std::vector<dof_id_type> & di) const
2160 {
2161   // We allow node==nullptr to request just SCALAR dofs
2162   // libmesh_assert(elem);
2163 
2164   LOG_SCOPE("dof_indices(Node)", "DofMap");
2165 
2166   // Clear the DOF indices vector
2167   di.clear();
2168 
2169   const unsigned int n_var_groups  = this->n_variable_groups();
2170   const unsigned int sys_num = this->sys_number();
2171 
2172   // Get the dof numbers
2173   for (unsigned int vg=0; vg<n_var_groups; vg++)
2174     {
2175       const VariableGroup & var = this->variable_group(vg);
2176       const unsigned int vars_in_group = var.n_variables();
2177 
2178       if (var.type().family == SCALAR)
2179         {
2180           for (unsigned int vig=0; vig != vars_in_group; ++vig)
2181             {
2182               std::vector<dof_id_type> di_new;
2183               this->SCALAR_dof_indices(di_new,var.number(vig));
2184               di.insert( di.end(), di_new.begin(), di_new.end());
2185             }
2186         }
2187       else
2188         {
2189           const int n_comp = node->n_comp_group(sys_num,vg);
2190           for (unsigned int vig=0; vig != vars_in_group; ++vig)
2191             {
2192               for (int i=0; i != n_comp; ++i)
2193                 {
2194                   const dof_id_type d =
2195                     node->dof_number(sys_num, vg, vig, i, n_comp);
2196                   libmesh_assert_not_equal_to
2197                     (d, DofObject::invalid_id);
2198                   di.push_back(d);
2199                 }
2200             }
2201         }
2202     }
2203 }
2204 
2205 
dof_indices(const Node * const node,std::vector<dof_id_type> & di,const unsigned int vn)2206 void DofMap::dof_indices (const Node * const node,
2207                           std::vector<dof_id_type> & di,
2208                           const unsigned int vn) const
2209 {
2210   if (vn == libMesh::invalid_uint)
2211     {
2212       this->dof_indices(node, di);
2213       return;
2214     }
2215 
2216   // We allow node==nullptr to request just SCALAR dofs
2217   // libmesh_assert(elem);
2218 
2219   LOG_SCOPE("dof_indices(Node)", "DofMap");
2220 
2221   // Clear the DOF indices vector
2222   di.clear();
2223 
2224   const unsigned int sys_num = this->sys_number();
2225 
2226   // Get the dof numbers
2227   const unsigned int vg = this->_variable_group_numbers[vn];
2228   const VariableGroup & var = this->variable_group(vg);
2229 
2230   if (var.type().family == SCALAR)
2231     {
2232       std::vector<dof_id_type> di_new;
2233       this->SCALAR_dof_indices(di_new,vn);
2234       di.insert( di.end(), di_new.begin(), di_new.end());
2235     }
2236   else
2237     {
2238       const unsigned int vig = vn - var.number();
2239       const int n_comp = node->n_comp_group(sys_num,vg);
2240       for (int i=0; i != n_comp; ++i)
2241         {
2242           const dof_id_type d =
2243             node->dof_number(sys_num, vg, vig, i, n_comp);
2244           libmesh_assert_not_equal_to
2245             (d, DofObject::invalid_id);
2246           di.push_back(d);
2247         }
2248     }
2249 }
2250 
2251 
dof_indices(const Elem & elem,unsigned int n,std::vector<dof_id_type> & di,const unsigned int vn)2252 void DofMap::dof_indices (const Elem & elem,
2253                           unsigned int n,
2254                           std::vector<dof_id_type> & di,
2255                           const unsigned int vn) const
2256 {
2257   this->_node_dof_indices(elem, n, elem.node_ref(n), di, vn);
2258 }
2259 
2260 
2261 
2262 #ifdef LIBMESH_ENABLE_AMR
2263 
old_dof_indices(const Elem & elem,unsigned int n,std::vector<dof_id_type> & di,const unsigned int vn)2264 void DofMap::old_dof_indices (const Elem & elem,
2265                               unsigned int n,
2266                               std::vector<dof_id_type> & di,
2267                               const unsigned int vn) const
2268 {
2269   const DofObject * old_obj = elem.node_ref(n).old_dof_object;
2270   libmesh_assert(old_obj);
2271   this->_node_dof_indices(elem, n, *old_obj, di, vn);
2272 }
2273 
2274 #endif // LIBMESH_ENABLE_AMR
2275 
2276 
2277 
_node_dof_indices(const Elem & elem,unsigned int n,const DofObject & obj,std::vector<dof_id_type> & di,const unsigned int vn)2278 void DofMap::_node_dof_indices (const Elem & elem,
2279                                 unsigned int n,
2280                                 const DofObject & obj,
2281                                 std::vector<dof_id_type> & di,
2282                                 const unsigned int vn) const
2283 {
2284   // Half of this is a cut and paste of _dof_indices code below, but
2285   // duplication actually seems cleaner than creating a helper
2286   // function with a million arguments and hoping the compiler inlines
2287   // it properly into one of our most highly trafficked functions.
2288 
2289   LOG_SCOPE("_node_dof_indices()", "DofMap");
2290 
2291   const unsigned int sys_num = this->sys_number();
2292   const std::pair<unsigned int, unsigned int>
2293     vg_and_offset = obj.var_to_vg_and_offset(sys_num,vn);
2294   const unsigned int vg = vg_and_offset.first;
2295   const unsigned int vig = vg_and_offset.second;
2296   const unsigned int n_comp = obj.n_comp_group(sys_num,vg);
2297 
2298   const VariableGroup & var = this->variable_group(vg);
2299   FEType fe_type = var.type();
2300   const bool extra_hanging_dofs =
2301     FEInterface::extra_hanging_dofs(fe_type);
2302 
2303   // There is a potential problem with h refinement.  Imagine a
2304   // quad9 that has a linear FE on it.  Then, on the hanging side,
2305   // it can falsely identify a DOF at the mid-edge node. This is why
2306   // we go through FEInterface instead of obj->n_comp() directly.
2307   const unsigned int nc =
2308     FEInterface::n_dofs_at_node(fe_type, &elem, n);
2309 
2310   // If this is a non-vertex on a hanging node with extra
2311   // degrees of freedom, we use the non-vertex dofs (which
2312   // come in reverse order starting from the end, to
2313   // simplify p refinement)
2314   if (extra_hanging_dofs && nc && !elem.is_vertex(n))
2315     {
2316       const int dof_offset = n_comp - nc;
2317 
2318       // We should never have fewer dofs than necessary on a
2319       // node unless we're getting indices on a parent element,
2320       // and we should never need the indices on such a node
2321       if (dof_offset < 0)
2322         {
2323           libmesh_assert(!elem.active());
2324           di.resize(di.size() + nc, DofObject::invalid_id);
2325         }
2326       else
2327         for (unsigned int i = dof_offset; i != n_comp; ++i)
2328           {
2329             const dof_id_type d =
2330               obj.dof_number(sys_num, vg, vig, i, n_comp);
2331             libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2332             di.push_back(d);
2333           }
2334     }
2335   // If this is a vertex or an element without extra hanging
2336   // dofs, our dofs come in forward order coming from the
2337   // beginning.  But we still might not have all those dofs, in cases
2338   // where a subdomain-restricted variable just had its subdomain
2339   // expanded.
2340   else
2341     {
2342       const unsigned int good_nc =
2343         std::min(static_cast<unsigned int>(n_comp), nc);
2344       for (unsigned int i=0; i != good_nc; ++i)
2345         {
2346           const dof_id_type d =
2347             obj.dof_number(sys_num, vg, vig, i, n_comp);
2348           libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2349           di.push_back(d);
2350         }
2351       for (unsigned int i=good_nc; i != nc; ++i)
2352         di.push_back(DofObject::invalid_id);
2353     }
2354 }
2355 
2356 
2357 
_dof_indices(const Elem & elem,int p_level,std::vector<dof_id_type> & di,const unsigned int vg,const unsigned int vig,const Node * const * nodes,unsigned int n_nodes,const unsigned int v,std::size_t & tot_size)2358 void DofMap::_dof_indices (const Elem & elem,
2359                            int p_level,
2360                            std::vector<dof_id_type> & di,
2361                            const unsigned int vg,
2362                            const unsigned int vig,
2363                            const Node * const * nodes,
2364                            unsigned int       n_nodes
2365 #ifdef DEBUG
2366                            ,
2367                            const unsigned int v,
2368                            std::size_t & tot_size
2369 #endif
2370                            ) const
2371 {
2372   const VariableGroup & var = this->variable_group(vg);
2373 
2374   if (var.active_on_subdomain(elem.subdomain_id()))
2375     {
2376       const ElemType type        = elem.type();
2377       const unsigned int sys_num = this->sys_number();
2378 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2379       const bool is_inf          = elem.infinite();
2380 #endif
2381 
2382       const bool extra_hanging_dofs =
2383         FEInterface::extra_hanging_dofs(var.type());
2384 
2385       FEType fe_type = var.type();
2386 
2387 #ifdef DEBUG
2388       // The number of dofs per element is non-static for subdivision FE
2389       if (var.type().family == SUBDIVISION)
2390         tot_size += n_nodes;
2391       else
2392         // FIXME: Is the passed-in p_level just elem.p_level()? If so,
2393         // this seems redundant.
2394         tot_size += FEInterface::n_dofs(fe_type, p_level, &elem);
2395 #endif
2396 
2397       // The total Order is not required when getting the function
2398       // pointer, it is only needed when the function is called (see
2399       // below).
2400       const FEInterface::n_dofs_at_node_ptr ndan =
2401         FEInterface::n_dofs_at_node_function(fe_type, &elem);
2402 
2403       // Get the node-based DOF numbers
2404       for (unsigned int n=0; n != n_nodes; n++)
2405         {
2406           const Node & node = *nodes[n];
2407 
2408           // Cache the intermediate lookups that are common to every
2409           // component
2410 #ifdef DEBUG
2411           const std::pair<unsigned int, unsigned int>
2412             vg_and_offset = node.var_to_vg_and_offset(sys_num,v);
2413           libmesh_assert_equal_to (vg, vg_and_offset.first);
2414           libmesh_assert_equal_to (vig, vg_and_offset.second);
2415 #endif
2416           const unsigned int n_comp = node.n_comp_group(sys_num,vg);
2417 
2418           // There is a potential problem with h refinement.  Imagine a
2419           // quad9 that has a linear FE on it.  Then, on the hanging side,
2420           // it can falsely identify a DOF at the mid-edge node. This is why
2421           // we go through FEInterface instead of node.n_comp() directly.
2422           const unsigned int nc =
2423 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2424             is_inf ?
2425             FEInterface::n_dofs_at_node(fe_type, p_level, &elem, n) :
2426 #endif
2427             ndan (type, static_cast<Order>(fe_type.order + p_level), n);
2428 
2429           // If this is a non-vertex on a hanging node with extra
2430           // degrees of freedom, we use the non-vertex dofs (which
2431           // come in reverse order starting from the end, to
2432           // simplify p refinement)
2433           if (extra_hanging_dofs && !elem.is_vertex(n))
2434             {
2435               const int dof_offset = n_comp - nc;
2436 
2437               // We should never have fewer dofs than necessary on a
2438               // node unless we're getting indices on a parent element,
2439               // and we should never need the indices on such a node
2440               if (dof_offset < 0)
2441                 {
2442                   libmesh_assert(!elem.active());
2443                   di.resize(di.size() + nc, DofObject::invalid_id);
2444                 }
2445               else
2446                 for (int i=n_comp-1; i>=dof_offset; i--)
2447                   {
2448                     const dof_id_type d =
2449                       node.dof_number(sys_num, vg, vig, i, n_comp);
2450                     libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2451                     di.push_back(d);
2452                   }
2453             }
2454           // If this is a vertex or an element without extra hanging
2455           // dofs, our dofs come in forward order coming from the
2456           // beginning
2457           else
2458             {
2459               // We have a good component index only if it's being
2460               // used on this FE type (nc) *and* it's available on
2461               // this DofObject (n_comp).
2462               const unsigned int good_nc = std::min(n_comp, nc);
2463               for (unsigned int i=0; i!=good_nc; ++i)
2464                 {
2465                   const dof_id_type d =
2466                     node.dof_number(sys_num, vg, vig, i, n_comp);
2467                   libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2468                   libmesh_assert_less (d, this->n_dofs());
2469                   di.push_back(d);
2470                 }
2471 
2472               // With fewer good component indices than we need, e.g.
2473               // due to subdomain expansion, the remaining expected
2474               // indices are marked invalid.
2475               if (n_comp < nc)
2476                 for (unsigned int i=n_comp; i!=nc; ++i)
2477                   di.push_back(DofObject::invalid_id);
2478             }
2479         }
2480 
2481       // If there are any element-based DOF numbers, get them
2482       const unsigned int nc = FEInterface::n_dofs_per_elem(fe_type, p_level, &elem);
2483 
2484       // We should never have fewer dofs than necessary on an
2485       // element unless we're getting indices on a parent element
2486       // (and we should never need those indices) or off-domain for a
2487       // subdomain-restricted variable (where invalid_id is the
2488       // correct thing to return)
2489       if (nc != 0)
2490         {
2491           const unsigned int n_comp = elem.n_comp_group(sys_num,vg);
2492           if (elem.n_systems() > sys_num && nc <= n_comp)
2493             {
2494               for (unsigned int i=0; i<nc; i++)
2495                 {
2496                   const dof_id_type d =
2497                     elem.dof_number(sys_num, vg, vig, i, n_comp);
2498                   libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2499 
2500                   di.push_back(d);
2501                 }
2502             }
2503           else
2504             {
2505               libmesh_assert(!elem.active() || fe_type.family == LAGRANGE || fe_type.family == SUBDIVISION);
2506               di.resize(di.size() + nc, DofObject::invalid_id);
2507             }
2508         }
2509     }
2510 }
2511 
2512 
2513 
SCALAR_dof_indices(std::vector<dof_id_type> & di,const unsigned int vn,const bool old_dofs)2514 void DofMap::SCALAR_dof_indices (std::vector<dof_id_type> & di,
2515                                  const unsigned int vn,
2516 #ifdef LIBMESH_ENABLE_AMR
2517                                  const bool old_dofs
2518 #else
2519                                  const bool
2520 #endif
2521                                  ) const
2522 {
2523   LOG_SCOPE("SCALAR_dof_indices()", "DofMap");
2524 
2525   libmesh_assert(this->variable(vn).type().family == SCALAR);
2526 
2527 #ifdef LIBMESH_ENABLE_AMR
2528   // If we're asking for old dofs then we'd better have some
2529   if (old_dofs)
2530     libmesh_assert_greater_equal(n_old_dofs(), n_SCALAR_dofs());
2531 
2532   dof_id_type my_idx = old_dofs ?
2533     this->_first_old_scalar_df[vn] : this->_first_scalar_df[vn];
2534 #else
2535   dof_id_type my_idx = this->_first_scalar_df[vn];
2536 #endif
2537 
2538   libmesh_assert_not_equal_to(my_idx, DofObject::invalid_id);
2539 
2540   // The number of SCALAR dofs comes from the variable order
2541   const int n_dofs_vn = this->variable(vn).type().order.get_order();
2542 
2543   di.resize(n_dofs_vn);
2544   for (int i = 0; i != n_dofs_vn; ++i)
2545     di[i] = my_idx++;
2546 }
2547 
2548 
2549 
semilocal_index(dof_id_type dof_index)2550 bool DofMap::semilocal_index (dof_id_type dof_index) const
2551 {
2552   // If it's not in the local indices
2553   if (!this->local_index(dof_index))
2554     {
2555       // and if it's not in the ghost indices, then we're not
2556       // semilocal
2557       if (!std::binary_search(_send_list.begin(), _send_list.end(), dof_index))
2558         return false;
2559     }
2560 
2561   return true;
2562 }
2563 
2564 
2565 
all_semilocal_indices(const std::vector<dof_id_type> & dof_indices_in)2566 bool DofMap::all_semilocal_indices (const std::vector<dof_id_type> & dof_indices_in) const
2567 {
2568   // We're all semilocal unless we find a counterexample
2569   for (const auto & di : dof_indices_in)
2570     if (!this->semilocal_index(di))
2571       return false;
2572 
2573   return true;
2574 }
2575 
2576 
2577 
2578 template <typename DofObjectSubclass>
is_evaluable(const DofObjectSubclass & obj,unsigned int var_num)2579 bool DofMap::is_evaluable(const DofObjectSubclass & obj,
2580                           unsigned int var_num) const
2581 {
2582   // Everything is evaluable on a local object
2583   if (obj.processor_id() == this->processor_id())
2584     return true;
2585 
2586   std::vector<dof_id_type> di;
2587 
2588   if (var_num == libMesh::invalid_uint)
2589     this->dof_indices(&obj, di);
2590   else
2591     this->dof_indices(&obj, di, var_num);
2592 
2593   return this->all_semilocal_indices(di);
2594 }
2595 
2596 
2597 
2598 #ifdef LIBMESH_ENABLE_AMR
2599 
old_dof_indices(const Elem * const elem,std::vector<dof_id_type> & di,const unsigned int vn)2600 void DofMap::old_dof_indices (const Elem * const elem,
2601                               std::vector<dof_id_type> & di,
2602                               const unsigned int vn) const
2603 {
2604   LOG_SCOPE("old_dof_indices()", "DofMap");
2605 
2606   libmesh_assert(elem);
2607 
2608   const ElemType type              = elem->type();
2609   const unsigned int sys_num       = this->sys_number();
2610   const unsigned int n_var_groups  = this->n_variable_groups();
2611 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2612   const bool is_inf                = elem->infinite();
2613 #endif
2614 
2615   // If we have dof indices stored on the elem, and there's no chance
2616   // that we only have those indices because we were just p refined,
2617   // then we should have old dof indices too.
2618   libmesh_assert(!elem->has_dofs(sys_num) ||
2619                  elem->p_refinement_flag() == Elem::JUST_REFINED ||
2620                  elem->old_dof_object);
2621 
2622   // Clear the DOF indices vector.
2623   di.clear();
2624 
2625   // Determine the nodes contributing to element elem
2626   std::vector<const Node *> elem_nodes;
2627   const Node * const * nodes_ptr;
2628   unsigned int n_nodes;
2629   if (elem->type() == TRI3SUBDIVISION)
2630     {
2631       // Subdivision surface FE require the 1-ring around elem
2632       const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2633       MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2634       nodes_ptr = elem_nodes.data();
2635       n_nodes = cast_int<unsigned int>(elem_nodes.size());
2636     }
2637   else
2638     {
2639       // All other FE use only the nodes of elem itself
2640       nodes_ptr = elem->get_nodes();
2641       n_nodes = elem->n_nodes();
2642     }
2643 
2644   // Get the dof numbers
2645   for (unsigned int vg=0; vg<n_var_groups; vg++)
2646     {
2647       const VariableGroup & var = this->variable_group(vg);
2648       const unsigned int vars_in_group = var.n_variables();
2649 
2650       for (unsigned int vig=0; vig<vars_in_group; vig++)
2651         {
2652           const unsigned int v = var.number(vig);
2653           if ((vn == v) || (vn == libMesh::invalid_uint))
2654             {
2655               if (var.type().family == SCALAR &&
2656                   (!elem ||
2657                    var.active_on_subdomain(elem->subdomain_id())))
2658                 {
2659                   // We asked for this variable, so add it to the vector.
2660                   std::vector<dof_id_type> di_new;
2661                   this->SCALAR_dof_indices(di_new,v,true);
2662                   di.insert( di.end(), di_new.begin(), di_new.end());
2663                 }
2664               else
2665                 if (var.active_on_subdomain(elem->subdomain_id()))
2666                   { // Do this for all the variables if one was not specified
2667                     // or just for the specified variable
2668 
2669                     // Increase the polynomial order on p refined elements,
2670                     // but make sure you get the right polynomial order for
2671                     // the OLD degrees of freedom
2672                     int p_adjustment = 0;
2673                     if (elem->p_refinement_flag() == Elem::JUST_REFINED)
2674                       {
2675                         libmesh_assert_greater (elem->p_level(), 0);
2676                         p_adjustment = -1;
2677                       }
2678                     else if (elem->p_refinement_flag() == Elem::JUST_COARSENED)
2679                       {
2680                         p_adjustment = 1;
2681                       }
2682 
2683                     // Compute the net amount of "extra" order, including Elem::p_level()
2684                     int extra_order = elem->p_level() + p_adjustment;
2685 
2686                     FEType fe_type = var.type();
2687 
2688                     const bool extra_hanging_dofs =
2689                       FEInterface::extra_hanging_dofs(fe_type);
2690 
2691                     const FEInterface::n_dofs_at_node_ptr ndan =
2692                       FEInterface::n_dofs_at_node_function(fe_type, elem);
2693 
2694                     // Get the node-based DOF numbers
2695                     for (unsigned int n=0; n<n_nodes; n++)
2696                       {
2697                         const Node * node = nodes_ptr[n];
2698                         const DofObject * old_dof_obj = node->old_dof_object;
2699                         libmesh_assert(old_dof_obj);
2700 
2701                         // There is a potential problem with h refinement.  Imagine a
2702                         // quad9 that has a linear FE on it.  Then, on the hanging side,
2703                         // it can falsely identify a DOF at the mid-edge node. This is why
2704                         // we call FEInterface instead of node->n_comp() directly.
2705                         const unsigned int nc =
2706 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2707                           is_inf ?
2708                           FEInterface::n_dofs_at_node(var.type(), extra_order, elem, n) :
2709 #endif
2710                           ndan (type, static_cast<Order>(var.type().order + extra_order), n);
2711 
2712                         const int n_comp = old_dof_obj->n_comp_group(sys_num,vg);
2713 
2714                         // If this is a non-vertex on a hanging node with extra
2715                         // degrees of freedom, we use the non-vertex dofs (which
2716                         // come in reverse order starting from the end, to
2717                         // simplify p refinement)
2718                         if (extra_hanging_dofs && !elem->is_vertex(n))
2719                           {
2720                             const int dof_offset = n_comp - nc;
2721 
2722                             // We should never have fewer dofs than necessary on a
2723                             // node unless we're getting indices on a parent element
2724                             // or a just-coarsened element
2725                             if (dof_offset < 0)
2726                               {
2727                                 libmesh_assert(!elem->active() || elem->refinement_flag() ==
2728                                                Elem::JUST_COARSENED);
2729                                 di.resize(di.size() + nc, DofObject::invalid_id);
2730                               }
2731                             else
2732                               for (int i=n_comp-1; i>=dof_offset; i--)
2733                                 {
2734                                   const dof_id_type d =
2735                                     old_dof_obj->dof_number(sys_num, vg, vig, i, n_comp);
2736 
2737                                   // On a newly-expanded subdomain, we
2738                                   // may have some DoFs that didn't
2739                                   // exist in the old system, in which
2740                                   // case we can't assert this:
2741                                   // libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2742 
2743                                   di.push_back(d);
2744                                 }
2745                           }
2746                         // If this is a vertex or an element without extra hanging
2747                         // dofs, our dofs come in forward order coming from the
2748                         // beginning.  But we still might not have all
2749                         // those dofs on the old_dof_obj, in cases
2750                         // where a subdomain-restricted variable just
2751                         // had its subdomain expanded.
2752                         else
2753                           {
2754                             const unsigned int old_nc =
2755                               std::min(static_cast<unsigned int>(n_comp), nc);
2756                             for (unsigned int i=0; i != old_nc; ++i)
2757                               {
2758                                 const dof_id_type d =
2759                                   old_dof_obj->dof_number(sys_num, vg, vig, i, n_comp);
2760 
2761                                 libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2762 
2763                                 di.push_back(d);
2764                               }
2765                             for (unsigned int i=old_nc; i != nc; ++i)
2766                               di.push_back(DofObject::invalid_id);
2767                           }
2768                       }
2769 
2770                     // If there are any element-based DOF numbers, get them
2771                     const unsigned int nc =
2772                       FEInterface::n_dofs_per_elem(fe_type, extra_order, elem);
2773 
2774                     if (nc != 0)
2775                       {
2776                         const DofObject * old_dof_obj = elem->old_dof_object;
2777                         libmesh_assert(old_dof_obj);
2778 
2779                         const unsigned int n_comp =
2780                           old_dof_obj->n_comp_group(sys_num,vg);
2781 
2782                         if (old_dof_obj->n_systems() > sys_num &&
2783                             nc <= n_comp)
2784                           {
2785 
2786                             for (unsigned int i=0; i<nc; i++)
2787                               {
2788                                 const dof_id_type d =
2789                                   old_dof_obj->dof_number(sys_num, vg, vig, i, n_comp);
2790 
2791                                 di.push_back(d);
2792                               }
2793                           }
2794                         else
2795                           {
2796                             // We should never have fewer dofs than
2797                             // necessary on an element unless we're
2798                             // getting indices on a parent element, a
2799                             // just-coarsened element ... or a
2800                             // subdomain-restricted variable with a
2801                             // just-expanded subdomain
2802                             // libmesh_assert(!elem->active() || fe_type.family == LAGRANGE ||
2803                             //                 elem->refinement_flag() == Elem::JUST_COARSENED);
2804                             di.resize(di.size() + nc, DofObject::invalid_id);
2805                           }
2806                       }
2807                   }
2808             }
2809         } // end loop over variables within group
2810     } // end loop over variable groups
2811 }
2812 
2813 #endif // LIBMESH_ENABLE_AMR
2814 
2815 
2816 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2817 
find_connected_dofs(std::vector<dof_id_type> & elem_dofs)2818 void DofMap::find_connected_dofs (std::vector<dof_id_type> & elem_dofs) const
2819 {
2820   typedef std::set<dof_id_type> RCSet;
2821 
2822   // First insert the DOFS we already depend on into the set.
2823   RCSet dof_set (elem_dofs.begin(), elem_dofs.end());
2824 
2825   bool done = true;
2826 
2827   // Next insert any dofs those might be constrained in terms
2828   // of.  Note that in this case we may not be done:  Those may
2829   // in turn depend on others.  So, we need to repeat this process
2830   // in that case until the system depends only on unconstrained
2831   // degrees of freedom.
2832   for (const auto & dof : elem_dofs)
2833     if (this->is_constrained_dof(dof))
2834       {
2835         // If the DOF is constrained
2836         DofConstraints::const_iterator
2837           pos = _dof_constraints.find(dof);
2838 
2839         libmesh_assert (pos != _dof_constraints.end());
2840 
2841         const DofConstraintRow & constraint_row = pos->second;
2842 
2843         // adaptive p refinement currently gives us lots of empty constraint
2844         // rows - we should optimize those DoFs away in the future.  [RHS]
2845         //libmesh_assert (!constraint_row.empty());
2846 
2847         // Add the DOFs this dof is constrained in terms of.
2848         // note that these dofs might also be constrained, so
2849         // we will need to call this function recursively.
2850         for (const auto & pr : constraint_row)
2851           if (!dof_set.count (pr.first))
2852             {
2853               dof_set.insert (pr.first);
2854               done = false;
2855             }
2856       }
2857 
2858 
2859   // If not done then we need to do more work
2860   // (obviously :-) )!
2861   if (!done)
2862     {
2863       // Fill the vector with the contents of the set
2864       elem_dofs.clear();
2865       elem_dofs.insert (elem_dofs.end(),
2866                         dof_set.begin(), dof_set.end());
2867 
2868 
2869       // May need to do this recursively.  It is possible
2870       // that we just replaced a constrained DOF with another
2871       // constrained DOF.
2872       this->find_connected_dofs (elem_dofs);
2873 
2874     } // end if (!done)
2875 }
2876 
2877 #endif // LIBMESH_ENABLE_CONSTRAINTS
2878 
2879 
2880 
print_info(std::ostream & os)2881 void DofMap::print_info(std::ostream & os) const
2882 {
2883   os << this->get_info();
2884 }
2885 
2886 
2887 
get_info()2888 std::string DofMap::get_info() const
2889 {
2890   std::ostringstream os;
2891 
2892   // If we didn't calculate the exact sparsity pattern, the threaded
2893   // sparsity pattern assembly may have just given us an upper bound
2894   // on sparsity.
2895   const char * may_equal = " <= ";
2896 
2897   // If we calculated the exact sparsity pattern, then we can report
2898   // exact bandwidth figures:
2899   for (const auto & mat : _matrices)
2900     if (mat->need_full_sparsity_pattern())
2901       may_equal = " = ";
2902 
2903   dof_id_type max_n_nz = 0, max_n_oz = 0;
2904   long double avg_n_nz = 0, avg_n_oz = 0;
2905 
2906   if (_n_nz)
2907     {
2908       for (const auto & val : *_n_nz)
2909         {
2910           max_n_nz = std::max(max_n_nz, val);
2911           avg_n_nz += val;
2912         }
2913 
2914       std::size_t n_nz_size = _n_nz->size();
2915 
2916       this->comm().max(max_n_nz);
2917       this->comm().sum(avg_n_nz);
2918       this->comm().sum(n_nz_size);
2919 
2920       avg_n_nz /= std::max(n_nz_size,std::size_t(1));
2921 
2922       libmesh_assert(_n_oz);
2923 
2924       for (const auto & val : *_n_oz)
2925         {
2926           max_n_oz = std::max(max_n_oz, val);
2927           avg_n_oz += val;
2928         }
2929 
2930       std::size_t n_oz_size = _n_oz->size();
2931 
2932       this->comm().max(max_n_oz);
2933       this->comm().sum(avg_n_oz);
2934       this->comm().sum(n_oz_size);
2935 
2936       avg_n_oz /= std::max(n_oz_size,std::size_t(1));
2937     }
2938 
2939   os << "    DofMap Sparsity\n      Average  On-Processor Bandwidth"
2940      << may_equal << avg_n_nz << '\n';
2941 
2942   os << "      Average Off-Processor Bandwidth"
2943      << may_equal << avg_n_oz << '\n';
2944 
2945   os << "      Maximum  On-Processor Bandwidth"
2946      << may_equal << max_n_nz << '\n';
2947 
2948   os << "      Maximum Off-Processor Bandwidth"
2949      << may_equal << max_n_oz << std::endl;
2950 
2951 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2952 
2953   std::size_t n_constraints = 0, max_constraint_length = 0,
2954     n_rhss = 0;
2955   long double avg_constraint_length = 0.;
2956 
2957   for (const auto & pr : _dof_constraints)
2958     {
2959       // Only count local constraints, then sum later
2960       const dof_id_type constrained_dof = pr.first;
2961       if (!this->local_index(constrained_dof))
2962         continue;
2963 
2964       const DofConstraintRow & row = pr.second;
2965       std::size_t rowsize = row.size();
2966 
2967       max_constraint_length = std::max(max_constraint_length,
2968                                        rowsize);
2969       avg_constraint_length += rowsize;
2970       n_constraints++;
2971 
2972       if (_primal_constraint_values.count(constrained_dof))
2973         n_rhss++;
2974     }
2975 
2976   this->comm().sum(n_constraints);
2977   this->comm().sum(n_rhss);
2978   this->comm().sum(avg_constraint_length);
2979   this->comm().max(max_constraint_length);
2980 
2981   os << "    DofMap Constraints\n      Number of DoF Constraints = "
2982      << n_constraints;
2983   if (n_rhss)
2984     os << '\n'
2985        << "      Number of Heterogenous Constraints= " << n_rhss;
2986   if (n_constraints)
2987     {
2988       avg_constraint_length /= n_constraints;
2989 
2990       os << '\n'
2991          << "      Average DoF Constraint Length= " << avg_constraint_length;
2992     }
2993 
2994 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2995   std::size_t n_node_constraints = 0, max_node_constraint_length = 0,
2996     n_node_rhss = 0;
2997   long double avg_node_constraint_length = 0.;
2998 
2999   for (const auto & pr : _node_constraints)
3000     {
3001       // Only count local constraints, then sum later
3002       const Node * node = pr.first;
3003       if (node->processor_id() != this->processor_id())
3004         continue;
3005 
3006       const NodeConstraintRow & row = pr.second.first;
3007       std::size_t rowsize = row.size();
3008 
3009       max_node_constraint_length = std::max(max_node_constraint_length,
3010                                             rowsize);
3011       avg_node_constraint_length += rowsize;
3012       n_node_constraints++;
3013 
3014       if (pr.second.second != Point(0))
3015         n_node_rhss++;
3016     }
3017 
3018   this->comm().sum(n_node_constraints);
3019   this->comm().sum(n_node_rhss);
3020   this->comm().sum(avg_node_constraint_length);
3021   this->comm().max(max_node_constraint_length);
3022 
3023   os << "\n      Number of Node Constraints = " << n_node_constraints;
3024   if (n_node_rhss)
3025     os << '\n'
3026        << "      Number of Heterogenous Node Constraints= " << n_node_rhss;
3027   if (n_node_constraints)
3028     {
3029       avg_node_constraint_length /= n_node_constraints;
3030       os << "\n      Maximum Node Constraint Length= " << max_node_constraint_length
3031          << '\n'
3032          << "      Average Node Constraint Length= " << avg_node_constraint_length;
3033     }
3034 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
3035 
3036   os << std::endl;
3037 
3038 #endif // LIBMESH_ENABLE_CONSTRAINTS
3039 
3040   return os.str();
3041 }
3042 
3043 
3044 template bool DofMap::is_evaluable<Elem>(const Elem &, unsigned int) const;
3045 template bool DofMap::is_evaluable<Node>(const Node &, unsigned int) const;
3046 
3047 } // namespace libMesh
3048