libMesh
dof_map.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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
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>
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,
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
112  {
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 
122  (sp->sparsity_pattern, sp->n_nz,
123  sp->n_oz, _extra_sparsity_context);
124  }
125 
128  (sp->sparsity_pattern, sp->n_nz, sp->n_oz);
129 
130  return std::unique_ptr<SparsityPattern::Build>(sp.release());
131 }
132 
133 
134 
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 
201 }
202 
203 
204 
205 // Destructor
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.
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 
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 
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 
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 
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 : IntRange<unsigned int>(0, 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 
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 &&
295  {
296  // We'd better have already computed the full sparsity pattern
297  // if we need it here
299  libmesh_assert(_sp.get());
300 
301  matrix.update_sparsity_pattern (_sp->sparsity_pattern);
302  }
303 
304  if (matrix.need_full_sparsity_pattern())
306 }
307 
308 
309 
311 {
312  return (std::find(_matrices.begin(), _matrices.end(),
313  &matrix) != _matrices.end());
314 }
315 
316 
317 
319 {
320  return mesh.node_ptr(i);
321 }
322 
323 
324 
326 {
327  return mesh.elem_ptr(i);
328 }
329 
330 
331 
332 template <typename iterator_type>
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;
377  requested_ids[obj->processor_id()].push_back(obj->id());
378  }
379 #ifdef DEBUG
380  for (auto p : IntRange<processor_id_type>(0, 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]
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 
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  const ElemType type = elem->type();
614  const unsigned int dim = elem->dim();
615 
616  FEType fe_type = base_fe_type;
617 
618 #ifdef LIBMESH_ENABLE_AMR
619  // Make sure we haven't done more p refinement than we can
620  // handle
621  if (elem->p_level() + base_fe_type.order >
622  FEInterface::max_order(base_fe_type, type))
623  {
624  libmesh_assert_less_msg(base_fe_type.order.get_order(),
625  FEInterface::max_order(base_fe_type,type),
626  "ERROR: Finite element "
627  << Utility::enum_to_string(base_fe_type.family)
628  << " on geometric element "
629  << Utility::enum_to_string(type)
630  << "\nonly supports FEInterface::max_order = "
631  << FEInterface::max_order(base_fe_type,type)
632  << ", not fe_type.order = "
633  << base_fe_type.order);
634 
635 # ifdef DEBUG
636  libMesh::err << "WARNING: Finite element "
637  << Utility::enum_to_string(base_fe_type.family)
638  << " on geometric element "
639  << Utility::enum_to_string(type) << std::endl
640  << "could not be p refined past FEInterface::max_order = "
641  << FEInterface::max_order(base_fe_type,type)
642  << std::endl;
643 # endif
644  elem->set_p_level(FEInterface::max_order(base_fe_type,type)
645  - base_fe_type.order);
646  }
647 #endif
648 
649  fe_type.order = static_cast<Order>(fe_type.order +
650  elem->p_level());
651 
652  // Allocate the vertex DOFs
653  for (auto n : elem->node_index_range())
654  {
655  Node & node = elem->node_ref(n);
656 
657  if (elem->is_vertex(n))
658  {
659  const unsigned int old_node_dofs =
660  node.n_comp_group(sys_num, vg);
661 
662  const unsigned int vertex_dofs =
663  std::max(FEInterface::n_dofs_at_node(dim, fe_type,
664  type, n),
665  old_node_dofs);
666 
667  // Some discontinuous FEs have no vertex dofs
668  if (vertex_dofs > old_node_dofs)
669  {
670  node.set_n_comp_group(sys_num, vg,
671  vertex_dofs);
672 
673  // Abusing dof_number to set a "this is a
674  // vertex" flag
675  node.set_vg_dof_base(sys_num, vg,
676  vertex_dofs);
677 
678  // libMesh::out << "sys_num,vg,old_node_dofs,vertex_dofs="
679  // << sys_num << ","
680  // << vg << ","
681  // << old_node_dofs << ","
682  // << vertex_dofs << '\n',
683  // node.debug_buffer();
684 
685  // libmesh_assert_equal_to (vertex_dofs, node.n_comp(sys_num, vg));
686  // libmesh_assert_equal_to (vertex_dofs, node.vg_dof_base(sys_num, vg));
687  }
688  }
689  }
690  } // done counting vertex dofs
691 
692  // count edge & face dofs next
693  for (auto & elem : mesh.active_element_ptr_range())
694  {
695  libmesh_assert(elem);
696 
697  // Skip the numbering if this variable is
698  // not active on this element's subdomain
699  if (!vg_description.active_on_subdomain(elem->subdomain_id()))
700  continue;
701 
702  const ElemType type = elem->type();
703  const unsigned int dim = elem->dim();
704 
705  FEType fe_type = base_fe_type;
706  fe_type.order = static_cast<Order>(fe_type.order +
707  elem->p_level());
708 
709  // Allocate the edge and face DOFs
710  for (auto n : elem->node_index_range())
711  {
712  Node & node = elem->node_ref(n);
713 
714  const unsigned int old_node_dofs =
715  node.n_comp_group(sys_num, vg);
716 
717  const unsigned int vertex_dofs = old_node_dofs?
718  cast_int<unsigned int>(node.vg_dof_base (sys_num,vg)):0;
719 
720  const unsigned int new_node_dofs =
721  FEInterface::n_dofs_at_node(dim, fe_type, type, n);
722 
723  // We've already allocated vertex DOFs
724  if (elem->is_vertex(n))
725  {
726  libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
727  // //if (vertex_dofs < new_node_dofs)
728  // libMesh::out << "sys_num,vg,old_node_dofs,vertex_dofs,new_node_dofs="
729  // << sys_num << ","
730  // << vg << ","
731  // << old_node_dofs << ","
732  // << vertex_dofs << ","
733  // << new_node_dofs << '\n',
734  // node.debug_buffer();
735 
736  libmesh_assert_greater_equal (vertex_dofs, new_node_dofs);
737  }
738  // We need to allocate the rest
739  else
740  {
741  // If this has no dofs yet, it needs no vertex
742  // dofs, so we just give it edge or face dofs
743  if (!old_node_dofs)
744  {
745  node.set_n_comp_group(sys_num, vg,
746  new_node_dofs);
747  // Abusing dof_number to set a "this has no
748  // vertex dofs" flag
749  if (new_node_dofs)
750  node.set_vg_dof_base(sys_num, vg, 0);
751  }
752 
753  // If this has dofs, but has no vertex dofs,
754  // it may still need more edge or face dofs if
755  // we're p-refined.
756  else if (vertex_dofs == 0)
757  {
758  if (new_node_dofs > old_node_dofs)
759  {
760  node.set_n_comp_group(sys_num, vg,
761  new_node_dofs);
762 
763  node.set_vg_dof_base(sys_num, vg,
764  vertex_dofs);
765  }
766  }
767  // If this is another element's vertex,
768  // add more (non-overlapping) edge/face dofs if
769  // necessary
770  else if (extra_hanging_dofs)
771  {
772  if (new_node_dofs > old_node_dofs - vertex_dofs)
773  {
774  node.set_n_comp_group(sys_num, vg,
775  vertex_dofs + new_node_dofs);
776 
777  node.set_vg_dof_base(sys_num, vg,
778  vertex_dofs);
779  }
780  }
781  // If this is another element's vertex, add any
782  // (overlapping) edge/face dofs if necessary
783  else
784  {
785  libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
786  if (new_node_dofs > old_node_dofs)
787  {
788  node.set_n_comp_group(sys_num, vg,
789  new_node_dofs);
790 
791  node.set_vg_dof_base (sys_num, vg,
792  vertex_dofs);
793  }
794  }
795  }
796  }
797  // Allocate the element DOFs
798  const unsigned int dofs_per_elem =
800  type);
801 
802  elem->set_n_comp_group(sys_num, vg, dofs_per_elem);
803 
804  }
805  } // end loop over variable groups
806 
807  // Calling DofMap::reinit() by itself makes little sense,
808  // so we won't bother with nonlocal DofObjects.
809  // Those will be fixed by distribute_dofs
810 
811  //------------------------------------------------------------
812  // Finally, clear all the current DOF indices
813  // (distribute_dofs expects them cleared!)
814  this->invalidate_dofs(mesh);
815 }
816 
817 
818 
820 {
821  const unsigned int sys_num = this->sys_number();
822 
823  // All the nodes
824  for (auto & node : mesh.node_ptr_range())
825  node->invalidate_dofs(sys_num);
826 
827  // All the active elements.
828  for (auto & elem : mesh.active_element_ptr_range())
829  elem->invalidate_dofs(sys_num);
830 }
831 
832 
833 
835 {
836  // we don't want to clear
837  // the coupling matrix!
838  // It should not change...
839  //_dof_coupling->clear();
840  //
841  // But it would be inconsistent to leave our coupling settings
842  // through a clear()...
843  _dof_coupling = nullptr;
844 
845  // Reset ghosting functor statuses
846  {
847  for (const auto & gf : _coupling_functors)
848  {
849  libmesh_assert(gf);
851  }
852  this->_coupling_functors.clear();
853 
854  // Go back to default coupling
855 
856  _default_coupling->set_dof_coupling(this->_dof_coupling);
857  _default_coupling->set_n_levels(this->use_coupled_neighbor_dofs(this->_mesh));
858 
860  }
861 
862 
863  {
864  for (const auto & gf : _algebraic_ghosting_functors)
865  {
866  libmesh_assert(gf);
868  }
869  this->_algebraic_ghosting_functors.clear();
870 
871  // Go back to default send_list generation
872 
873  // _default_evaluating->set_dof_coupling(this->_dof_coupling);
874  _default_evaluating->set_n_levels(1);
876  }
877 
878  this->_shared_functors.clear();
879 
880  _variables.clear();
881  _variable_groups.clear();
882  _variable_group_numbers.clear();
883  _first_df.clear();
884  _end_df.clear();
885  _first_scalar_df.clear();
886  this->clear_send_list();
887  this->clear_sparsity();
889 
890 #ifdef LIBMESH_ENABLE_AMR
891 
892  _dof_constraints.clear();
893  _stashed_dof_constraints.clear();
896  _n_old_dfs = 0;
897  _first_old_df.clear();
898  _end_old_df.clear();
899  _first_old_scalar_df.clear();
900 
901 #endif
902 
903  _matrices.clear();
904 
905  _n_dfs = 0;
906 }
907 
908 
909 
911 {
912  // This function must be run on all processors at once
913  parallel_object_only();
914 
915  // Log how long it takes to distribute the degrees of freedom
916  LOG_SCOPE("distribute_dofs()", "DofMap");
917 
918  libmesh_assert (mesh.is_prepared());
919 
920  const processor_id_type proc_id = this->processor_id();
921  const processor_id_type n_proc = this->n_processors();
922 
923  // libmesh_assert_greater (this->n_variables(), 0);
924  libmesh_assert_less (proc_id, n_proc);
925 
926  // re-init in case the mesh has changed
927  this->reinit(mesh);
928 
929  // By default distribute variables in a
930  // var-major fashion, but allow run-time
931  // specification
932  bool node_major_dofs = libMesh::on_command_line ("--node-major-dofs");
933 
934  // The DOF counter, will be incremented as we encounter
935  // new degrees of freedom
936  dof_id_type next_free_dof = 0;
937 
938  // Clear the send list before we rebuild it
939  this->clear_send_list();
940 
941  // Set temporary DOF indices on this processor
942  if (node_major_dofs)
943  this->distribute_local_dofs_node_major (next_free_dof, mesh);
944  else
945  this->distribute_local_dofs_var_major (next_free_dof, mesh);
946 
947  // Get DOF counts on all processors
948  std::vector<dof_id_type> dofs_on_proc(n_proc, 0);
949  this->comm().allgather(next_free_dof, dofs_on_proc);
950 
951  // Resize and fill the _first_df and _end_df arrays
952 #ifdef LIBMESH_ENABLE_AMR
955 #endif
956 
957  _first_df.resize(n_proc);
958  _end_df.resize (n_proc);
959 
960  // Get DOF offsets
961  _first_df[0] = 0;
962  for (processor_id_type i=1; i < n_proc; ++i)
963  _first_df[i] = _end_df[i-1] = _first_df[i-1] + dofs_on_proc[i-1];
964  _end_df[n_proc-1] = _first_df[n_proc-1] + dofs_on_proc[n_proc-1];
965 
966  // Clear all the current DOF indices
967  // (distribute_dofs expects them cleared!)
968  this->invalidate_dofs(mesh);
969 
970  next_free_dof = _first_df[proc_id];
971 
972  // Set permanent DOF indices on this processor
973  if (node_major_dofs)
974  this->distribute_local_dofs_node_major (next_free_dof, mesh);
975  else
976  this->distribute_local_dofs_var_major (next_free_dof, mesh);
977 
978  libmesh_assert_equal_to (next_free_dof, _end_df[proc_id]);
979 
980  //------------------------------------------------------------
981  // At this point, all n_comp and dof_number values on local
982  // DofObjects should be correct, but a DistributedMesh might have
983  // incorrect values on non-local DofObjects. Let's request the
984  // correct values from each other processor.
985 
986  if (this->n_processors() > 1)
987  {
989  mesh.nodes_end(),
991 
993  mesh.elements_end(),
995  }
996 
997 #ifdef DEBUG
998  {
999  const unsigned int
1000  sys_num = this->sys_number();
1001 
1002  // Processors should all agree on DoF ids for the newly numbered
1003  // system.
1005 
1006  // DoF processor ids should match DofObject processor ids
1007  for (auto & node : mesh.node_ptr_range())
1008  {
1009  DofObject const * const dofobj = node;
1010  const processor_id_type obj_proc_id = dofobj->processor_id();
1011 
1012  for (auto v : IntRange<unsigned int>(0, dofobj->n_vars(sys_num)))
1013  for (auto c : IntRange<unsigned int>(0, dofobj->n_comp(sys_num,v)))
1014  {
1015  const dof_id_type dofid = dofobj->dof_number(sys_num,v,c);
1016  libmesh_assert_greater_equal (dofid, this->first_dof(obj_proc_id));
1017  libmesh_assert_less (dofid, this->end_dof(obj_proc_id));
1018  }
1019  }
1020 
1021  for (auto & elem : mesh.element_ptr_range())
1022  {
1023  DofObject const * const dofobj = elem;
1024  const processor_id_type obj_proc_id = dofobj->processor_id();
1025 
1026  for (auto v : IntRange<unsigned int>(0, dofobj->n_vars(sys_num)))
1027  for (auto c : IntRange<unsigned int>(0, dofobj->n_comp(sys_num,v)))
1028  {
1029  const dof_id_type dofid = dofobj->dof_number(sys_num,v,c);
1030  libmesh_assert_greater_equal (dofid, this->first_dof(obj_proc_id));
1031  libmesh_assert_less (dofid, this->end_dof(obj_proc_id));
1032  }
1033  }
1034  }
1035 #endif
1036 
1037  // Set the total number of degrees of freedom, then start finding
1038  // SCALAR degrees of freedom
1039 #ifdef LIBMESH_ENABLE_AMR
1040  _n_old_dfs = _n_dfs;
1042 #endif
1043  _n_dfs = _end_df[n_proc-1];
1044  _first_scalar_df.clear();
1046  dof_id_type current_SCALAR_dof_index = n_dofs() - n_SCALAR_dofs();
1047 
1048  // Calculate and cache the initial DoF indices for SCALAR variables.
1049  // This is an O(N_vars) calculation so we want to do it once per
1050  // renumbering rather than once per SCALAR_dof_indices() call
1051 
1052  for (auto v : IntRange<unsigned int>(0, this->n_variables()))
1053  if (this->variable(v).type().family == SCALAR)
1054  {
1055  _first_scalar_df[v] = current_SCALAR_dof_index;
1056  current_SCALAR_dof_index += this->variable(v).type().order.get_order();
1057  }
1058 
1059  // Allow our GhostingFunctor objects to reinit if necessary
1060  for (const auto & gf : _algebraic_ghosting_functors)
1061  {
1062  libmesh_assert(gf);
1063  gf->dofmap_reinit();
1064  }
1065 
1066  for (const auto & gf : _coupling_functors)
1067  {
1068  libmesh_assert(gf);
1069  gf->dofmap_reinit();
1070  }
1071 
1072  // Note that in the add_neighbors_to_send_list nodes on processor
1073  // boundaries that are shared by multiple elements are added for
1074  // each element.
1075  this->add_neighbors_to_send_list(mesh);
1076 
1077  // Here we used to clean up that data structure; now System and
1078  // EquationSystems call that for us, after we've added constraint
1079  // dependencies to the send_list too.
1080  // this->sort_send_list ();
1081 }
1082 
1083 
1084 void DofMap::local_variable_indices(std::vector<dof_id_type> & idx,
1085  const MeshBase & mesh,
1086  unsigned int var_num) const
1087 {
1088  // Count dofs in the *exact* order that distribute_dofs numbered
1089  // them, so that we can assume ascending indices and use push_back
1090  // instead of find+insert.
1091 
1092  const unsigned int sys_num = this->sys_number();
1093 
1094  // If this isn't a SCALAR variable, we need to find all its field
1095  // dofs on the mesh
1096  if (this->variable_type(var_num).family != SCALAR)
1097  {
1098  const Variable & var(this->variable(var_num));
1099 
1100  for (auto & elem : mesh.active_local_element_ptr_range())
1101  {
1102  if (!var.active_on_subdomain(elem->subdomain_id()))
1103  continue;
1104 
1105  // Only count dofs connected to active
1106  // elements on this processor.
1107  const unsigned int n_nodes = elem->n_nodes();
1108 
1109  // First get any new nodal DOFS
1110  for (unsigned int n=0; n<n_nodes; n++)
1111  {
1112  Node & node = elem->node_ref(n);
1113 
1114  if (node.processor_id() != this->processor_id())
1115  continue;
1116 
1117  const unsigned int n_comp = node.n_comp(sys_num, var_num);
1118  for(unsigned int i=0; i<n_comp; i++)
1119  {
1120  const dof_id_type index = node.dof_number(sys_num,var_num,i);
1121  libmesh_assert (this->local_index(index));
1122 
1123  if (idx.empty() || index > idx.back())
1124  idx.push_back(index);
1125  }
1126  }
1127 
1128  // Next get any new element DOFS
1129  const unsigned int n_comp = elem->n_comp(sys_num, var_num);
1130  for (unsigned int i=0; i<n_comp; i++)
1131  {
1132  const dof_id_type index = elem->dof_number(sys_num,var_num,i);
1133  if (idx.empty() || index > idx.back())
1134  idx.push_back(index);
1135  }
1136  } // done looping over elements
1137 
1138 
1139  // we may have missed assigning DOFs to nodes that we own
1140  // but to which we have no connected elements matching our
1141  // variable restriction criterion. this will happen, for example,
1142  // if variable V is restricted to subdomain S. We may not own
1143  // any elements which live in S, but we may own nodes which are
1144  // *connected* to elements which do. in this scenario these nodes
1145  // will presently have unnumbered DOFs. we need to take care of
1146  // them here since we own them and no other processor will touch them.
1147  for (const auto & node : mesh.local_node_ptr_range())
1148  {
1149  libmesh_assert(node);
1150 
1151  const unsigned int n_comp = node->n_comp(sys_num, var_num);
1152  for (unsigned int i=0; i<n_comp; i++)
1153  {
1154  const dof_id_type index = node->dof_number(sys_num,var_num,i);
1155  if (idx.empty() || index > idx.back())
1156  idx.push_back(index);
1157  }
1158  }
1159  }
1160  // Otherwise, count up the SCALAR dofs, if we're on the processor
1161  // that holds this SCALAR variable
1162  else if (this->processor_id() == (this->n_processors()-1))
1163  {
1164  std::vector<dof_id_type> di_scalar;
1165  this->SCALAR_dof_indices(di_scalar,var_num);
1166  idx.insert( idx.end(), di_scalar.begin(), di_scalar.end());
1167  }
1168 }
1169 
1170 
1172  MeshBase & mesh)
1173 {
1174  const unsigned int sys_num = this->sys_number();
1175  const unsigned int n_var_groups = this->n_variable_groups();
1176 
1177  // Our numbering here must be kept consistent with the numbering
1178  // scheme assumed by DofMap::local_variable_indices!
1179 
1180  //-------------------------------------------------------------------------
1181  // First count and assign temporary numbers to local dofs
1182  for (auto & elem : mesh.active_local_element_ptr_range())
1183  {
1184  // Only number dofs connected to active
1185  // elements on this processor.
1186  const unsigned int n_nodes = elem->n_nodes();
1187 
1188  // First number the nodal DOFS
1189  for (unsigned int n=0; n<n_nodes; n++)
1190  {
1191  Node & node = elem->node_ref(n);
1192 
1193  for (unsigned vg=0; vg<n_var_groups; vg++)
1194  {
1195  const VariableGroup & vg_description(this->variable_group(vg));
1196 
1197  if ((vg_description.type().family != SCALAR) &&
1198  (vg_description.active_on_subdomain(elem->subdomain_id())))
1199  {
1200  // assign dof numbers (all at once) if this is
1201  // our node and if they aren't already there
1202  if ((node.n_comp_group(sys_num,vg) > 0) &&
1203  (node.processor_id() == this->processor_id()) &&
1204  (node.vg_dof_base(sys_num,vg) ==
1206  {
1207  node.set_vg_dof_base(sys_num, vg,
1208  next_free_dof);
1209  next_free_dof += (vg_description.n_variables()*
1210  node.n_comp_group(sys_num,vg));
1211  //node.debug_buffer();
1212  }
1213  }
1214  }
1215  }
1216 
1217  // Now number the element DOFS
1218  for (unsigned vg=0; vg<n_var_groups; vg++)
1219  {
1220  const VariableGroup & vg_description(this->variable_group(vg));
1221 
1222  if ((vg_description.type().family != SCALAR) &&
1223  (vg_description.active_on_subdomain(elem->subdomain_id())))
1224  if (elem->n_comp_group(sys_num,vg) > 0)
1225  {
1226  libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
1228 
1229  elem->set_vg_dof_base(sys_num,
1230  vg,
1231  next_free_dof);
1232 
1233  next_free_dof += (vg_description.n_variables()*
1234  elem->n_comp(sys_num,vg));
1235  }
1236  }
1237  } // done looping over elements
1238 
1239 
1240  // we may have missed assigning DOFs to nodes that we own
1241  // but to which we have no connected elements matching our
1242  // variable restriction criterion. this will happen, for example,
1243  // if variable V is restricted to subdomain S. We may not own
1244  // any elements which live in S, but we may own nodes which are
1245  // *connected* to elements which do. in this scenario these nodes
1246  // will presently have unnumbered DOFs. we need to take care of
1247  // them here since we own them and no other processor will touch them.
1248  for (auto & node : mesh.local_node_ptr_range())
1249  for (unsigned vg=0; vg<n_var_groups; vg++)
1250  {
1251  const VariableGroup & vg_description(this->variable_group(vg));
1252 
1253  if (node->n_comp_group(sys_num,vg))
1254  if (node->vg_dof_base(sys_num,vg) == DofObject::invalid_id)
1255  {
1256  node->set_vg_dof_base (sys_num,
1257  vg,
1258  next_free_dof);
1259 
1260  next_free_dof += (vg_description.n_variables()*
1261  node->n_comp(sys_num,vg));
1262  }
1263  }
1264 
1265  // Finally, count up the SCALAR dofs
1266  this->_n_SCALAR_dofs = 0;
1267  for (unsigned vg=0; vg<n_var_groups; vg++)
1268  {
1269  const VariableGroup & vg_description(this->variable_group(vg));
1270 
1271  if (vg_description.type().family == SCALAR)
1272  {
1273  this->_n_SCALAR_dofs += (vg_description.n_variables()*
1274  vg_description.type().order.get_order());
1275  continue;
1276  }
1277  }
1278 
1279  // Only increment next_free_dof if we're on the processor
1280  // that holds this SCALAR variable
1281  if (this->processor_id() == (this->n_processors()-1))
1282  next_free_dof += _n_SCALAR_dofs;
1283 
1284 #ifdef DEBUG
1285  {
1286  // libMesh::out << "next_free_dof=" << next_free_dof << std::endl
1287  // << "_n_SCALAR_dofs=" << _n_SCALAR_dofs << std::endl;
1288 
1289  // Make sure we didn't miss any nodes
1290  MeshTools::libmesh_assert_valid_procids<Node>(mesh);
1291 
1292  for (auto & node : mesh.local_node_ptr_range())
1293  {
1294  unsigned int n_var_g = node->n_var_groups(this->sys_number());
1295  for (unsigned int vg=0; vg != n_var_g; ++vg)
1296  {
1297  unsigned int n_comp_g =
1298  node->n_comp_group(this->sys_number(), vg);
1299  dof_id_type my_first_dof = n_comp_g ?
1300  node->vg_dof_base(this->sys_number(), vg) : 0;
1301  libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
1302  }
1303  }
1304  }
1305 #endif // DEBUG
1306 }
1307 
1308 
1309 
1311  MeshBase & mesh)
1312 {
1313  const unsigned int sys_num = this->sys_number();
1314  const unsigned int n_var_groups = this->n_variable_groups();
1315 
1316  // Our numbering here must be kept consistent with the numbering
1317  // scheme assumed by DofMap::local_variable_indices!
1318 
1319  //-------------------------------------------------------------------------
1320  // First count and assign temporary numbers to local dofs
1321  for (unsigned vg=0; vg<n_var_groups; vg++)
1322  {
1323  const VariableGroup & vg_description(this->variable_group(vg));
1324 
1325  const unsigned int n_vars_in_group = vg_description.n_variables();
1326 
1327  // Skip the SCALAR dofs
1328  if (vg_description.type().family == SCALAR)
1329  continue;
1330 
1331  for (auto & elem : mesh.active_local_element_ptr_range())
1332  {
1333  // Only number dofs connected to active elements on this
1334  // processor and only variables which are active on on this
1335  // element's subdomain.
1336  if (!vg_description.active_on_subdomain(elem->subdomain_id()))
1337  continue;
1338 
1339  const unsigned int n_nodes = elem->n_nodes();
1340 
1341  // First number the nodal DOFS
1342  for (unsigned int n=0; n<n_nodes; n++)
1343  {
1344  Node & node = elem->node_ref(n);
1345 
1346  // assign dof numbers (all at once) if this is
1347  // our node and if they aren't already there
1348  if ((node.n_comp_group(sys_num,vg) > 0) &&
1349  (node.processor_id() == this->processor_id()) &&
1350  (node.vg_dof_base(sys_num,vg) ==
1352  {
1353  node.set_vg_dof_base(sys_num, vg, next_free_dof);
1354 
1355  next_free_dof += (n_vars_in_group*
1356  node.n_comp_group(sys_num,vg));
1357  }
1358  }
1359 
1360  // Now number the element DOFS
1361  if (elem->n_comp_group(sys_num,vg) > 0)
1362  {
1363  libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
1365 
1366  elem->set_vg_dof_base(sys_num,
1367  vg,
1368  next_free_dof);
1369 
1370  next_free_dof += (n_vars_in_group*
1371  elem->n_comp_group(sys_num,vg));
1372  }
1373  } // end loop on elements
1374 
1375  // we may have missed assigning DOFs to nodes that we own
1376  // but to which we have no connected elements matching our
1377  // variable restriction criterion. this will happen, for example,
1378  // if variable V is restricted to subdomain S. We may not own
1379  // any elements which live in S, but we may own nodes which are
1380  // *connected* to elements which do. in this scenario these nodes
1381  // will presently have unnumbered DOFs. we need to take care of
1382  // them here since we own them and no other processor will touch them.
1383  for (auto & node : mesh.local_node_ptr_range())
1384  if (node->n_comp_group(sys_num,vg))
1385  if (node->vg_dof_base(sys_num,vg) == DofObject::invalid_id)
1386  {
1387  node->set_vg_dof_base (sys_num,
1388  vg,
1389  next_free_dof);
1390 
1391  next_free_dof += (n_vars_in_group*
1392  node->n_comp_group(sys_num,vg));
1393  }
1394  } // end loop on variable groups
1395 
1396  // Finally, count up the SCALAR dofs
1397  this->_n_SCALAR_dofs = 0;
1398  for (unsigned vg=0; vg<n_var_groups; vg++)
1399  {
1400  const VariableGroup & vg_description(this->variable_group(vg));
1401 
1402  if (vg_description.type().family == SCALAR)
1403  {
1404  this->_n_SCALAR_dofs += (vg_description.n_variables()*
1405  vg_description.type().order.get_order());
1406  continue;
1407  }
1408  }
1409 
1410  // Only increment next_free_dof if we're on the processor
1411  // that holds this SCALAR variable
1412  if (this->processor_id() == (this->n_processors()-1))
1413  next_free_dof += _n_SCALAR_dofs;
1414 
1415 #ifdef DEBUG
1416  {
1417  // Make sure we didn't miss any nodes
1418  MeshTools::libmesh_assert_valid_procids<Node>(mesh);
1419 
1420  for (auto & node : mesh.local_node_ptr_range())
1421  {
1422  unsigned int n_var_g = node->n_var_groups(this->sys_number());
1423  for (unsigned int vg=0; vg != n_var_g; ++vg)
1424  {
1425  unsigned int n_comp_g =
1426  node->n_comp_group(this->sys_number(), vg);
1427  dof_id_type my_first_dof = n_comp_g ?
1428  node->vg_dof_base(this->sys_number(), vg) : 0;
1429  libmesh_assert_not_equal_to (my_first_dof, DofObject::invalid_id);
1430  }
1431  }
1432  }
1433 #endif // DEBUG
1434 }
1435 
1436 
1437 
1438 void
1439 DofMap::
1441  std::set<CouplingMatrix *> & temporary_coupling_matrices,
1442  const std::set<GhostingFunctor *>::iterator & gf_begin,
1443  const std::set<GhostingFunctor *>::iterator & gf_end,
1444  const MeshBase::const_element_iterator & elems_begin,
1445  const MeshBase::const_element_iterator & elems_end,
1447 {
1448  for (const auto & gf : as_range(gf_begin, gf_end))
1449  {
1450  GhostingFunctor::map_type more_elements_to_ghost;
1451 
1452  libmesh_assert(gf);
1453  (*gf)(elems_begin, elems_end, p, more_elements_to_ghost);
1454 
1455  for (const auto & pr : more_elements_to_ghost)
1456  {
1457  GhostingFunctor::map_type::iterator existing_it =
1458  elements_to_ghost.find (pr.first);
1459  if (existing_it == elements_to_ghost.end())
1460  elements_to_ghost.insert(pr);
1461  else
1462  {
1463  if (existing_it->second)
1464  {
1465  if (pr.second)
1466  {
1467  // If this isn't already a temporary
1468  // then we need to make one so we'll
1469  // have a non-const matrix to merge
1470  if (temporary_coupling_matrices.empty() ||
1471  temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second)) == temporary_coupling_matrices.end())
1472  {
1473  CouplingMatrix * cm = new CouplingMatrix(*existing_it->second);
1474  temporary_coupling_matrices.insert(cm);
1475  existing_it->second = cm;
1476  }
1477  const_cast<CouplingMatrix &>(*existing_it->second) &= *pr.second;
1478  }
1479  else
1480  {
1481  // Any existing_it matrix merged with a full
1482  // matrix (symbolized as nullptr) gives another
1483  // full matrix (symbolizable as nullptr).
1484 
1485  // So if existing_it->second is a temporary then
1486  // we don't need it anymore; we might as well
1487  // remove it to keep the set of temporaries
1488  // small.
1489  std::set<CouplingMatrix *>::iterator temp_it =
1490  temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second));
1491  if (temp_it != temporary_coupling_matrices.end())
1492  temporary_coupling_matrices.erase(temp_it);
1493 
1494  existing_it->second = nullptr;
1495  }
1496  }
1497  // else we have a nullptr already, then we have a full
1498  // coupling matrix, already, and merging with anything
1499  // else won't change that, so we're done.
1500  }
1501  }
1502  }
1503 }
1504 
1505 
1506 
1508 {
1509  LOG_SCOPE("add_neighbors_to_send_list()", "DofMap");
1510 
1511  // Return immediately if there's no ghost data
1512  if (this->n_processors() == 1)
1513  return;
1514 
1515  const unsigned int n_var = this->n_variables();
1516 
1517  MeshBase::const_element_iterator local_elem_it
1518  = mesh.active_local_elements_begin();
1519  const MeshBase::const_element_iterator local_elem_end
1520  = mesh.active_local_elements_end();
1521 
1522  GhostingFunctor::map_type elements_to_send;
1523 
1524  // Man, I wish we had guaranteed unique_ptr availability...
1525  std::set<CouplingMatrix *> temporary_coupling_matrices;
1526 
1527  // We need to add dofs to the send list if they've been directly
1528  // requested by an algebraic ghosting functor or they've been
1529  // indirectly requested by a coupling functor.
1530  this->merge_ghost_functor_outputs(elements_to_send,
1531  temporary_coupling_matrices,
1534  local_elem_it, local_elem_end, mesh.processor_id());
1535 
1536  this->merge_ghost_functor_outputs(elements_to_send,
1537  temporary_coupling_matrices,
1538  this->coupling_functors_begin(),
1539  this->coupling_functors_end(),
1540  local_elem_it, local_elem_end, mesh.processor_id());
1541 
1542  // Making a list of non-zero coupling matrix columns is an
1543  // O(N_var^2) operation. We cache it so we only have to do it once
1544  // per CouplingMatrix and not once per element.
1545  std::map<const CouplingMatrix *, std::vector<unsigned int>>
1546  column_variable_lists;
1547 
1548  for (auto & pr : elements_to_send)
1549  {
1550  const Elem * const partner = pr.first;
1551 
1552  // We asked ghosting functors not to give us local elements
1553  libmesh_assert_not_equal_to
1554  (partner->processor_id(), this->processor_id());
1555 
1556  const CouplingMatrix * ghost_coupling = pr.second;
1557 
1558  // Loop over any present coupling matrix column variables if we
1559  // have a coupling matrix, or just add all variables to
1560  // send_list if not.
1561  if (ghost_coupling)
1562  {
1563  libmesh_assert_equal_to (ghost_coupling->size(), n_var);
1564 
1565  // Try to find a cached list of column variables.
1566  std::map<const CouplingMatrix *, std::vector<unsigned int>>::const_iterator
1567  column_variable_list = column_variable_lists.find(ghost_coupling);
1568 
1569  // If we didn't find it, then we need to create it.
1570  if (column_variable_list == column_variable_lists.end())
1571  {
1572  std::pair<std::map<const CouplingMatrix *, std::vector<unsigned int>>::iterator, bool>
1573  inserted_variable_list_pair = column_variable_lists.insert(std::make_pair(ghost_coupling,
1574  std::vector<unsigned int>()));
1575  column_variable_list = inserted_variable_list_pair.first;
1576 
1577  std::vector<unsigned int> & new_variable_list =
1578  inserted_variable_list_pair.first->second;
1579 
1580  std::vector<unsigned char> has_variable(n_var, false);
1581 
1582  for (unsigned int vi = 0; vi != n_var; ++vi)
1583  {
1584  ConstCouplingRow ccr(vi, *ghost_coupling);
1585 
1586  for (const auto & vj : ccr)
1587  has_variable[vj] = true;
1588  }
1589  for (unsigned int vj = 0; vj != n_var; ++vj)
1590  {
1591  if (has_variable[vj])
1592  new_variable_list.push_back(vj);
1593  }
1594  }
1595 
1596  const std::vector<unsigned int> & variable_list =
1597  column_variable_list->second;
1598 
1599  for (const auto & vj : variable_list)
1600  {
1601  std::vector<dof_id_type> di;
1602  this->dof_indices (partner, di, vj);
1603 
1604  // Insert the remote DOF indices into the send list
1605  for (auto d : di)
1606  if (!this->local_index(d))
1607  _send_list.push_back(d);
1608  }
1609  }
1610  else
1611  {
1612  std::vector<dof_id_type> di;
1613  this->dof_indices (partner, di);
1614 
1615  // Insert the remote DOF indices into the send list
1616  for (const auto & dof : di)
1617  if (!this->local_index(dof))
1618  _send_list.push_back(dof);
1619  }
1620 
1621  }
1622 
1623  // We're now done with any merged coupling matrices we had to create.
1624  for (auto & mat : temporary_coupling_matrices)
1625  delete mat;
1626 
1627  //-------------------------------------------------------------------------
1628  // Our coupling functors added dofs from neighboring elements to the
1629  // send list, but we may still need to add non-local dofs from local
1630  // elements.
1631  //-------------------------------------------------------------------------
1632 
1633  // Loop over the active local elements, adding all active elements
1634  // that neighbor an active local element to the send list.
1635  for ( ; local_elem_it != local_elem_end; ++local_elem_it)
1636  {
1637  const Elem * elem = *local_elem_it;
1638 
1639  std::vector<dof_id_type> di;
1640  this->dof_indices (elem, di);
1641 
1642  // Insert the remote DOF indices into the send list
1643  for (const auto & dof : di)
1644  if (!this->local_index(dof))
1645  _send_list.push_back(dof);
1646  }
1647 }
1648 
1649 
1650 
1652 {
1653  LOG_SCOPE("prepare_send_list()", "DofMap");
1654 
1655  // Return immediately if there's no ghost data
1656  if (this->n_processors() == 1)
1657  return;
1658 
1659  // Check to see if we have any extra stuff to add to the send_list
1661  {
1662  if (_augment_send_list)
1663  {
1664  libmesh_here();
1665  libMesh::out << "WARNING: You have specified both an extra send list function and object.\n"
1666  << " Are you sure this is what you meant to do??"
1667  << std::endl;
1668  }
1669 
1671  }
1672 
1673  if (_augment_send_list)
1675 
1676  // First sort the send list. After this
1677  // duplicated elements will be adjacent in the
1678  // vector
1679  std::sort(_send_list.begin(), _send_list.end());
1680 
1681  // Now use std::unique to remove duplicate entries
1682  std::vector<dof_id_type>::iterator new_end =
1683  std::unique (_send_list.begin(), _send_list.end());
1684 
1685  // Remove the end of the send_list. Use the "swap trick"
1686  // from Effective STL
1687  std::vector<dof_id_type> (_send_list.begin(), new_end).swap (_send_list);
1688 }
1689 
1691 {
1692  this->clear_send_list();
1693  this->add_neighbors_to_send_list(mesh);
1694 
1695 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1696  // This is assuming that we only need to recommunicate
1697  // the constraints and no new ones have been added since
1698  // a previous call to reinit_constraints.
1699  this->process_constraints(mesh);
1700 #endif
1701  this->prepare_send_list();
1702 }
1703 
1704 void DofMap::set_implicit_neighbor_dofs(bool implicit_neighbor_dofs)
1705 {
1707  _implicit_neighbor_dofs = implicit_neighbor_dofs;
1708 }
1709 
1710 
1712 {
1713  // If we were asked on the command line, then we need to
1714  // include sensitivities between neighbor degrees of freedom
1715  bool implicit_neighbor_dofs =
1716  libMesh::on_command_line ("--implicit-neighbor-dofs");
1717 
1718  // If the user specifies --implicit-neighbor-dofs 0, then
1719  // presumably he knows what he is doing and we won't try to
1720  // automatically turn it on even when all the variables are
1721  // discontinuous.
1722  if (implicit_neighbor_dofs)
1723  {
1724  // No flag provided defaults to 'true'
1725  int flag = 1;
1726  flag = libMesh::command_line_next ("--implicit-neighbor-dofs", flag);
1727 
1728  if (!flag)
1729  {
1730  // The user said --implicit-neighbor-dofs 0, so he knows
1731  // what he is doing and really doesn't want it.
1732  return false;
1733  }
1734  }
1735 
1736  // Possibly override the commandline option, if set_implicit_neighbor_dofs
1737  // has been called.
1739  {
1740  implicit_neighbor_dofs = _implicit_neighbor_dofs;
1741 
1742  // Again, if the user explicitly says implicit_neighbor_dofs = false,
1743  // then we return here.
1744  if (!implicit_neighbor_dofs)
1745  return false;
1746  }
1747 
1748  // Look at all the variables in this system. If every one is
1749  // discontinuous then the user must be doing DG/FVM, so be nice
1750  // and force implicit_neighbor_dofs=true.
1751  {
1752  bool all_discontinuous_dofs = true;
1753 
1754  for (auto var : IntRange<unsigned int>(0, this->n_variables()))
1755  if (FEAbstract::build (mesh.mesh_dimension(),
1756  this->variable_type(var))->get_continuity() != DISCONTINUOUS)
1757  all_discontinuous_dofs = false;
1758 
1759  if (all_discontinuous_dofs)
1760  implicit_neighbor_dofs = true;
1761  }
1762 
1763  return implicit_neighbor_dofs;
1764 }
1765 
1766 
1767 
1769 {
1770  _sp = this->build_sparsity(mesh);
1771 
1772  // It is possible that some \p SparseMatrix implementations want to
1773  // see the sparsity pattern before we throw it away. If so, we
1774  // share a view of its arrays, and we pass it in to the matrices.
1776  {
1777  _n_nz = &_sp->n_nz;
1778  _n_oz = &_sp->n_oz;
1779 
1780  for (const auto & mat : _matrices)
1781  mat->update_sparsity_pattern (_sp->sparsity_pattern);
1782  }
1783  // If we don't need the full sparsity pattern anymore, steal the
1784  // arrays we do need and free the rest of the memory
1785  else
1786  {
1787  if (!_n_nz)
1788  _n_nz = new std::vector<dof_id_type>();
1789  _n_nz->swap(_sp->n_nz);
1790  if (!_n_oz)
1791  _n_oz = new std::vector<dof_id_type>();
1792  _n_oz->swap(_sp->n_oz);
1793 
1794  _sp.reset();
1795  }
1796 }
1797 
1798 
1799 
1801 {
1803  {
1804  libmesh_assert(_sp.get());
1805  libmesh_assert(!_n_nz || _n_nz == &_sp->n_nz);
1806  libmesh_assert(!_n_oz || _n_oz == &_sp->n_oz);
1807  _sp.reset();
1808  }
1809  else
1810  {
1811  libmesh_assert(!_sp.get());
1812  delete _n_nz;
1813  delete _n_oz;
1814  }
1815  _n_nz = nullptr;
1816  _n_oz = nullptr;
1817 }
1818 
1819 
1820 
1822 {
1825 }
1826 
1827 
1828 
1830 {
1831  this->add_coupling_functor(this->default_coupling());
1833 }
1834 
1835 
1836 
1837 void
1839  bool to_mesh)
1840 {
1841  _coupling_functors.insert(&coupling_functor);
1842  if (to_mesh)
1843  _mesh.add_ghosting_functor(coupling_functor);
1844 }
1845 
1846 
1847 
1848 void
1850 {
1851  _coupling_functors.erase(&coupling_functor);
1852  _mesh.remove_ghosting_functor(coupling_functor);
1853 
1854  auto it = _shared_functors.find(&coupling_functor);
1855  if (it != _shared_functors.end())
1856  _shared_functors.erase(it);
1857 }
1858 
1859 
1860 
1861 void
1863  bool to_mesh)
1864 {
1865  _algebraic_ghosting_functors.insert(&evaluable_functor);
1866  if (to_mesh)
1867  _mesh.add_ghosting_functor(evaluable_functor);
1868 }
1869 
1870 
1871 
1872 void
1874 {
1875  _algebraic_ghosting_functors.erase(&evaluable_functor);
1876  _mesh.remove_ghosting_functor(evaluable_functor);
1877 
1878  auto it = _shared_functors.find(&evaluable_functor);
1879  if (it != _shared_functors.end())
1880  _shared_functors.erase(it);
1881 }
1882 
1883 
1884 
1886  const std::vector<dof_id_type> & dof_indices_in,
1887  DenseVectorBase<Number> & Ue) const
1888 {
1889  const unsigned int n_original_dofs = dof_indices_in.size();
1890 
1891 #ifdef LIBMESH_ENABLE_AMR
1892 
1893  // Trivial mapping
1894  libmesh_assert_equal_to (dof_indices_in.size(), Ue.size());
1895  bool has_constrained_dofs = false;
1896 
1897  for (unsigned int il=0; il != n_original_dofs; ++il)
1898  {
1899  const dof_id_type ig = dof_indices_in[il];
1900 
1901  if (this->is_constrained_dof (ig)) has_constrained_dofs = true;
1902 
1903  libmesh_assert_less (ig, Ug.size());
1904 
1905  Ue.el(il) = Ug(ig);
1906  }
1907 
1908  // If the element has any constrained DOFs then we need
1909  // to account for them in the mapping. This will handle
1910  // the case that the input vector is not constrained.
1911  if (has_constrained_dofs)
1912  {
1913  // Copy the input DOF indices.
1914  std::vector<dof_id_type> constrained_dof_indices(dof_indices_in);
1915 
1918 
1919  this->build_constraint_matrix_and_vector (C, H, constrained_dof_indices);
1920 
1921  libmesh_assert_equal_to (dof_indices_in.size(), C.m());
1922  libmesh_assert_equal_to (constrained_dof_indices.size(), C.n());
1923 
1924  // zero-out Ue
1925  Ue.zero();
1926 
1927  // compute Ue = C Ug, with proper mapping.
1928  for (unsigned int i=0; i != n_original_dofs; i++)
1929  {
1930  Ue.el(i) = H(i);
1931 
1932  const unsigned int n_constrained =
1933  cast_int<unsigned int>(constrained_dof_indices.size());
1934  for (unsigned int j=0; j<n_constrained; j++)
1935  {
1936  const dof_id_type jg = constrained_dof_indices[j];
1937 
1938  // If Ug is a serial or ghosted vector, then this assert is
1939  // overzealous. If Ug is a parallel vector, then this assert
1940  // is redundant.
1941  // libmesh_assert ((jg >= Ug.first_local_index()) &&
1942  // (jg < Ug.last_local_index()));
1943 
1944  Ue.el(i) += C(i,j)*Ug(jg);
1945  }
1946  }
1947  }
1948 
1949 #else
1950 
1951  // Trivial mapping
1952 
1953  libmesh_assert_equal_to (n_original_dofs, Ue.size());
1954 
1955  for (unsigned int il=0; il<n_original_dofs; il++)
1956  {
1957  const dof_id_type ig = dof_indices_in[il];
1958 
1959  libmesh_assert ((ig >= Ug.first_local_index()) && (ig < Ug.last_local_index()));
1960 
1961  Ue.el(il) = Ug(ig);
1962  }
1963 
1964 #endif
1965 }
1966 
1967 void DofMap::dof_indices (const Elem * const elem,
1968  std::vector<dof_id_type> & di) const
1969 {
1970  // We now allow elem==nullptr to request just SCALAR dofs
1971  // libmesh_assert(elem);
1972 
1973  // If we are asking for current indices on an element, it ought to
1974  // be an active element (or a Side proxy, which also thinks it's
1975  // active)
1976  libmesh_assert(!elem || elem->active());
1977 
1978  LOG_SCOPE("dof_indices()", "DofMap");
1979 
1980  // Clear the DOF indices vector
1981  di.clear();
1982 
1983  const unsigned int n_var_groups = this->n_variable_groups();
1984 
1985 #ifdef DEBUG
1986  // Check that sizes match in DEBUG mode
1987  std::size_t tot_size = 0;
1988 #endif
1989 
1990  if (elem && elem->type() == TRI3SUBDIVISION)
1991  {
1992  // Subdivision surface FE require the 1-ring around elem
1993  const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
1994 
1995  // Ghost subdivision elements have no real dofs
1996  if (!sd_elem->is_ghost())
1997  {
1998  // Determine the nodes contributing to element elem
1999  std::vector<const Node *> elem_nodes;
2000  MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2001 
2002  // Get the dof numbers
2003  for (unsigned int vg=0; vg<n_var_groups; vg++)
2004  {
2005  const VariableGroup & var = this->variable_group(vg);
2006  const unsigned int vars_in_group = var.n_variables();
2007 
2008  if (var.type().family == SCALAR &&
2009  var.active_on_subdomain(elem->subdomain_id()))
2010  {
2011  for (unsigned int vig=0; vig != vars_in_group; ++vig)
2012  {
2013 #ifdef DEBUG
2014  tot_size += var.type().order;
2015 #endif
2016  std::vector<dof_id_type> di_new;
2017  this->SCALAR_dof_indices(di_new,var.number(vig));
2018  di.insert( di.end(), di_new.begin(), di_new.end());
2019  }
2020  }
2021  else
2022  for (unsigned int vig=0; vig != vars_in_group; ++vig)
2023  {
2024  _dof_indices(*elem, elem->p_level(), di, vg, vig,
2025  elem_nodes.data(),
2026  cast_int<unsigned int>(elem_nodes.size())
2027 #ifdef DEBUG
2028  , var.number(vig), tot_size
2029 #endif
2030  );
2031  }
2032  }
2033  }
2034 
2035  return;
2036  }
2037 
2038  // Get the dof numbers for each variable
2039  const unsigned int n_nodes = elem ? elem->n_nodes() : 0;
2040  for (unsigned int vg=0; vg<n_var_groups; vg++)
2041  {
2042  const VariableGroup & var = this->variable_group(vg);
2043  const unsigned int vars_in_group = var.n_variables();
2044 
2045  if (var.type().family == SCALAR &&
2046  (!elem ||
2047  var.active_on_subdomain(elem->subdomain_id())))
2048  {
2049  for (unsigned int vig=0; vig != vars_in_group; ++vig)
2050  {
2051 #ifdef DEBUG
2052  tot_size += var.type().order;
2053 #endif
2054  std::vector<dof_id_type> di_new;
2055  this->SCALAR_dof_indices(di_new,var.number(vig));
2056  di.insert( di.end(), di_new.begin(), di_new.end());
2057  }
2058  }
2059  else if (elem)
2060  for (unsigned int vig=0; vig != vars_in_group; ++vig)
2061  {
2062  _dof_indices(*elem, elem->p_level(), di, vg, vig,
2063  elem->get_nodes(), n_nodes
2064 #ifdef DEBUG
2065  , var.number(vig), tot_size
2066 #endif
2067  );
2068  }
2069  }
2070 
2071 #ifdef DEBUG
2072  libmesh_assert_equal_to (tot_size, di.size());
2073 #endif
2074 }
2075 
2076 
2077 void DofMap::dof_indices (const Elem * const elem,
2078  std::vector<dof_id_type> & di,
2079  const unsigned int vn,
2080  int p_level) const
2081 {
2082  // We now allow elem==nullptr to request just SCALAR dofs
2083  // libmesh_assert(elem);
2084 
2085  LOG_SCOPE("dof_indices()", "DofMap");
2086 
2087  // Clear the DOF indices vector
2088  di.clear();
2089 
2090  // Use the default p refinement level?
2091  if (p_level == -12345)
2092  p_level = elem ? elem->p_level() : 0;
2093 
2094  const unsigned int vg = this->_variable_group_numbers[vn];
2095  const VariableGroup & var = this->variable_group(vg);
2096  const unsigned int vig = vn - var.number();
2097 
2098 #ifdef DEBUG
2099  // Check that sizes match in DEBUG mode
2100  std::size_t tot_size = 0;
2101 #endif
2102 
2103  if (elem && elem->type() == TRI3SUBDIVISION)
2104  {
2105  // Subdivision surface FE require the 1-ring around elem
2106  const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2107 
2108  // Ghost subdivision elements have no real dofs
2109  if (!sd_elem->is_ghost())
2110  {
2111  // Determine the nodes contributing to element elem
2112  std::vector<const Node *> elem_nodes;
2113  MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2114 
2115  _dof_indices(*elem, p_level, di, vg, vig, elem_nodes.data(),
2116  cast_int<unsigned int>(elem_nodes.size())
2117 #ifdef DEBUG
2118  , vn, tot_size
2119 #endif
2120  );
2121  }
2122 
2123  return;
2124  }
2125 
2126  // Get the dof numbers
2127  if (var.type().family == SCALAR &&
2128  (!elem ||
2129  var.active_on_subdomain(elem->subdomain_id())))
2130  {
2131 #ifdef DEBUG
2132  tot_size += var.type().order;
2133 #endif
2134  std::vector<dof_id_type> di_new;
2135  this->SCALAR_dof_indices(di_new,vn);
2136  di.insert( di.end(), di_new.begin(), di_new.end());
2137  }
2138  else if (elem)
2139  _dof_indices(*elem, p_level, di, vg, vig, elem->get_nodes(),
2140  elem->n_nodes()
2141 #ifdef DEBUG
2142  , vn, tot_size
2143 #endif
2144  );
2145 
2146 #ifdef DEBUG
2147  libmesh_assert_equal_to (tot_size, di.size());
2148 #endif
2149 }
2150 
2151 
2152 void DofMap::dof_indices (const Node * const node,
2153  std::vector<dof_id_type> & di) const
2154 {
2155  // We allow node==nullptr to request just SCALAR dofs
2156  // libmesh_assert(elem);
2157 
2158  LOG_SCOPE("dof_indices(Node)", "DofMap");
2159 
2160  // Clear the DOF indices vector
2161  di.clear();
2162 
2163  const unsigned int n_var_groups = this->n_variable_groups();
2164  const unsigned int sys_num = this->sys_number();
2165 
2166  // Get the dof numbers
2167  for (unsigned int vg=0; vg<n_var_groups; vg++)
2168  {
2169  const VariableGroup & var = this->variable_group(vg);
2170  const unsigned int vars_in_group = var.n_variables();
2171 
2172  if (var.type().family == SCALAR)
2173  {
2174  for (unsigned int vig=0; vig != vars_in_group; ++vig)
2175  {
2176  std::vector<dof_id_type> di_new;
2177  this->SCALAR_dof_indices(di_new,var.number(vig));
2178  di.insert( di.end(), di_new.begin(), di_new.end());
2179  }
2180  }
2181  else
2182  {
2183  const int n_comp = node->n_comp_group(sys_num,vg);
2184  for (unsigned int vig=0; vig != vars_in_group; ++vig)
2185  {
2186  for (int i=0; i != n_comp; ++i)
2187  {
2188  const dof_id_type d =
2189  node->dof_number(sys_num, vg, vig, i, n_comp);
2190  libmesh_assert_not_equal_to
2191  (d, DofObject::invalid_id);
2192  di.push_back(d);
2193  }
2194  }
2195  }
2196  }
2197 }
2198 
2199 
2200 void DofMap::dof_indices (const Node * const node,
2201  std::vector<dof_id_type> & di,
2202  const unsigned int vn) const
2203 {
2204  if (vn == libMesh::invalid_uint)
2205  {
2206  this->dof_indices(node, di);
2207  return;
2208  }
2209 
2210  // We allow node==nullptr to request just SCALAR dofs
2211  // libmesh_assert(elem);
2212 
2213  LOG_SCOPE("dof_indices(Node)", "DofMap");
2214 
2215  // Clear the DOF indices vector
2216  di.clear();
2217 
2218  const unsigned int sys_num = this->sys_number();
2219 
2220  // Get the dof numbers
2221  const unsigned int vg = this->_variable_group_numbers[vn];
2222  const VariableGroup & var = this->variable_group(vg);
2223 
2224  if (var.type().family == SCALAR)
2225  {
2226  std::vector<dof_id_type> di_new;
2227  this->SCALAR_dof_indices(di_new,vn);
2228  di.insert( di.end(), di_new.begin(), di_new.end());
2229  }
2230  else
2231  {
2232  const unsigned int vig = vn - var.number();
2233  const int n_comp = node->n_comp_group(sys_num,vg);
2234  for (int i=0; i != n_comp; ++i)
2235  {
2236  const dof_id_type d =
2237  node->dof_number(sys_num, vg, vig, i, n_comp);
2238  libmesh_assert_not_equal_to
2239  (d, DofObject::invalid_id);
2240  di.push_back(d);
2241  }
2242  }
2243 }
2244 
2245 
2246 void DofMap::dof_indices (const Elem & elem,
2247  unsigned int n,
2248  std::vector<dof_id_type> & di,
2249  const unsigned int vn) const
2250 {
2251  this->_node_dof_indices(elem, n, elem.node_ref(n), di, vn);
2252 }
2253 
2254 
2255 
2256 #ifdef LIBMESH_ENABLE_AMR
2257 
2258 void DofMap::old_dof_indices (const Elem & elem,
2259  unsigned int n,
2260  std::vector<dof_id_type> & di,
2261  const unsigned int vn) const
2262 {
2263  const DofObject * old_obj = elem.node_ref(n).old_dof_object;
2264  libmesh_assert(old_obj);
2265  this->_node_dof_indices(elem, n, *old_obj, di, vn);
2266 }
2267 
2268 #endif // LIBMESH_ENABLE_AMR
2269 
2270 
2271 
2272 void DofMap::_node_dof_indices (const Elem & elem,
2273  unsigned int n,
2274  const DofObject & obj,
2275  std::vector<dof_id_type> & di,
2276  const unsigned int vn) const
2277 {
2278  // Half of this is a cut and paste of _dof_indices code below, but
2279  // duplication actually seems cleaner than creating a helper
2280  // function with a million arguments and hoping the compiler inlines
2281  // it properly into one of our most highly trafficked functions.
2282 
2283  LOG_SCOPE("_node_dof_indices()", "DofMap");
2284 
2285  const ElemType type = elem.type();
2286  const unsigned int dim = elem.dim();
2287 
2288  const unsigned int sys_num = this->sys_number();
2289  const std::pair<unsigned int, unsigned int>
2290  vg_and_offset = obj.var_to_vg_and_offset(sys_num,vn);
2291  const unsigned int vg = vg_and_offset.first;
2292  const unsigned int vig = vg_and_offset.second;
2293  const unsigned int n_comp = obj.n_comp_group(sys_num,vg);
2294 
2295  const VariableGroup & var = this->variable_group(vg);
2296  FEType fe_type = var.type();
2297  fe_type.order = static_cast<Order>(fe_type.order +
2298  elem.p_level());
2299  const bool extra_hanging_dofs =
2301 
2302  // There is a potential problem with h refinement. Imagine a
2303  // quad9 that has a linear FE on it. Then, on the hanging side,
2304  // it can falsely identify a DOF at the mid-edge node. This is why
2305  // we go through FEInterface instead of obj->n_comp() directly.
2306  const unsigned int nc =
2307  FEInterface::n_dofs_at_node(dim, fe_type, type, n);
2308 
2309  // If this is a non-vertex on a hanging node with extra
2310  // degrees of freedom, we use the non-vertex dofs (which
2311  // come in reverse order starting from the end, to
2312  // simplify p refinement)
2313  if (extra_hanging_dofs && nc && !elem.is_vertex(n))
2314  {
2315  const int dof_offset = n_comp - nc;
2316 
2317  // We should never have fewer dofs than necessary on a
2318  // node unless we're getting indices on a parent element,
2319  // and we should never need the indices on such a node
2320  if (dof_offset < 0)
2321  {
2322  libmesh_assert(!elem.active());
2323  di.resize(di.size() + nc, DofObject::invalid_id);
2324  }
2325  else
2326  for (unsigned int i = dof_offset; i != n_comp; ++i)
2327  {
2328  const dof_id_type d =
2329  obj.dof_number(sys_num, vg, vig, i, n_comp);
2330  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2331  di.push_back(d);
2332  }
2333  }
2334  // If this is a vertex or an element without extra hanging
2335  // dofs, our dofs come in forward order coming from the
2336  // beginning
2337  else
2338  for (unsigned int i=0; i<nc; i++)
2339  {
2340  const dof_id_type d =
2341  obj.dof_number(sys_num, vg, vig, i, n_comp);
2342  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2343  di.push_back(d);
2344  }
2345 }
2346 
2347 
2348 
2349 void DofMap::_dof_indices (const Elem & elem,
2350  int p_level,
2351  std::vector<dof_id_type> & di,
2352  const unsigned int vg,
2353  const unsigned int vig,
2354  const Node * const * nodes,
2355  unsigned int n_nodes
2356 #ifdef DEBUG
2357  ,
2358  const unsigned int v,
2359  std::size_t & tot_size
2360 #endif
2361  ) const
2362 {
2363  const VariableGroup & var = this->variable_group(vg);
2364 
2365  if (var.active_on_subdomain(elem.subdomain_id()))
2366  {
2367  const ElemType type = elem.type();
2368  const unsigned int sys_num = this->sys_number();
2369  const unsigned int dim = elem.dim();
2370 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2371  const bool is_inf = elem.infinite();
2372 #endif
2373 
2374  // Increase the polynomial order on p refined elements
2375  FEType fe_type = var.type();
2376  fe_type.order = static_cast<Order>(fe_type.order + p_level);
2377 
2378  const bool extra_hanging_dofs =
2380 
2381 #ifdef DEBUG
2382  // The number of dofs per element is non-static for subdivision FE
2383  if (fe_type.family == SUBDIVISION)
2384  tot_size += n_nodes;
2385  else
2386  tot_size += FEInterface::n_dofs(dim,fe_type,type);
2387 #endif
2388 
2389  const FEInterface::n_dofs_at_node_ptr ndan =
2391 
2392  // Get the node-based DOF numbers
2393  for (unsigned int n=0; n != n_nodes; n++)
2394  {
2395  const Node & node = *nodes[n];
2396 
2397  // Cache the intermediate lookups that are common to every
2398  // component
2399 #ifdef DEBUG
2400  const std::pair<unsigned int, unsigned int>
2401  vg_and_offset = node.var_to_vg_and_offset(sys_num,v);
2402  libmesh_assert_equal_to (vg, vg_and_offset.first);
2403  libmesh_assert_equal_to (vig, vg_and_offset.second);
2404 #endif
2405  const unsigned int n_comp = node.n_comp_group(sys_num,vg);
2406 
2407  // There is a potential problem with h refinement. Imagine a
2408  // quad9 that has a linear FE on it. Then, on the hanging side,
2409  // it can falsely identify a DOF at the mid-edge node. This is why
2410  // we go through FEInterface instead of node.n_comp() directly.
2411  const unsigned int nc =
2412 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2413  is_inf ?
2414  FEInterface::n_dofs_at_node(dim, fe_type, type, n) :
2415 #endif
2416  ndan (type, fe_type.order, n);
2417 
2418  // If this is a non-vertex on a hanging node with extra
2419  // degrees of freedom, we use the non-vertex dofs (which
2420  // come in reverse order starting from the end, to
2421  // simplify p refinement)
2422  if (extra_hanging_dofs && !elem.is_vertex(n))
2423  {
2424  const int dof_offset = n_comp - nc;
2425 
2426  // We should never have fewer dofs than necessary on a
2427  // node unless we're getting indices on a parent element,
2428  // and we should never need the indices on such a node
2429  if (dof_offset < 0)
2430  {
2431  libmesh_assert(!elem.active());
2432  di.resize(di.size() + nc, DofObject::invalid_id);
2433  }
2434  else
2435  for (int i=n_comp-1; i>=dof_offset; i--)
2436  {
2437  const dof_id_type d =
2438  node.dof_number(sys_num, vg, vig, i, n_comp);
2439  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2440  di.push_back(d);
2441  }
2442  }
2443  // If this is a vertex or an element without extra hanging
2444  // dofs, our dofs come in forward order coming from the
2445  // beginning
2446  else
2447  for (unsigned int i=0; i<nc; i++)
2448  {
2449  const dof_id_type d =
2450  node.dof_number(sys_num, vg, vig, i, n_comp);
2451  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2452  di.push_back(d);
2453  }
2454  }
2455 
2456  // If there are any element-based DOF numbers, get them
2457  const unsigned int nc = FEInterface::n_dofs_per_elem(dim,
2458  fe_type,
2459  type);
2460  // We should never have fewer dofs than necessary on an
2461  // element unless we're getting indices on a parent element,
2462  // and we should never need those indices
2463  if (nc != 0)
2464  {
2465  const unsigned int n_comp = elem.n_comp_group(sys_num,vg);
2466  if (elem.n_systems() > sys_num && nc <= n_comp)
2467  {
2468  for (unsigned int i=0; i<nc; i++)
2469  {
2470  const dof_id_type d =
2471  elem.dof_number(sys_num, vg, vig, i, n_comp);
2472  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2473 
2474  di.push_back(d);
2475  }
2476  }
2477  else
2478  {
2479  libmesh_assert(!elem.active() || fe_type.family == LAGRANGE || fe_type.family == SUBDIVISION);
2480  di.resize(di.size() + nc, DofObject::invalid_id);
2481  }
2482  }
2483  }
2484 }
2485 
2486 
2487 
2488 void DofMap::SCALAR_dof_indices (std::vector<dof_id_type> & di,
2489  const unsigned int vn,
2490 #ifdef LIBMESH_ENABLE_AMR
2491  const bool old_dofs
2492 #else
2493  const bool
2494 #endif
2495  ) const
2496 {
2497  LOG_SCOPE("SCALAR_dof_indices()", "DofMap");
2498 
2499  libmesh_assert(this->variable(vn).type().family == SCALAR);
2500 
2501 #ifdef LIBMESH_ENABLE_AMR
2502  // If we're asking for old dofs then we'd better have some
2503  if (old_dofs)
2504  libmesh_assert_greater_equal(n_old_dofs(), n_SCALAR_dofs());
2505 
2506  dof_id_type my_idx = old_dofs ?
2507  this->_first_old_scalar_df[vn] : this->_first_scalar_df[vn];
2508 #else
2509  dof_id_type my_idx = this->_first_scalar_df[vn];
2510 #endif
2511 
2512  libmesh_assert_not_equal_to(my_idx, DofObject::invalid_id);
2513 
2514  // The number of SCALAR dofs comes from the variable order
2515  const int n_dofs_vn = this->variable(vn).type().order.get_order();
2516 
2517  di.resize(n_dofs_vn);
2518  for (int i = 0; i != n_dofs_vn; ++i)
2519  di[i] = my_idx++;
2520 }
2521 
2522 
2523 
2524 bool DofMap::semilocal_index (dof_id_type dof_index) const
2525 {
2526  // If it's not in the local indices
2527  if (!this->local_index(dof_index))
2528  {
2529  // and if it's not in the ghost indices, then we're not
2530  // semilocal
2531  if (!std::binary_search(_send_list.begin(), _send_list.end(), dof_index))
2532  return false;
2533  }
2534 
2535  return true;
2536 }
2537 
2538 
2539 
2540 bool DofMap::all_semilocal_indices (const std::vector<dof_id_type> & dof_indices_in) const
2541 {
2542  // We're all semilocal unless we find a counterexample
2543  for (const auto & di : dof_indices_in)
2544  if (!this->semilocal_index(di))
2545  return false;
2546 
2547  return true;
2548 }
2549 
2550 
2551 
2552 template <typename DofObjectSubclass>
2553 bool DofMap::is_evaluable(const DofObjectSubclass & obj,
2554  unsigned int var_num) const
2555 {
2556  // Everything is evaluable on a local object
2557  if (obj.processor_id() == this->processor_id())
2558  return true;
2559 
2560  std::vector<dof_id_type> di;
2561 
2562  if (var_num == libMesh::invalid_uint)
2563  this->dof_indices(&obj, di);
2564  else
2565  this->dof_indices(&obj, di, var_num);
2566 
2567  return this->all_semilocal_indices(di);
2568 }
2569 
2570 
2571 
2572 #ifdef LIBMESH_ENABLE_AMR
2573 
2574 void DofMap::old_dof_indices (const Elem * const elem,
2575  std::vector<dof_id_type> & di,
2576  const unsigned int vn) const
2577 {
2578  LOG_SCOPE("old_dof_indices()", "DofMap");
2579 
2580  libmesh_assert(elem);
2581 
2582  const ElemType type = elem->type();
2583  const unsigned int sys_num = this->sys_number();
2584  const unsigned int n_var_groups = this->n_variable_groups();
2585  const unsigned int dim = elem->dim();
2586 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2587  const bool is_inf = elem->infinite();
2588 #endif
2589 
2590  // If we have dof indices stored on the elem, and there's no chance
2591  // that we only have those indices because we were just p refined,
2592  // then we should have old dof indices too.
2593  libmesh_assert(!elem->has_dofs(sys_num) ||
2595  elem->old_dof_object);
2596 
2597  // Clear the DOF indices vector.
2598  di.clear();
2599 
2600  // Determine the nodes contributing to element elem
2601  std::vector<const Node *> elem_nodes;
2602  const Node * const * nodes_ptr;
2603  unsigned int n_nodes;
2604  if (elem->type() == TRI3SUBDIVISION)
2605  {
2606  // Subdivision surface FE require the 1-ring around elem
2607  const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
2608  MeshTools::Subdivision::find_one_ring(sd_elem, elem_nodes);
2609  nodes_ptr = elem_nodes.data();
2610  n_nodes = cast_int<unsigned int>(elem_nodes.size());
2611  }
2612  else
2613  {
2614  // All other FE use only the nodes of elem itself
2615  nodes_ptr = elem->get_nodes();
2616  n_nodes = elem->n_nodes();
2617  }
2618 
2619  // Get the dof numbers
2620  for (unsigned int vg=0; vg<n_var_groups; vg++)
2621  {
2622  const VariableGroup & var = this->variable_group(vg);
2623  const unsigned int vars_in_group = var.n_variables();
2624 
2625  for (unsigned int vig=0; vig<vars_in_group; vig++)
2626  {
2627  const unsigned int v = var.number(vig);
2628  if ((vn == v) || (vn == libMesh::invalid_uint))
2629  {
2630  if (var.type().family == SCALAR &&
2631  (!elem ||
2632  var.active_on_subdomain(elem->subdomain_id())))
2633  {
2634  // We asked for this variable, so add it to the vector.
2635  std::vector<dof_id_type> di_new;
2636  this->SCALAR_dof_indices(di_new,v,true);
2637  di.insert( di.end(), di_new.begin(), di_new.end());
2638  }
2639  else
2640  if (var.active_on_subdomain(elem->subdomain_id()))
2641  { // Do this for all the variables if one was not specified
2642  // or just for the specified variable
2643 
2644  // Increase the polynomial order on p refined elements,
2645  // but make sure you get the right polynomial order for
2646  // the OLD degrees of freedom
2647  int p_adjustment = 0;
2648  if (elem->p_refinement_flag() == Elem::JUST_REFINED)
2649  {
2650  libmesh_assert_greater (elem->p_level(), 0);
2651  p_adjustment = -1;
2652  }
2653  else if (elem->p_refinement_flag() == Elem::JUST_COARSENED)
2654  {
2655  p_adjustment = 1;
2656  }
2657  FEType fe_type = var.type();
2658  fe_type.order = static_cast<Order>(fe_type.order +
2659  elem->p_level() +
2660  p_adjustment);
2661 
2662  const bool extra_hanging_dofs =
2664 
2665  const FEInterface::n_dofs_at_node_ptr ndan =
2667 
2668  // Get the node-based DOF numbers
2669  for (unsigned int n=0; n<n_nodes; n++)
2670  {
2671  const Node * node = nodes_ptr[n];
2672  const DofObject * old_dof_obj = node->old_dof_object;
2673  libmesh_assert(old_dof_obj);
2674 
2675  // There is a potential problem with h refinement. Imagine a
2676  // quad9 that has a linear FE on it. Then, on the hanging side,
2677  // it can falsely identify a DOF at the mid-edge node. This is why
2678  // we call FEInterface instead of node->n_comp() directly.
2679  const unsigned int nc =
2680 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2681  is_inf ?
2682  FEInterface::n_dofs_at_node(dim, fe_type, type, n) :
2683 #endif
2684  ndan (type, fe_type.order, n);
2685 
2686  const int n_comp = old_dof_obj->n_comp_group(sys_num,vg);
2687 
2688  // If this is a non-vertex on a hanging node with extra
2689  // degrees of freedom, we use the non-vertex dofs (which
2690  // come in reverse order starting from the end, to
2691  // simplify p refinement)
2692  if (extra_hanging_dofs && !elem->is_vertex(n))
2693  {
2694  const int dof_offset = n_comp - nc;
2695 
2696  // We should never have fewer dofs than necessary on a
2697  // node unless we're getting indices on a parent element
2698  // or a just-coarsened element
2699  if (dof_offset < 0)
2700  {
2701  libmesh_assert(!elem->active() || elem->refinement_flag() ==
2703  di.resize(di.size() + nc, DofObject::invalid_id);
2704  }
2705  else
2706  for (int i=n_comp-1; i>=dof_offset; i--)
2707  {
2708  const dof_id_type d =
2709  old_dof_obj->dof_number(sys_num, vg, vig, i, n_comp);
2710  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2711  di.push_back(d);
2712  }
2713  }
2714  // If this is a vertex or an element without extra hanging
2715  // dofs, our dofs come in forward order coming from the
2716  // beginning
2717  else
2718  for (unsigned int i=0; i<nc; i++)
2719  {
2720  const dof_id_type d =
2721  old_dof_obj->dof_number(sys_num, vg, vig, i, n_comp);
2722  libmesh_assert_not_equal_to (d, DofObject::invalid_id);
2723  di.push_back(d);
2724  }
2725  }
2726 
2727  // If there are any element-based DOF numbers, get them
2728  const unsigned int nc = FEInterface::n_dofs_per_elem(dim,
2729  fe_type,
2730  type);
2731 
2732  // We should never have fewer dofs than necessary on an
2733  // element unless we're getting indices on a parent element
2734  // or a just-coarsened element
2735  if (nc != 0)
2736  {
2737  const DofObject * old_dof_obj = elem->old_dof_object;
2738  libmesh_assert(old_dof_obj);
2739 
2740  const unsigned int n_comp =
2741  old_dof_obj->n_comp_group(sys_num,vg);
2742 
2743  if (old_dof_obj->n_systems() > sys_num &&
2744  nc <= n_comp)
2745  {
2746 
2747  for (unsigned int i=0; i<nc; i++)
2748  {
2749  const dof_id_type d =
2750  old_dof_obj->dof_number(sys_num, vg, vig, i, n_comp);
2751 
2752  di.push_back(d);
2753  }
2754  }
2755  else
2756  {
2757  libmesh_assert(!elem->active() || fe_type.family == LAGRANGE ||
2759  di.resize(di.size() + nc, DofObject::invalid_id);
2760  }
2761  }
2762  }
2763  }
2764  } // end loop over variables within group
2765  } // end loop over variable groups
2766 }
2767 
2768 #endif // LIBMESH_ENABLE_AMR
2769 
2770 
2771 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2772 
2773 void DofMap::find_connected_dofs (std::vector<dof_id_type> & elem_dofs) const
2774 {
2775  typedef std::set<dof_id_type> RCSet;
2776 
2777  // First insert the DOFS we already depend on into the set.
2778  RCSet dof_set (elem_dofs.begin(), elem_dofs.end());
2779 
2780  bool done = true;
2781 
2782  // Next insert any dofs those might be constrained in terms
2783  // of. Note that in this case we may not be done: Those may
2784  // in turn depend on others. So, we need to repeat this process
2785  // in that case until the system depends only on unconstrained
2786  // degrees of freedom.
2787  for (const auto & dof : elem_dofs)
2788  if (this->is_constrained_dof(dof))
2789  {
2790  // If the DOF is constrained
2791  DofConstraints::const_iterator
2792  pos = _dof_constraints.find(dof);
2793 
2794  libmesh_assert (pos != _dof_constraints.end());
2795 
2796  const DofConstraintRow & constraint_row = pos->second;
2797 
2798  // adaptive p refinement currently gives us lots of empty constraint
2799  // rows - we should optimize those DoFs away in the future. [RHS]
2800  //libmesh_assert (!constraint_row.empty());
2801 
2802  // Add the DOFs this dof is constrained in terms of.
2803  // note that these dofs might also be constrained, so
2804  // we will need to call this function recursively.
2805  for (const auto & pr : constraint_row)
2806  if (!dof_set.count (pr.first))
2807  {
2808  dof_set.insert (pr.first);
2809  done = false;
2810  }
2811  }
2812 
2813 
2814  // If not done then we need to do more work
2815  // (obviously :-) )!
2816  if (!done)
2817  {
2818  // Fill the vector with the contents of the set
2819  elem_dofs.clear();
2820  elem_dofs.insert (elem_dofs.end(),
2821  dof_set.begin(), dof_set.end());
2822 
2823 
2824  // May need to do this recursively. It is possible
2825  // that we just replaced a constrained DOF with another
2826  // constrained DOF.
2827  this->find_connected_dofs (elem_dofs);
2828 
2829  } // end if (!done)
2830 }
2831 
2832 #endif // LIBMESH_ENABLE_CONSTRAINTS
2833 
2834 
2835 
2836 void DofMap::print_info(std::ostream & os) const
2837 {
2838  os << this->get_info();
2839 }
2840 
2841 
2842 
2843 std::string DofMap::get_info() const
2844 {
2845  std::ostringstream os;
2846 
2847  // If we didn't calculate the exact sparsity pattern, the threaded
2848  // sparsity pattern assembly may have just given us an upper bound
2849  // on sparsity.
2850  const char * may_equal = " <= ";
2851 
2852  // If we calculated the exact sparsity pattern, then we can report
2853  // exact bandwidth figures:
2854  for (const auto & mat : _matrices)
2855  if (mat->need_full_sparsity_pattern())
2856  may_equal = " = ";
2857 
2858  dof_id_type max_n_nz = 0, max_n_oz = 0;
2859  long double avg_n_nz = 0, avg_n_oz = 0;
2860 
2861  if (_n_nz)
2862  {
2863  for (const auto & val : *_n_nz)
2864  {
2865  max_n_nz = std::max(max_n_nz, val);
2866  avg_n_nz += val;
2867  }
2868 
2869  std::size_t n_nz_size = _n_nz->size();
2870 
2871  this->comm().max(max_n_nz);
2872  this->comm().sum(avg_n_nz);
2873  this->comm().sum(n_nz_size);
2874 
2875  avg_n_nz /= std::max(n_nz_size,std::size_t(1));
2876 
2878 
2879  for (const auto & val : *_n_oz)
2880  {
2881  max_n_oz = std::max(max_n_oz, val);
2882  avg_n_oz += val;
2883  }
2884 
2885  std::size_t n_oz_size = _n_oz->size();
2886 
2887  this->comm().max(max_n_oz);
2888  this->comm().sum(avg_n_oz);
2889  this->comm().sum(n_oz_size);
2890 
2891  avg_n_oz /= std::max(n_oz_size,std::size_t(1));
2892  }
2893 
2894  os << " DofMap Sparsity\n Average On-Processor Bandwidth"
2895  << may_equal << avg_n_nz << '\n';
2896 
2897  os << " Average Off-Processor Bandwidth"
2898  << may_equal << avg_n_oz << '\n';
2899 
2900  os << " Maximum On-Processor Bandwidth"
2901  << may_equal << max_n_nz << '\n';
2902 
2903  os << " Maximum Off-Processor Bandwidth"
2904  << may_equal << max_n_oz << std::endl;
2905 
2906 #ifdef LIBMESH_ENABLE_CONSTRAINTS
2907 
2908  std::size_t n_constraints = 0, max_constraint_length = 0,
2909  n_rhss = 0;
2910  long double avg_constraint_length = 0.;
2911 
2912  for (const auto & pr : _dof_constraints)
2913  {
2914  // Only count local constraints, then sum later
2915  const dof_id_type constrained_dof = pr.first;
2916  if (!this->local_index(constrained_dof))
2917  continue;
2918 
2919  const DofConstraintRow & row = pr.second;
2920  std::size_t rowsize = row.size();
2921 
2922  max_constraint_length = std::max(max_constraint_length,
2923  rowsize);
2924  avg_constraint_length += rowsize;
2925  n_constraints++;
2926 
2927  if (_primal_constraint_values.count(constrained_dof))
2928  n_rhss++;
2929  }
2930 
2931  this->comm().sum(n_constraints);
2932  this->comm().sum(n_rhss);
2933  this->comm().sum(avg_constraint_length);
2934  this->comm().max(max_constraint_length);
2935 
2936  os << " DofMap Constraints\n Number of DoF Constraints = "
2937  << n_constraints;
2938  if (n_rhss)
2939  os << '\n'
2940  << " Number of Heterogenous Constraints= " << n_rhss;
2941  if (n_constraints)
2942  {
2943  avg_constraint_length /= n_constraints;
2944 
2945  os << '\n'
2946  << " Average DoF Constraint Length= " << avg_constraint_length;
2947  }
2948 
2949 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
2950  std::size_t n_node_constraints = 0, max_node_constraint_length = 0,
2951  n_node_rhss = 0;
2952  long double avg_node_constraint_length = 0.;
2953 
2954  for (const auto & pr : _node_constraints)
2955  {
2956  // Only count local constraints, then sum later
2957  const Node * node = pr.first;
2958  if (node->processor_id() != this->processor_id())
2959  continue;
2960 
2961  const NodeConstraintRow & row = pr.second.first;
2962  std::size_t rowsize = row.size();
2963 
2964  max_node_constraint_length = std::max(max_node_constraint_length,
2965  rowsize);
2966  avg_node_constraint_length += rowsize;
2967  n_node_constraints++;
2968 
2969  if (pr.second.second != Point(0))
2970  n_node_rhss++;
2971  }
2972 
2973  this->comm().sum(n_node_constraints);
2974  this->comm().sum(n_node_rhss);
2975  this->comm().sum(avg_node_constraint_length);
2976  this->comm().max(max_node_constraint_length);
2977 
2978  os << "\n Number of Node Constraints = " << n_node_constraints;
2979  if (n_node_rhss)
2980  os << '\n'
2981  << " Number of Heterogenous Node Constraints= " << n_node_rhss;
2982  if (n_node_constraints)
2983  {
2984  avg_node_constraint_length /= n_node_constraints;
2985  os << "\n Maximum Node Constraint Length= " << max_node_constraint_length
2986  << '\n'
2987  << " Average Node Constraint Length= " << avg_node_constraint_length;
2988  }
2989 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
2990 
2991  os << std::endl;
2992 
2993 #endif // LIBMESH_ENABLE_CONSTRAINTS
2994 
2995  return os.str();
2996 }
2997 
2998 
2999 template bool DofMap::is_evaluable<Elem>(const Elem &, unsigned int) const;
3000 template bool DofMap::is_evaluable<Node>(const Node &, unsigned int) const;
3001 
3002 } // namespace libMesh
libMesh::DofMap::_adjoint_constraint_values
AdjointDofConstraintValues _adjoint_constraint_values
Definition: dof_map.h:1823
libMesh::DofMap::add_algebraic_ghosting_functor
void add_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor, bool to_mesh=true)
Adds a functor which can specify algebraic ghosting requirements for use with distributed vectors.
Definition: dof_map.C:1862
libMesh::DofMap::_primal_constraint_values
DofConstraintValueMap _primal_constraint_values
Definition: dof_map.h:1821
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::DofMap::compute_sparsity
void compute_sparsity(const MeshBase &)
Computes the sparsity pattern for the matrices corresponding to proc_id and sends that data to Linear...
Definition: dof_map.C:1768
libMesh::DofMap::elem_ptr
DofObject * elem_ptr(MeshBase &mesh, dof_id_type i) const
Definition: dof_map.C:325
libMesh::PeriodicBoundaries
We're using a class instead of a typedef to allow forward declarations and future flexibility.
Definition: periodic_boundaries.h:51
libMesh::DofMap::add_variable_group
void add_variable_group(const VariableGroup &var_group)
Definition: dof_map.C:257
libMesh::DofMap::is_evaluable
bool is_evaluable(const DofObjectSubclass &obj, unsigned int var_num=libMesh::invalid_uint) const
Definition: dof_map.C:2553
libMesh::OrderWrapper::get_order
int get_order() const
Explicitly request the order as an int.
Definition: fe_type.h:77
libMesh::invalid_uint
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:249
libMesh::DofMap::is_attached
bool is_attached(SparseMatrix< Number > &matrix)
Matrices should not be attached more than once.
Definition: dof_map.C:310
libMesh::DofMap::n_old_dofs
dof_id_type n_old_dofs() const
Definition: dof_map.h:1384
libMesh::FEInterface::extra_hanging_dofs
static bool extra_hanging_dofs(const FEType &fe_t)
Definition: fe_interface.C:1656
libMesh::DofMap::_mesh
MeshBase & _mesh
The mesh that system uses.
Definition: dof_map.h:1644
libMesh::VariableGroup::number
unsigned int number(unsigned int v) const
Definition: variable.h:255
libMesh::DofObject::n_systems
unsigned int n_systems() const
Definition: dof_object.h:861
libMesh::DofMap::_stashed_dof_constraints
DofConstraints _stashed_dof_constraints
Definition: dof_map.h:1819
libMesh::command_line_next
T command_line_next(std::string name, T default_value)
Use GetPot's search()/next() functions to get following arguments from the command line.
Definition: libmesh.C:963
libMesh::DofMap::prepare_send_list
void prepare_send_list()
Takes the _send_list vector (which may have duplicate entries) and sorts it.
Definition: dof_map.C:1651
libMesh::DofMap::_default_coupling
std::unique_ptr< DefaultCoupling > _default_coupling
The default coupling GhostingFunctor, used to implement standard libMesh sparsity pattern constructio...
Definition: dof_map.h:1713
libMesh::Elem::JUST_REFINED
Definition: elem.h:1172
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::DofObject::old_dof_object
DofObject * old_dof_object
This object on the last mesh.
Definition: dof_object.h:81
libMesh::DofMap::dof_indices
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1967
libMesh::MeshTools::Subdivision::find_one_ring
void find_one_ring(const Tri3Subdivision *elem, std::vector< const Node * > &nodes)
Determines the 1-ring of element elem, and writes it to the nodes vector.
Definition: mesh_subdivision_support.C:31
libMesh::DofMap::_periodic_boundaries
std::unique_ptr< PeriodicBoundaries > _periodic_boundaries
Data structure containing periodic boundaries.
Definition: dof_map.h:1839
libMesh::DofMap::_node_constraints
NodeConstraints _node_constraints
Data structure containing DofObject constraints.
Definition: dof_map.h:1830
libMesh::Elem::n_nodes
virtual unsigned int n_nodes() const =0
libMesh::DofMap::_first_old_df
std::vector< dof_id_type > _first_old_df
First old DOF index on processor p.
Definition: dof_map.h:1799
libMesh::DofMap::extract_local_vector
void extract_local_vector(const NumericVector< Number > &Ug, const std::vector< dof_id_type > &dof_indices, DenseVectorBase< Number > &Ue) const
Builds the local element vector Ue from the global vector Ug, accounting for any constrained degrees ...
Definition: dof_map.C:1885
libMesh::DofMap::n_dofs
dof_id_type n_dofs() const
Definition: dof_map.h:625
libMesh::DofMap::_first_old_scalar_df
std::vector< dof_id_type > _first_old_scalar_df
First old DOF index for SCALAR variable v, or garbage for non-SCALAR variable v.
Definition: dof_map.h:1810
libMesh::FEInterface::n_dofs_at_node_function
static n_dofs_at_node_ptr n_dofs_at_node_function(const unsigned int dim, const FEType &fe_t)
Definition: fe_interface.C:523
libMesh::DofMap::_extra_sparsity_function
void(* _extra_sparsity_function)(SparsityPattern::Graph &, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz, void *)
A function pointer to a function to call to add extra entries to the sparsity pattern.
Definition: dof_map.h:1683
libMesh::DofMap::~DofMap
~DofMap()
Destructor.
Definition: dof_map.C:206
libMesh::DofMap::need_full_sparsity_pattern
bool need_full_sparsity_pattern
Default false; set to true if any attached matrix requires a full sparsity pattern.
Definition: dof_map.h:1756
libMesh::DofMap::get_info
std::string get_info() const
Gets summary info about the sparsity bandwidth and constraints.
Definition: dof_map.C:2843
libMesh::FEInterface::n_dofs_at_node
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:503
libMesh::DofMap::build_sparsity
std::unique_ptr< SparsityPattern::Build > build_sparsity(const MeshBase &mesh) const
Builds a sparsity pattern.
Definition: dof_map.C:63
libMesh::DofMap::end_dof
dof_id_type end_dof() const
Definition: dof_map.h:695
libMesh::NumericVector::last_local_index
virtual numeric_index_type last_local_index() const =0
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::DofMap::_n_nz
std::vector< dof_id_type > * _n_nz
The number of on-processor nonzeros in my portion of the global matrix.
Definition: dof_map.h:1770
libMesh::DofMap::DofMap
DofMap(const unsigned int sys_number, MeshBase &mesh)
Constructor.
Definition: dof_map.C:135
libMesh::DofMap::_end_old_df
std::vector< dof_id_type > _end_old_df
Last old DOF index (plus 1) on processor p.
Definition: dof_map.h:1804
libMesh::DofMap::remove_coupling_functor
void remove_coupling_functor(GhostingFunctor &coupling_functor)
Removes a functor which was previously added to the set of coupling functors, from both this DofMap a...
Definition: dof_map.C:1849
libMesh::DofMap::distribute_local_dofs_var_major
void distribute_local_dofs_var_major(dof_id_type &next_free_dof, MeshBase &mesh)
Distributes the global degrees of freedom, for dofs on this processor.
Definition: dof_map.C:1310
libMesh::DofMap::use_coupled_neighbor_dofs
bool use_coupled_neighbor_dofs(const MeshBase &mesh) const
Tells other library functions whether or not this problem includes coupling between dofs in neighbori...
Definition: dof_map.C:1711
libMesh::DofMap::semilocal_index
bool semilocal_index(dof_id_type dof_index) const
Definition: dof_map.C:2524
libMesh::DofMap::set_error_on_cyclic_constraint
void set_error_on_cyclic_constraint(bool error_on_cyclic_constraint)
Specify whether or not we perform an extra (opt-mode enabled) check for constraint loops.
Definition: dof_map.C:243
libMesh::Elem::infinite
virtual bool infinite() const =0
libMesh::DofMap::coupling_functors_begin
std::set< GhostingFunctor * >::const_iterator coupling_functors_begin() const
Beginning of range of coupling functors.
Definition: dof_map.h:329
libMesh::DofMap::clear_send_list
void clear_send_list()
Clears the _send_list vector.
Definition: dof_map.h:470
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::DenseMatrix< Number >
libMesh::DofObject::set_n_comp_group
void set_n_comp_group(const unsigned int s, const unsigned int vg, const unsigned int ncomp)
Sets the number of components for VariableGroup vg of system s associated with this DofObject.
Definition: dof_object.C:410
libMesh::DofMap::n_variables
unsigned int n_variables() const
Definition: dof_map.h:592
libMesh::DofMap::_n_SCALAR_dofs
dof_id_type _n_SCALAR_dofs
The total number of SCALAR dofs associated to all SCALAR variables.
Definition: dof_map.h:1787
libMesh::DofObject::n_var_groups
unsigned int n_var_groups(const unsigned int s) const
Definition: dof_object.h:881
libMesh::GhostingFunctor::map_type
std::unordered_map< const Elem *, const CouplingMatrix * > map_type
What elements do we care about and what variables do we care about on each element?
Definition: ghosting_functor.h:171
libMesh::DenseVectorBase
Defines an abstract dense vector base class for use in Finite Element-type computations.
Definition: dof_map.h:71
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::DofMap::process_constraints
void process_constraints(MeshBase &)
Postprocesses any constrained degrees of freedom to be constrained only in terms of unconstrained dof...
Definition: dof_map_constraints.C:3334
libMesh::DofMap::_dof_constraints
DofConstraints _dof_constraints
Data structure containing DOF constraints.
Definition: dof_map.h:1819
libMesh::DenseMatrixBase::n
unsigned int n() const
Definition: dense_matrix_base.h:109
libMesh::Variable::type
const FEType & type() const
Definition: variable.h:119
libMesh::DofObject::set_vg_dof_base
void set_vg_dof_base(const unsigned int s, const unsigned int vg, const dof_id_type db)
VariableGroup DoF indices are indexed as id = base + var_in_vg*ncomp + comp This method allows for di...
Definition: dof_object.h:1209
libMesh::DofMap::_adjoint_dirichlet_boundaries
std::vector< DirichletBoundaries * > _adjoint_dirichlet_boundaries
Data structure containing Dirichlet functions.
Definition: dof_map.h:1859
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::DofMap::default_coupling
DefaultCoupling & default_coupling()
Default coupling functor.
Definition: dof_map.h:341
libMesh::DenseVectorBase::el
virtual T el(const unsigned int i) const =0
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::DofMap::_variable_group_numbers
std::vector< unsigned int > _variable_group_numbers
The variable group number for each variable.
Definition: dof_map.h:1634
libMesh::DofMap::_first_scalar_df
std::vector< dof_id_type > _first_scalar_df
First DOF index for SCALAR variable v, or garbage for non-SCALAR variable v.
Definition: dof_map.h:1667
libMesh::DofMap::n_variable_groups
unsigned int n_variable_groups() const
Definition: dof_map.h:584
libMesh::DofMap::invalidate_dofs
void invalidate_dofs(MeshBase &mesh) const
Invalidates all active DofObject dofs for this system.
Definition: dof_map.C:819
libMesh::FEInterface::n_dofs_at_node_ptr
unsigned int(* n_dofs_at_node_ptr)(const ElemType, const Order, const unsigned int)
Definition: fe_interface.h:126
libMesh::DofMap::is_constrained_dof
bool is_constrained_dof(const dof_id_type dof) const
Definition: dof_map.h:1962
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
libMesh::DofObject::processor_id
processor_id_type processor_id() const
Definition: dof_object.h:829
libMesh::MeshBase::elements_begin
virtual element_iterator elements_begin()=0
Iterate over all the elements in the Mesh.
libMesh::Elem::active
bool active() const
Definition: elem.h:2345
parallel_sync.h
libMesh::DofObject
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
Definition: dof_object.h:53
libMesh::SparseMatrix< Number >
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DirichletBoundaries
We're using a class instead of a typedef to allow forward declarations and future flexibility.
Definition: dirichlet_boundaries.h:207
libMesh::DofMap::_send_list
std::vector< dof_id_type > _send_list
A list containing all the global DOF indices that affect the solution on my processor.
Definition: dof_map.h:1673
libMesh::DofMap::_first_df
std::vector< dof_id_type > _first_df
First DOF index on processor p.
Definition: dof_map.h:1656
libMesh::SparseMatrix::need_full_sparsity_pattern
virtual bool need_full_sparsity_pattern() const
Definition: sparse_matrix.h:120
libMesh::NumericVector::size
virtual numeric_index_type size() const =0
libMesh::DofMap::old_dof_indices
void old_dof_indices(const Elem &elem, unsigned int n, std::vector< dof_id_type > &di, const unsigned int vn) const
Appends to the vector di the old global degree of freedom indices for elem.node_ref(n),...
Definition: dof_map.C:2258
libMesh::NumericVector< Number >
libMesh::TypeVector::size
auto size() const -> decltype(std::norm(T()))
Definition: type_vector.h:944
libMesh::DofMap::set_error_on_constraint_loop
void set_error_on_constraint_loop(bool error_on_constraint_loop)
Definition: dof_map.C:250
libMesh::CouplingMatrix::empty
bool empty() const
Definition: coupling_matrix.h:631
libMesh::DofMap::_n_oz
std::vector< dof_id_type > * _n_oz
The number of off-processor nonzeros in my portion of the global matrix; allocated similar to _n_nz.
Definition: dof_map.h:1776
libMesh::DenseMatrixBase::m
unsigned int m() const
Definition: dense_matrix_base.h:104
libMesh::DofMap::add_default_ghosting
void add_default_ghosting()
Add the default functor(s) for coupling and algebraic ghosting.
Definition: dof_map.C:1829
libMesh::DofObject::dof_number
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
Definition: dof_object.h:956
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::DofMap::reinit
void reinit(MeshBase &mesh)
Reinitialize the underlying data structures conformal to the current mesh.
Definition: dof_map.C:484
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::DofMap::distribute_dofs
void distribute_dofs(MeshBase &)
Distribute dofs on the current mesh.
Definition: dof_map.C:910
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::DofMap::AugmentSparsityPattern::augment_sparsity_pattern
virtual void augment_sparsity_pattern(SparsityPattern::Graph &sparsity, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz)=0
User-defined function to augment the sparsity pattern.
libMesh::MeshBase::remove_ghosting_functor
void remove_ghosting_functor(GhostingFunctor &ghosting_functor)
Removes a functor which was previously added to the set of ghosting functors.
Definition: mesh_base.C:450
libMesh::DofMap::local_variable_indices
void local_variable_indices(std::vector< dof_id_type > &idx, const MeshBase &mesh, unsigned int var_num) const
Fills an array of those dof indices which belong to the given variable number and live on the current...
Definition: dof_map.C:1084
libMesh::FEInterface::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:534
libMesh::DofMap::sys_number
unsigned int sys_number() const
Definition: dof_map.h:1876
libMesh::ParallelObject::n_processors
processor_id_type n_processors() const
Definition: parallel_object.h:100
libMesh::Elem::p_refinement_flag
RefinementState p_refinement_flag() const
Definition: elem.h:2630
libMesh::Tri3Subdivision::is_ghost
bool is_ghost() const
Definition: face_tri3_subdivision.h:116
libMesh::GhostingFunctor
This abstract base class defines the interface by which library code and user code can report associa...
Definition: ghosting_functor.h:153
libMesh::DofObject::n_comp
unsigned int n_comp(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:926
libMesh::Elem::get_nodes
const Node *const * get_nodes() const
Definition: elem.h:2001
libMesh::DofMap::_augment_sparsity_pattern
AugmentSparsityPattern * _augment_sparsity_pattern
Function object to call to add extra entries to the sparsity pattern.
Definition: dof_map.h:1678
libMesh::DofMap::set_implicit_neighbor_dofs
void set_implicit_neighbor_dofs(bool implicit_neighbor_dofs)
Allow the implicit_neighbor_dofs flag to be set programmatically.
Definition: dof_map.C:1704
libMesh::DofMap::_extra_sparsity_context
void * _extra_sparsity_context
A pointer associated with the extra sparsity that can optionally be passed in.
Definition: dof_map.h:1690
libMesh::DofMap::variable
const Variable & variable(const unsigned int c) const
Definition: dof_map.h:1894
libMesh::DofMap::clear
void clear()
Free all new memory associated with the object, but restore its original state, with the mesh pointer...
Definition: dof_map.C:834
libMesh::FEAbstract::build
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
Definition: fe_abstract.C:72
libMesh::ParallelObject::processor_id
processor_id_type processor_id() const
Definition: parallel_object.h:106
libMesh::DofObject::invalid_id
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:421
libMesh::DofMap::_n_dfs
dof_id_type _n_dfs
Total number of degrees of freedom.
Definition: dof_map.h:1781
parallel_implementation.h
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::DofMap::SCALAR_dof_indices
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
Fills the vector di with the global degree of freedom indices corresponding to the SCALAR variable vn...
Definition: dof_map.C:2488
libMesh::Variable
This class defines the notion of a variable in the system.
Definition: variable.h:49
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::DISCONTINUOUS
Definition: enum_fe_family.h:78
libMesh::DofObject::var_to_vg_and_offset
std::pair< unsigned int, unsigned int > var_to_vg_and_offset(const unsigned int s, const unsigned int var) const
Definition: dof_object.h:1112
libMesh::DofMap::_n_old_dfs
dof_id_type _n_old_dfs
Total number of degrees of freedom on old dof objects.
Definition: dof_map.h:1794
libMesh::DofMap::set_nonlocal_dof_objects
void set_nonlocal_dof_objects(iterator_type objects_begin, iterator_type objects_end, MeshBase &mesh, dofobject_accessor objects)
Helper function for distributing dofs in parallel.
Definition: dof_map.C:333
libMesh::DofMap::find_connected_dofs
void find_connected_dofs(std::vector< dof_id_type > &elem_dofs) const
Finds all the DOFS associated with the element DOFs elem_dofs.
Definition: dof_map.C:2773
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::DofMap::_coupling_functors
std::set< GhostingFunctor * > _coupling_functors
The list of all GhostingFunctor objects to be used when coupling degrees of freedom in matrix sparsit...
Definition: dof_map.h:1744
libMesh::MeshTools::Generation::Private::idx
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
A useful inline function which replaces the macros used previously.
Definition: mesh_generation.C:72
libMesh::DofMap::build_constraint_matrix_and_vector
void build_constraint_matrix_and_vector(DenseMatrix< Number > &C, DenseVector< Number > &H, std::vector< dof_id_type > &elem_dofs, int qoi_index=-1, const bool called_recursively=false) const
Build the constraint matrix C and the forcing vector H associated with the element degree of freedom ...
Definition: dof_map_constraints.C:2561
libMesh::Threads::parallel_reduce
void parallel_reduce(const Range &range, Body &body)
Execute the provided reduction operation in parallel on the specified range.
Definition: threads_none.h:101
libMesh::FEInterface::n_dofs
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:472
libMesh::DofMap::variable_type
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:1924
libMesh::as_range
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
libMesh::DofMap::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints summary info about the sparsity bandwidth and constraints.
Definition: dof_map.C:2836
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::Elem::JUST_COARSENED
Definition: elem.h:1173
libMesh::DofMap::algebraic_ghosting_functors_end
std::set< GhostingFunctor * >::const_iterator algebraic_ghosting_functors_end() const
End of range of algebraic ghosting functors.
Definition: dof_map.h:397
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::DofMap::algebraic_ghosting_functors_begin
std::set< GhostingFunctor * >::const_iterator algebraic_ghosting_functors_begin() const
Beginning of range of algebraic ghosting functors.
Definition: dof_map.h:391
libMesh::MeshBase::const_element_iterator
The definition of the const_element_iterator struct.
Definition: mesh_base.h:1891
libMesh::DofObject::n_vars
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:891
libMesh::DenseVectorBase::size
virtual unsigned int size() const =0
libMesh::DofMap::_augment_send_list
AugmentSendList * _augment_send_list
Function object to call to add extra entries to the send list.
Definition: dof_map.h:1695
libMesh::DofMap::_variable_groups
std::vector< VariableGroup > _variable_groups
The finite element type for each variable group.
Definition: dof_map.h:1629
libMesh::MeshTools::libmesh_assert_valid_dof_ids
void libmesh_assert_valid_dof_ids(const MeshBase &mesh, unsigned int sysnum=libMesh::invalid_uint)
A function for verifying that degree of freedom indexing matches across processors.
Definition: mesh_tools.C:1566
libMesh::DofMap::reinit_send_list
void reinit_send_list(MeshBase &mesh)
Clears the _send_list vector and then rebuilds it.
Definition: dof_map.C:1690
libMesh::DefaultCoupling
This class implements the default algebraic coupling in libMesh: elements couple to themselves,...
Definition: default_coupling.h:41
libMesh::DofMap::remove_algebraic_ghosting_functor
void remove_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor)
Removes a functor which was previously added to the set of algebraic ghosting functors,...
Definition: dof_map.C:1873
libMesh::DofMap::n_dofs_on_processor
dof_id_type n_dofs_on_processor(const processor_id_type proc) const
Definition: dof_map.h:641
libMesh::DofMap::_matrices
std::vector< SparseMatrix< Number > * > _matrices
Additional matrices handled by this object.
Definition: dof_map.h:1651
libMesh::MeshBase::add_ghosting_functor
void add_ghosting_functor(GhostingFunctor &ghosting_functor)
Adds a functor which can specify ghosting requirements for use on distributed meshes.
Definition: mesh_base.h:1089
libMesh::Elem::refinement_flag
RefinementState refinement_flag() const
Definition: elem.h:2614
libMesh::DofMap::attach_matrix
void attach_matrix(SparseMatrix< Number > &matrix)
Additional matrices may be handled with this DofMap.
Definition: dof_map.C:274
libMesh::DofMap::local_index
bool local_index(dof_id_type dof_index) const
Definition: dof_map.h:815
libMesh::DofMap::remove_default_ghosting
void remove_default_ghosting()
Remove any default ghosting functor(s).
Definition: dof_map.C:1821
swap
void swap(Iterator &lhs, Iterator &rhs)
swap, used to implement op=
Definition: variant_filter_iterator.h:478
libMesh::SparseMatrix::update_sparsity_pattern
virtual void update_sparsity_pattern(const SparsityPattern::Graph &)
Updates the matrix sparsity pattern.
Definition: sparse_matrix.h:128
libMesh::DofMap::coupling_functors_end
std::set< GhostingFunctor * >::const_iterator coupling_functors_end() const
End of range of coupling functors.
Definition: dof_map.h:335
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::MeshBase::nodes_begin
virtual node_iterator nodes_begin()=0
Iterate over all the nodes in the Mesh.
libMesh::FEType::order
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:197
libMesh::DofMap::_algebraic_ghosting_functors
std::set< GhostingFunctor * > _algebraic_ghosting_functors
The list of all GhostingFunctor objects to be used when distributing ghosted vectors.
Definition: dof_map.h:1731
libMesh::DofMap::all_semilocal_indices
bool all_semilocal_indices(const std::vector< dof_id_type > &dof_indices) const
Definition: dof_map.C:2540
libMesh::DofMap::clear_sparsity
void clear_sparsity()
Clears the sparsity pattern.
Definition: dof_map.C:1800
libMesh::Elem::subdomain_id
subdomain_id_type subdomain_id() const
Definition: elem.h:2069
libMesh::DofObject::id
dof_id_type id() const
Definition: dof_object.h:767
libMesh::DofMap::_end_df
std::vector< dof_id_type > _end_df
Last DOF index (plus 1) on processor p.
Definition: dof_map.h:1661
libMesh::StoredRange
The StoredRange class defines a contiguous, divisible set of objects.
Definition: stored_range.h:52
libMesh::Elem::is_vertex
virtual bool is_vertex(const unsigned int i) const =0
libMesh::SparseMatrix::attach_dof_map
void attach_dof_map(const DofMap &dof_map)
Get a pointer to the DofMap to use.
Definition: sparse_matrix.h:108
libMesh::DofObject::invalid_processor_id
static const processor_id_type invalid_processor_id
An invalid processor_id to distinguish DoFs that have not been assigned to a processor.
Definition: dof_object.h:432
libMesh::DofMap::_shared_functors
std::map< GhostingFunctor *, std::shared_ptr< GhostingFunctor > > _shared_functors
Hang on to references to any GhostingFunctor objects we were passed in shared_ptr form.
Definition: dof_map.h:1750
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::DenseVectorBase::zero
virtual void zero()=0
Set every element in the vector to 0.
libMesh::DofMap::node_ptr
DofObject * node_ptr(MeshBase &mesh, dof_id_type i) const
Definition: dof_map.C:318
libMesh::DofMap::distribute_local_dofs_node_major
void distribute_local_dofs_node_major(dof_id_type &next_free_dof, MeshBase &mesh)
Distributes the global degrees of freedom for dofs on this processor.
Definition: dof_map.C:1171
libMesh::on_command_line
bool on_command_line(std::string arg)
Definition: libmesh.C:898
libMesh::DofMap::AugmentSendList::augment_send_list
virtual void augment_send_list(std::vector< dof_id_type > &send_list)=0
User-defined function to augment the send list.
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::DofMap::is_periodic_boundary
bool is_periodic_boundary(const boundary_id_type boundaryid) const
Definition: dof_map.C:224
libMesh::DofObject::n_comp_group
unsigned int n_comp_group(const unsigned int s, const unsigned int vg) const
Definition: dof_object.h:939
libMesh::DofMap::default_algebraic_ghosting
DefaultCoupling & default_algebraic_ghosting()
Default algebraic ghosting functor.
Definition: dof_map.h:403
libMesh::DofMap::_default_evaluating
std::unique_ptr< DefaultCoupling > _default_evaluating
The default algebraic GhostingFunctor, used to implement standard libMesh send_list construction.
Definition: dof_map.h:1721
libMesh::ConstCouplingRow
This proxy class acts like a container of indices from a single coupling row.
Definition: coupling_matrix.h:337
libMesh::DofMap::_implicit_neighbor_dofs
bool _implicit_neighbor_dofs
Definition: dof_map.h:1869
libMesh::DofMap::_dof_coupling
CouplingMatrix * _dof_coupling
Degree of freedom coupling.
Definition: dof_map.h:1434
libMesh::err
OStreamProxy err
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::CouplingMatrix
This class defines a coupling matrix.
Definition: coupling_matrix.h:54
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::SUBDIVISION
Definition: enum_fe_family.h:55
libMesh::DofObject::has_dofs
bool has_dofs(const unsigned int s=libMesh::invalid_uint) const
Definition: dof_object.h:1138
libMesh::Elem::node_ref
const Node & node_ref(const unsigned int i) const
Definition: elem.h:2031
libMesh::TRI3SUBDIVISION
Definition: enum_elem_type.h:69
libMesh::VariableGroup::n_variables
unsigned int n_variables() const
Definition: variable.h:217
libMesh::Variable::active_on_subdomain
bool active_on_subdomain(subdomain_id_type sid) const
Definition: variable.h:136
libMesh::DofMap::n_SCALAR_dofs
dof_id_type n_SCALAR_dofs() const
Definition: dof_map.h:630
libMesh::ParallelObject
An object whose state is distributed along a set of processors.
Definition: parallel_object.h:55
libMesh::DofMap::_error_on_constraint_loop
bool _error_on_constraint_loop
This flag indicates whether or not we do an opt-mode check for the presence of constraint loops,...
Definition: dof_map.h:1619
libMesh::DofMap::_variables
std::vector< Variable > _variables
The finite element type for each variable.
Definition: dof_map.h:1624
libMesh::VariableGroup
This class defines a logically grouped set of variables in the system.
Definition: variable.h:172
libMesh::Tri3Subdivision
The Tri3Subdivision element is a three-noded subdivision surface shell element used in mechanics calc...
Definition: face_tri3_subdivision.h:40
libMesh::DofMap::_dof_indices
void _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 #ifdef DEBUG, const unsigned int v, std::size_t &tot_size #endif) const
Helper function that gets the dof indices on the current element for a non-SCALAR type variable,...
Definition: dof_map.C:2349
libMesh::out
OStreamProxy out
libMesh::DofMap::first_dof
dof_id_type first_dof() const
Definition: dof_map.h:653
libMesh::DofObject::vg_dof_base
dof_id_type vg_dof_base(const unsigned int s, const unsigned int vg) const
VariableGroup DoF indices are indexed as id = base + var_in_vg*ncomp + comp This method allows for di...
Definition: dof_object.h:1229
libMesh::DofConstraintRow
std::map< dof_id_type, Real, std::less< dof_id_type >, Threads::scalable_allocator< std::pair< const dof_id_type, Real > > > DofConstraintRow
A row of the Dof constraint matrix.
Definition: dof_map.h:97
libMesh::FEInterface::max_order
static unsigned int max_order(const FEType &fe_t, const ElemType &el_t)
Definition: fe_interface.C:1281
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::DofMap::_sp
std::unique_ptr< SparsityPattern::Build > _sp
The sparsity pattern of the global matrix, kept around if it might be needed by future additions of t...
Definition: dof_map.h:1762
libMesh::DofMap::_implicit_neighbor_dofs_initialized
bool _implicit_neighbor_dofs_initialized
Bools to indicate if we override the –implicit_neighbor_dofs commandline options.
Definition: dof_map.h:1868
libMesh::DenseVector< Number >
libMesh::DofMap::_node_dof_indices
void _node_dof_indices(const Elem &elem, unsigned int n, const DofObject &obj, std::vector< dof_id_type > &di, const unsigned int vn) const
Helper function that implements the element-nodal versions of dof_indices and old_dof_indices.
Definition: dof_map.C:2272
libMesh::DofMap::_extra_send_list_context
void * _extra_send_list_context
A pointer associated with the extra send list that can optionally be passed in.
Definition: dof_map.h:1705
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::DofMap::variable_group
const VariableGroup & variable_group(const unsigned int c) const
Definition: dof_map.h:1884
libMesh::DofMap::_extra_send_list_function
void(* _extra_send_list_function)(std::vector< dof_id_type > &, void *)
A function pointer to a function to call to add extra entries to the send list.
Definition: dof_map.h:1700
libMesh::DofMap::merge_ghost_functor_outputs
static void 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)
Definition: dof_map.C:1440
libMesh::NumericVector::first_local_index
virtual numeric_index_type first_local_index() const =0
libMesh::NodeConstraintRow
std::map< const Node *, Real, std::less< const Node * >, Threads::scalable_allocator< std::pair< const Node *const, Real > > > NodeConstraintRow
A row of the Node constraint mapping.
Definition: dof_map.h:145
libMesh::DofMap::add_neighbors_to_send_list
void add_neighbors_to_send_list(MeshBase &mesh)
Adds entries to the _send_list vector corresponding to DoFs on elements neighboring the current proce...
Definition: dof_map.C:1507
libMesh::DofMap::add_coupling_functor
void add_coupling_functor(GhostingFunctor &coupling_functor, bool to_mesh=true)
Adds a functor which can specify coupling requirements for creation of sparse matrices.
Definition: dof_map.C:1838