21 #include "libmesh/dof_map.h" 
   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"  
   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"  
   55 #include <unordered_map> 
   62 std::unique_ptr<SparsityPattern::Build>
 
   67   LOG_SCOPE(
"build_sparsity()", 
"DofMap");
 
   90   auto sp = libmesh_make_unique<SparsityPattern::Build>
 
   95      implicit_neighbor_dofs,
 
   99                                             mesh.active_local_elements_end()), *sp);
 
  108   libmesh_assert_equal_to (sp->sparsity_pattern.size(), n_dofs_on_proc);
 
  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??" 
  122         (sp->sparsity_pattern, sp->n_nz,
 
  128       (sp->sparsity_pattern, sp->n_nz, sp->n_oz);
 
  130   return std::unique_ptr<SparsityPattern::Build>(sp.release());
 
  138   _dof_coupling(nullptr),
 
  139   _error_on_constraint_loop(false),
 
  142   _variable_group_numbers(),
 
  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),
 
  158   need_full_sparsity_pattern(false),
 
  163 #ifdef LIBMESH_ENABLE_AMR
 
  167   _first_old_scalar_df()
 
  169 #ifdef LIBMESH_ENABLE_CONSTRAINTS
 
  171   , _stashed_dof_constraints()
 
  172   , _primal_constraint_values()
 
  173   , _adjoint_constraint_values()
 
  175 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
 
  176   , _node_constraints()
 
  178 #ifdef LIBMESH_ENABLE_PERIODIC
 
  181 #ifdef LIBMESH_ENABLE_DIRICHLET
 
  183   , _adjoint_dirichlet_boundaries()
 
  185   , _implicit_neighbor_dofs_initialized(false),
 
  186   _implicit_neighbor_dofs(false)
 
  194 #ifdef LIBMESH_ENABLE_PERIODIC 
  215 #ifdef LIBMESH_ENABLE_DIRICHLET 
  222 #ifdef LIBMESH_ENABLE_PERIODIC 
  276   parallel_object_only();
 
  289   bool computed_sparsity_already =
 
  292   this->
comm().max(computed_sparsity_already);
 
  293   if (computed_sparsity_already &&
 
  320   return mesh.node_ptr(i);
 
  327   return mesh.elem_ptr(i);
 
  332 template <
typename iterator_type>
 
  334                                       iterator_type objects_end,
 
  336                                       dofobject_accessor objects)
 
  339   parallel_object_only();
 
  343   std::unordered_map<processor_id_type, dof_id_type> ghost_objects_from_proc;
 
  345   iterator_type it  = objects_begin;
 
  347   for (; it != objects_end; ++it)
 
  356           ghost_objects_from_proc[obj_procid]++;
 
  361   std::map<processor_id_type, std::vector<dof_id_type>>
 
  366   for (
auto pair : ghost_objects_from_proc)
 
  370         requested_ids[p].reserve(pair.second);
 
  373   for (it = objects_begin; it != objects_end; ++it)
 
  382       if (ghost_objects_from_proc.count(p))
 
  383         libmesh_assert_equal_to (requested_ids[p].size(), ghost_objects_from_proc[p]);
 
  389   typedef std::vector<dof_id_type> datum;
 
  391   auto gather_functor =
 
  392     [
this, &
mesh, &objects]
 
  394      const std::vector<dof_id_type> & ids,
 
  395      std::vector<datum> & 
data)
 
  402       const std::size_t query_size = ids.size();
 
  404       data.resize(query_size);
 
  405       for (
auto & d : 
data)
 
  406         d.resize(2 * n_var_groups);
 
  408       for (std::size_t i=0; i != query_size; ++i)
 
  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)
 
  416               unsigned int n_comp_g =
 
  418               data[i][vg] = n_comp_g;
 
  422               data[i][n_var_groups+vg] = my_first_dof;
 
  427   auto action_functor =
 
  428     [
this, &
mesh, &objects]
 
  430      const std::vector<dof_id_type> & ids,
 
  431      const std::vector<datum> & 
data)
 
  442           libmesh_assert_equal_to (requested->
processor_id(), pid);
 
  443           for (
unsigned int vg=0; vg != n_var_groups; ++vg)
 
  445               unsigned int n_comp_g =
 
  446                 cast_int<unsigned int>(
data[i][vg]);
 
  453                     (sys_num, vg, my_first_dof);
 
  459   datum * ex = 
nullptr;
 
  460   Parallel::pull_parallel_vector_data
 
  461     (this->
comm(), requested_ids, gather_functor, action_functor, ex);
 
  465   for (it = objects_begin; it != objects_end; ++it)
 
  470       for (
unsigned int v=0; v != num_variables; ++v)
 
  472           unsigned int n_comp =
 
  488   LOG_SCOPE(
"reinit()", 
"DofMap");
 
  503   unsigned int standard_n_levels =
 
  519   std::vector<unsigned int> n_vars_per_group;  n_vars_per_group.reserve (n_var_groups);
 
  521   for (
unsigned int vg=0; vg<n_var_groups; vg++)
 
  522     n_vars_per_group.push_back (this->variable_group(vg).n_variables());
 
  524 #ifdef LIBMESH_ENABLE_AMR 
  529   for (
auto & node : 
mesh.node_ptr_range())
 
  531       node->clear_old_dof_object();
 
  535   for (
auto & elem : 
mesh.element_ptr_range())
 
  537       elem->clear_old_dof_object();
 
  546   for (
auto & elem : 
mesh.element_ptr_range())
 
  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();
 
  559       if (elem->has_dofs(sys_num))
 
  560         elem->set_old_dof_object();
 
  563 #endif // #ifdef LIBMESH_ENABLE_AMR 
  572   for (
auto & node : 
mesh.node_ptr_range())
 
  573     node->set_n_vars_per_group(sys_num, n_vars_per_group);
 
  576   for (
auto & elem : 
mesh.element_ptr_range())
 
  577     elem->set_n_vars_per_group(sys_num, n_vars_per_group);
 
  584   for (
unsigned int vg=0; vg<n_var_groups; vg++)
 
  588       const unsigned int n_var_in_group = vg_description.
n_variables();
 
  589       const FEType & base_fe_type        = vg_description.
type();
 
  600       const bool extra_hanging_dofs =
 
  604       for (
auto & elem : 
mesh.active_element_ptr_range())
 
  614           const unsigned int dim = elem->dim();
 
  616           FEType fe_type = base_fe_type;
 
  618 #ifdef LIBMESH_ENABLE_AMR 
  621           if (elem->p_level() + base_fe_type.
order >
 
  626                                       "ERROR: Finite element " 
  628                                       << 
" on geometric element " 
  630                                       << 
"\nonly supports FEInterface::max_order = " 
  632                                       << 
", not fe_type.order = " 
  633                                       << base_fe_type.
order);
 
  638                            << 
" on geometric element " 
  640                            << 
"could not be p refined past FEInterface::max_order = " 
  645                                 - base_fe_type.
order);
 
  649           fe_type.
order = static_cast<Order>(fe_type.
order +
 
  653           for (
auto n : elem->node_index_range())
 
  655               Node & node = elem->node_ref(n);
 
  657               if (elem->is_vertex(n))
 
  659                   const unsigned int old_node_dofs =
 
  662                   const unsigned int vertex_dofs =
 
  668                   if (vertex_dofs > old_node_dofs)
 
  693       for (
auto & elem : 
mesh.active_element_ptr_range())
 
  703           const unsigned int dim = elem->dim();
 
  705           FEType fe_type = base_fe_type;
 
  706           fe_type.
order = static_cast<Order>(fe_type.
order +
 
  710           for (
auto n : elem->node_index_range())
 
  712               Node & node = elem->node_ref(n);
 
  714               const unsigned int old_node_dofs =
 
  717               const unsigned int vertex_dofs = old_node_dofs?
 
  718                 cast_int<unsigned int>(node.
vg_dof_base (sys_num,vg)):0;
 
  720               const unsigned int new_node_dofs =
 
  724               if (elem->is_vertex(n))
 
  726                   libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
 
  736                   libmesh_assert_greater_equal (vertex_dofs,   new_node_dofs);
 
  756                   else if (vertex_dofs == 0)
 
  758                       if (new_node_dofs > old_node_dofs)
 
  770                   else if (extra_hanging_dofs)
 
  772                       if (new_node_dofs > old_node_dofs - vertex_dofs)
 
  775                                                 vertex_dofs + new_node_dofs);
 
  785                       libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
 
  786                       if (new_node_dofs > old_node_dofs)
 
  798           const unsigned int dofs_per_elem =
 
  802           elem->set_n_comp_group(sys_num, vg, dofs_per_elem);
 
  821   const unsigned int sys_num = this->
sys_number();
 
  824   for (
auto & node : 
mesh.node_ptr_range())
 
  825     node->invalidate_dofs(sys_num);
 
  828   for (
auto & elem : 
mesh.active_element_ptr_range())
 
  829     elem->invalidate_dofs(sys_num);
 
  852     this->_coupling_functors.clear();
 
  869     this->_algebraic_ghosting_functors.clear();
 
  890 #ifdef LIBMESH_ENABLE_AMR 
  913   parallel_object_only();
 
  916   LOG_SCOPE(
"distribute_dofs()", 
"DofMap");
 
  924   libmesh_assert_less (proc_id, n_proc);
 
  948   std::vector<dof_id_type> dofs_on_proc(n_proc, 0);
 
  949   this->
comm().allgather(next_free_dof, dofs_on_proc);
 
  952 #ifdef LIBMESH_ENABLE_AMR 
  978   libmesh_assert_equal_to (next_free_dof, 
_end_df[proc_id]);
 
 1007     for (
auto & node : 
mesh.node_ptr_range())
 
 1016               libmesh_assert_greater_equal (dofid, this->
first_dof(obj_proc_id));
 
 1017               libmesh_assert_less (dofid, this->
end_dof(obj_proc_id));
 
 1021     for (
auto & elem : 
mesh.element_ptr_range())
 
 1030               libmesh_assert_greater_equal (dofid, this->
first_dof(obj_proc_id));
 
 1031               libmesh_assert_less (dofid, this->
end_dof(obj_proc_id));
 
 1039 #ifdef LIBMESH_ENABLE_AMR 
 1063       gf->dofmap_reinit();
 
 1069       gf->dofmap_reinit();
 
 1086                                     unsigned int var_num)
 const 
 1092   const unsigned int sys_num       = this->
sys_number();
 
 1100       for (
auto & elem : 
mesh.active_local_element_ptr_range())
 
 1107           const unsigned int n_nodes = elem->n_nodes();
 
 1110           for (
unsigned int n=0; n<
n_nodes; n++)
 
 1112               Node & node = elem->node_ref(n);
 
 1117               const unsigned int n_comp = node.
n_comp(sys_num, var_num);
 
 1118               for(
unsigned int i=0; i<n_comp; i++)
 
 1123                   if (
idx.empty() || index > 
idx.back())
 
 1124                     idx.push_back(index);
 
 1129           const unsigned int n_comp = elem->n_comp(sys_num, var_num);
 
 1130           for (
unsigned int i=0; i<n_comp; i++)
 
 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);
 
 1147       for (
const auto & node : 
mesh.local_node_ptr_range())
 
 1151           const unsigned int n_comp = node->n_comp(sys_num, var_num);
 
 1152           for (
unsigned int i=0; i<n_comp; i++)
 
 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);
 
 1164       std::vector<dof_id_type> di_scalar;
 
 1166       idx.insert( 
idx.end(), di_scalar.begin(), di_scalar.end());
 
 1174   const unsigned int sys_num       = this->
sys_number();
 
 1182   for (
auto & elem : 
mesh.active_local_element_ptr_range())
 
 1186       const unsigned int n_nodes = elem->n_nodes();
 
 1189       for (
unsigned int n=0; n<
n_nodes; n++)
 
 1191           Node & node = elem->node_ref(n);
 
 1193           for (
unsigned vg=0; vg<n_var_groups; vg++)
 
 1218       for (
unsigned vg=0; vg<n_var_groups; vg++)
 
 1224             if (elem->n_comp_group(sys_num,vg) > 0)
 
 1226                 libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
 
 1229                 elem->set_vg_dof_base(sys_num,
 
 1234                                   elem->n_comp(sys_num,vg));
 
 1248   for (
auto & node : 
mesh.local_node_ptr_range())
 
 1249     for (
unsigned vg=0; vg<n_var_groups; vg++)
 
 1253         if (node->n_comp_group(sys_num,vg))
 
 1256               node->set_vg_dof_base (sys_num,
 
 1261                                 node->n_comp(sys_num,vg));
 
 1267   for (
unsigned vg=0; vg<n_var_groups; vg++)
 
 1290     MeshTools::libmesh_assert_valid_procids<Node>(
mesh);
 
 1292     for (
auto & node : 
mesh.local_node_ptr_range())
 
 1294         unsigned int n_var_g = node->n_var_groups(this->
sys_number());
 
 1295         for (
unsigned int vg=0; vg != n_var_g; ++vg)
 
 1297             unsigned int n_comp_g =
 
 1300               node->vg_dof_base(this->
sys_number(), vg) : 0;
 
 1313   const unsigned int sys_num      = this->
sys_number();
 
 1321   for (
unsigned vg=0; vg<n_var_groups; vg++)
 
 1325       const unsigned int n_vars_in_group = vg_description.
n_variables();
 
 1331       for (
auto & elem : 
mesh.active_local_element_ptr_range())
 
 1339           const unsigned int n_nodes = elem->n_nodes();
 
 1342           for (
unsigned int n=0; n<
n_nodes; n++)
 
 1344               Node & node = elem->node_ref(n);
 
 1355                   next_free_dof += (n_vars_in_group*
 
 1361           if (elem->n_comp_group(sys_num,vg) > 0)
 
 1363               libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
 
 1366               elem->set_vg_dof_base(sys_num,
 
 1370               next_free_dof += (n_vars_in_group*
 
 1371                                 elem->n_comp_group(sys_num,vg));
 
 1383       for (
auto & node : 
mesh.local_node_ptr_range())
 
 1384         if (node->n_comp_group(sys_num,vg))
 
 1387               node->set_vg_dof_base (sys_num,
 
 1391               next_free_dof += (n_vars_in_group*
 
 1392                                 node->n_comp_group(sys_num,vg));
 
 1398   for (
unsigned vg=0; vg<n_var_groups; vg++)
 
 1418     MeshTools::libmesh_assert_valid_procids<Node>(
mesh);
 
 1420     for (
auto & node : 
mesh.local_node_ptr_range())
 
 1422         unsigned int n_var_g = node->n_var_groups(this->
sys_number());
 
 1423         for (
unsigned int vg=0; vg != n_var_g; ++vg)
 
 1425             unsigned int n_comp_g =
 
 1428               node->vg_dof_base(this->
sys_number(), vg) : 0;
 
 1441                             std::set<CouplingMatrix *> & temporary_coupling_matrices,
 
 1442                             const std::set<GhostingFunctor *>::iterator & gf_begin,
 
 1443                             const std::set<GhostingFunctor *>::iterator & gf_end,
 
 1448   for (
const auto & gf : 
as_range(gf_begin, gf_end))
 
 1453       (*gf)(elems_begin, elems_end, p, more_elements_to_ghost);
 
 1455       for (
const auto & pr : more_elements_to_ghost)
 
 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);
 
 1463               if (existing_it->second)
 
 1470                       if (temporary_coupling_matrices.empty() ||
 
 1471                           temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second)) == temporary_coupling_matrices.end())
 
 1474                           temporary_coupling_matrices.insert(cm);
 
 1475                           existing_it->second = cm;
 
 1477                       const_cast<CouplingMatrix &>(*existing_it->second) &= *pr.second;
 
 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);
 
 1494                       existing_it->second = 
nullptr;
 
 1509   LOG_SCOPE(
"add_neighbors_to_send_list()", 
"DofMap");
 
 1518     = 
mesh.active_local_elements_begin();
 
 1520     = 
mesh.active_local_elements_end();
 
 1525   std::set<CouplingMatrix *> temporary_coupling_matrices;
 
 1531                                     temporary_coupling_matrices,
 
 1534                                     local_elem_it, local_elem_end, 
mesh.processor_id());
 
 1537                                     temporary_coupling_matrices,
 
 1540                                     local_elem_it, local_elem_end, 
mesh.processor_id());
 
 1545   std::map<const CouplingMatrix *, std::vector<unsigned int>>
 
 1546     column_variable_lists;
 
 1548   for (
auto & pr : elements_to_send)
 
 1550       const Elem * 
const partner = pr.first;
 
 1553       libmesh_assert_not_equal_to
 
 1563           libmesh_assert_equal_to (ghost_coupling->size(), n_var);
 
 1566           std::map<const CouplingMatrix *, std::vector<unsigned int>>::const_iterator
 
 1567             column_variable_list = column_variable_lists.find(ghost_coupling);
 
 1570           if (column_variable_list == column_variable_lists.end())
 
 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;
 
 1577               std::vector<unsigned int> & new_variable_list =
 
 1578                 inserted_variable_list_pair.first->second;
 
 1580               std::vector<unsigned char> has_variable(n_var, 
false);
 
 1582               for (
unsigned int vi = 0; vi != n_var; ++vi)
 
 1586                   for (
const auto & vj : ccr)
 
 1587                     has_variable[vj] = 
true;
 
 1589               for (
unsigned int vj = 0; vj != n_var; ++vj)
 
 1591                   if (has_variable[vj])
 
 1592                     new_variable_list.push_back(vj);
 
 1596           const std::vector<unsigned int> & variable_list =
 
 1597             column_variable_list->second;
 
 1599           for (
const auto & vj : variable_list)
 
 1601               std::vector<dof_id_type> di;
 
 1612           std::vector<dof_id_type> di;
 
 1616           for (
const auto & dof : di)
 
 1624   for (
auto & mat : temporary_coupling_matrices)
 
 1635   for ( ; local_elem_it != local_elem_end; ++local_elem_it)
 
 1637       const Elem * elem = *local_elem_it;
 
 1639       std::vector<dof_id_type> di;
 
 1643       for (
const auto & dof : di)
 
 1653   LOG_SCOPE(
"prepare_send_list()", 
"DofMap");
 
 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??" 
 1682   std::vector<dof_id_type>::iterator new_end =
 
 1695 #ifdef LIBMESH_ENABLE_CONSTRAINTS 
 1715   bool implicit_neighbor_dofs =
 
 1722   if (implicit_neighbor_dofs)
 
 1744       if (!implicit_neighbor_dofs)
 
 1752     bool all_discontinuous_dofs = 
true;
 
 1757         all_discontinuous_dofs = 
false;
 
 1759     if (all_discontinuous_dofs)
 
 1760       implicit_neighbor_dofs = 
true;
 
 1763   return implicit_neighbor_dofs;
 
 1781         mat->update_sparsity_pattern (
_sp->sparsity_pattern);
 
 1788         _n_nz = 
new std::vector<dof_id_type>();
 
 1791         _n_oz = 
new std::vector<dof_id_type>();
 
 1886                                    const std::vector<dof_id_type> & dof_indices_in,
 
 1889   const unsigned int n_original_dofs = dof_indices_in.size();
 
 1891 #ifdef LIBMESH_ENABLE_AMR 
 1894   libmesh_assert_equal_to (dof_indices_in.size(), Ue.
size());
 
 1895   bool has_constrained_dofs = 
false;
 
 1897   for (
unsigned int il=0; il != n_original_dofs; ++il)
 
 1903       libmesh_assert_less (ig, Ug.
size());
 
 1911   if (has_constrained_dofs)
 
 1914       std::vector<dof_id_type> constrained_dof_indices(dof_indices_in);
 
 1921       libmesh_assert_equal_to (dof_indices_in.size(), C.
m());
 
 1922       libmesh_assert_equal_to (constrained_dof_indices.size(), C.
n());
 
 1928       for (
unsigned int i=0; i != n_original_dofs; i++)
 
 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++)
 
 1936               const dof_id_type jg = constrained_dof_indices[j];
 
 1944               Ue.
el(i) += C(i,j)*Ug(jg);
 
 1953   libmesh_assert_equal_to (n_original_dofs, Ue.
size());
 
 1955   for (
unsigned int il=0; il<n_original_dofs; il++)
 
 1968                           std::vector<dof_id_type> & di)
 const 
 1978   LOG_SCOPE(
"dof_indices()", 
"DofMap");
 
 1987   std::size_t tot_size = 0;
 
 1993       const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
 
 1999           std::vector<const Node *> elem_nodes;
 
 2003           for (
unsigned int vg=0; vg<n_var_groups; vg++)
 
 2006               const unsigned int vars_in_group = var.
n_variables();
 
 2011                   for (
unsigned int vig=0; vig != vars_in_group; ++vig)
 
 2016                       std::vector<dof_id_type> di_new;
 
 2018                       di.insert( di.end(), di_new.begin(), di_new.end());
 
 2022                 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
 
 2026                                  cast_int<unsigned int>(elem_nodes.size())
 
 2028                                  , var.
number(vig), tot_size
 
 2040   for (
unsigned int vg=0; vg<n_var_groups; vg++)
 
 2043       const unsigned int vars_in_group = var.
n_variables();
 
 2049           for (
unsigned int vig=0; vig != vars_in_group; ++vig)
 
 2054               std::vector<dof_id_type> di_new;
 
 2056               di.insert( di.end(), di_new.begin(), di_new.end());
 
 2060         for (
unsigned int vig=0; vig != vars_in_group; ++vig)
 
 2065                          , var.
number(vig), tot_size
 
 2072   libmesh_assert_equal_to (tot_size, di.size());
 
 2078                           std::vector<dof_id_type> & di,
 
 2079                           const unsigned int vn,
 
 2085   LOG_SCOPE(
"dof_indices()", 
"DofMap");
 
 2091   if (p_level == -12345)
 
 2092     p_level = elem ? elem->
p_level() : 0;
 
 2096   const unsigned int vig = vn - var.
number();
 
 2100   std::size_t tot_size = 0;
 
 2106       const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
 
 2112           std::vector<const Node *> elem_nodes;
 
 2115           _dof_indices(*elem, p_level, di, vg, vig, elem_nodes.data(),
 
 2116                        cast_int<unsigned int>(elem_nodes.size())
 
 2134       std::vector<dof_id_type> di_new;
 
 2136       di.insert( di.end(), di_new.begin(), di_new.end());
 
 2147   libmesh_assert_equal_to (tot_size, di.size());
 
 2153                           std::vector<dof_id_type> & di)
 const 
 2158   LOG_SCOPE(
"dof_indices(Node)", 
"DofMap");
 
 2164   const unsigned int sys_num = this->
sys_number();
 
 2167   for (
unsigned int vg=0; vg<n_var_groups; vg++)
 
 2170       const unsigned int vars_in_group = var.
n_variables();
 
 2174           for (
unsigned int vig=0; vig != vars_in_group; ++vig)
 
 2176               std::vector<dof_id_type> di_new;
 
 2178               di.insert( di.end(), di_new.begin(), di_new.end());
 
 2184           for (
unsigned int vig=0; vig != vars_in_group; ++vig)
 
 2186               for (
int i=0; i != n_comp; ++i)
 
 2189                     node->
dof_number(sys_num, vg, vig, i, n_comp);
 
 2190                   libmesh_assert_not_equal_to
 
 2201                           std::vector<dof_id_type> & di,
 
 2202                           const unsigned int vn)
 const 
 2213   LOG_SCOPE(
"dof_indices(Node)", 
"DofMap");
 
 2218   const unsigned int sys_num = this->
sys_number();
 
 2226       std::vector<dof_id_type> di_new;
 
 2228       di.insert( di.end(), di_new.begin(), di_new.end());
 
 2232       const unsigned int vig = vn - var.
number();
 
 2234       for (
int i=0; i != n_comp; ++i)
 
 2237             node->
dof_number(sys_num, vg, vig, i, n_comp);
 
 2238           libmesh_assert_not_equal_to
 
 2248                           std::vector<dof_id_type> & di,
 
 2249                           const unsigned int vn)
 const 
 2256 #ifdef LIBMESH_ENABLE_AMR 
 2260                               std::vector<dof_id_type> & di,
 
 2261                               const unsigned int vn)
 const 
 2268 #endif // LIBMESH_ENABLE_AMR 
 2275                                 std::vector<dof_id_type> & di,
 
 2276                                 const unsigned int vn)
 const 
 2283   LOG_SCOPE(
"_node_dof_indices()", 
"DofMap");
 
 2286   const unsigned int dim = elem.
dim();
 
 2288   const unsigned int sys_num = this->
sys_number();
 
 2289   const std::pair<unsigned int, unsigned int>
 
 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);
 
 2297   fe_type.
order = static_cast<Order>(fe_type.
order +
 
 2299   const bool extra_hanging_dofs =
 
 2306   const unsigned int nc =
 
 2313   if (extra_hanging_dofs && nc && !elem.
is_vertex(n))
 
 2315       const int dof_offset = n_comp - nc;
 
 2326         for (
unsigned int i = dof_offset; i != n_comp; ++i)
 
 2338     for (
unsigned int i=0; i<nc; i++)
 
 2351                            std::vector<dof_id_type> & di,
 
 2352                            const unsigned int vg,
 
 2353                            const unsigned int vig,
 
 2354                            const Node * 
const * nodes,
 
 2358                            const unsigned int v,
 
 2359                            std::size_t & tot_size
 
 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();
 
 2376       fe_type.
order = static_cast<Order>(fe_type.
order + p_level);
 
 2378       const bool extra_hanging_dofs =
 
 2393       for (
unsigned int n=0; n != 
n_nodes; n++)
 
 2395           const Node & node = *nodes[n];
 
 2400           const std::pair<unsigned int, unsigned int>
 
 2402           libmesh_assert_equal_to (vg, vg_and_offset.first);
 
 2403           libmesh_assert_equal_to (vig, vg_and_offset.second);
 
 2405           const unsigned int n_comp = node.
n_comp_group(sys_num,vg);
 
 2411           const unsigned int nc =
 
 2412 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
 2416             ndan (type, fe_type.
order, n);
 
 2422           if (extra_hanging_dofs && !elem.
is_vertex(n))
 
 2424               const int dof_offset = n_comp - nc;
 
 2435                 for (
int i=n_comp-1; i>=dof_offset; i--)
 
 2438                       node.
dof_number(sys_num, vg, vig, i, n_comp);
 
 2447             for (
unsigned int i=0; i<nc; i++)
 
 2450                   node.
dof_number(sys_num, vg, vig, i, n_comp);
 
 2465           const unsigned int n_comp = elem.
n_comp_group(sys_num,vg);
 
 2466           if (elem.
n_systems() > sys_num && nc <= n_comp)
 
 2468               for (
unsigned int i=0; i<nc; i++)
 
 2471                     elem.
dof_number(sys_num, vg, vig, i, n_comp);
 
 2489                                  const unsigned int vn,
 
 2490 #ifdef LIBMESH_ENABLE_AMR
 
 2497   LOG_SCOPE(
"SCALAR_dof_indices()", 
"DofMap");
 
 2501 #ifdef LIBMESH_ENABLE_AMR 
 2517   di.resize(n_dofs_vn);
 
 2518   for (
int i = 0; i != n_dofs_vn; ++i)
 
 2543   for (
const auto & di : dof_indices_in)
 
 2552 template <
typename DofObjectSub
class>
 
 2554                           unsigned int var_num)
 const 
 2560   std::vector<dof_id_type> di;
 
 2572 #ifdef LIBMESH_ENABLE_AMR 
 2575                               std::vector<dof_id_type> & di,
 
 2576                               const unsigned int vn)
 const 
 2578   LOG_SCOPE(
"old_dof_indices()", 
"DofMap");
 
 2583   const unsigned int sys_num       = this->
sys_number();
 
 2585   const unsigned int dim           = elem->
dim();
 
 2586 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
 2587   const bool is_inf                = elem->
infinite();
 
 2601   std::vector<const Node *> elem_nodes;
 
 2602   const Node * 
const * nodes_ptr;
 
 2607       const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
 
 2609       nodes_ptr = elem_nodes.data();
 
 2610       n_nodes = cast_int<unsigned int>(elem_nodes.size());
 
 2620   for (
unsigned int vg=0; vg<n_var_groups; vg++)
 
 2623       const unsigned int vars_in_group = var.
n_variables();
 
 2625       for (
unsigned int vig=0; vig<vars_in_group; vig++)
 
 2627           const unsigned int v = var.
number(vig);
 
 2635                   std::vector<dof_id_type> di_new;
 
 2637                   di.insert( di.end(), di_new.begin(), di_new.end());
 
 2647                     int p_adjustment = 0;
 
 2650                         libmesh_assert_greater (elem->
p_level(), 0);
 
 2658                     fe_type.
order = static_cast<Order>(fe_type.
order +
 
 2662                     const bool extra_hanging_dofs =
 
 2669                     for (
unsigned int n=0; n<
n_nodes; n++)
 
 2671                         const Node * node = nodes_ptr[n];
 
 2679                         const unsigned int nc =
 
 2680 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
 2684                           ndan (type, fe_type.
order, n);
 
 2686                         const int n_comp = old_dof_obj->
n_comp_group(sys_num,vg);
 
 2692                         if (extra_hanging_dofs && !elem->
is_vertex(n))
 
 2694                             const int dof_offset = n_comp - nc;
 
 2706                               for (
int i=n_comp-1; i>=dof_offset; i--)
 
 2709                                     old_dof_obj->
dof_number(sys_num, vg, vig, i, n_comp);
 
 2718                           for (
unsigned int i=0; i<nc; i++)
 
 2721                                 old_dof_obj->
dof_number(sys_num, vg, vig, i, n_comp);
 
 2740                         const unsigned int n_comp =
 
 2743                         if (old_dof_obj->
n_systems() > sys_num &&
 
 2747                             for (
unsigned int i=0; i<nc; i++)
 
 2750                                   old_dof_obj->
dof_number(sys_num, vg, vig, i, n_comp);
 
 2768 #endif // LIBMESH_ENABLE_AMR 
 2771 #ifdef LIBMESH_ENABLE_CONSTRAINTS 
 2775   typedef std::set<dof_id_type> RCSet;
 
 2778   RCSet dof_set (elem_dofs.begin(), elem_dofs.end());
 
 2787   for (
const auto & dof : elem_dofs)
 
 2791         DofConstraints::const_iterator
 
 2805         for (
const auto & pr : constraint_row)
 
 2806           if (!dof_set.count (pr.first))
 
 2808               dof_set.insert (pr.first);
 
 2820       elem_dofs.insert (elem_dofs.end(),
 
 2821                         dof_set.begin(), dof_set.end());
 
 2832 #endif // LIBMESH_ENABLE_CONSTRAINTS 
 2845   std::ostringstream os;
 
 2850   const char * may_equal = 
" <= ";
 
 2855     if (mat->need_full_sparsity_pattern())
 
 2859   long double avg_n_nz = 0, avg_n_oz = 0;
 
 2863       for (
const auto & val : *
_n_nz)
 
 2865           max_n_nz = std::max(max_n_nz, val);
 
 2869       std::size_t n_nz_size = 
_n_nz->size();
 
 2871       this->
comm().max(max_n_nz);
 
 2872       this->
comm().sum(avg_n_nz);
 
 2873       this->
comm().sum(n_nz_size);
 
 2875       avg_n_nz /= std::max(n_nz_size,std::size_t(1));
 
 2879       for (
const auto & val : *
_n_oz)
 
 2881           max_n_oz = std::max(max_n_oz, val);
 
 2885       std::size_t n_oz_size = 
_n_oz->size();
 
 2887       this->
comm().max(max_n_oz);
 
 2888       this->
comm().sum(avg_n_oz);
 
 2889       this->
comm().sum(n_oz_size);
 
 2891       avg_n_oz /= std::max(n_oz_size,std::size_t(1));
 
 2894   os << 
"    DofMap Sparsity\n      Average  On-Processor Bandwidth" 
 2895      << may_equal << avg_n_nz << 
'\n';
 
 2897   os << 
"      Average Off-Processor Bandwidth" 
 2898      << may_equal << avg_n_oz << 
'\n';
 
 2900   os << 
"      Maximum  On-Processor Bandwidth" 
 2901      << may_equal << max_n_nz << 
'\n';
 
 2903   os << 
"      Maximum Off-Processor Bandwidth" 
 2904      << may_equal << max_n_oz << std::endl;
 
 2906 #ifdef LIBMESH_ENABLE_CONSTRAINTS 
 2908   std::size_t n_constraints = 0, max_constraint_length = 0,
 
 2910   long double avg_constraint_length = 0.;
 
 2920       std::size_t rowsize = row.size();
 
 2922       max_constraint_length = std::max(max_constraint_length,
 
 2924       avg_constraint_length += rowsize;
 
 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);
 
 2936   os << 
"    DofMap Constraints\n      Number of DoF Constraints = " 
 2940        << 
"      Number of Heterogenous Constraints= " << n_rhss;
 
 2943       avg_constraint_length /= n_constraints;
 
 2946          << 
"      Average DoF Constraint Length= " << avg_constraint_length;
 
 2949 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 
 2950   std::size_t n_node_constraints = 0, max_node_constraint_length = 0,
 
 2952   long double avg_node_constraint_length = 0.;
 
 2957       const Node * node = pr.first;
 
 2962       std::size_t rowsize = row.
size();
 
 2964       max_node_constraint_length = std::max(max_node_constraint_length,
 
 2966       avg_node_constraint_length += rowsize;
 
 2967       n_node_constraints++;
 
 2969       if (pr.second.second != 
Point(0))
 
 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);
 
 2978   os << 
"\n      Number of Node Constraints = " << n_node_constraints;
 
 2981        << 
"      Number of Heterogenous Node Constraints= " << n_node_rhss;
 
 2982   if (n_node_constraints)
 
 2984       avg_node_constraint_length /= n_node_constraints;
 
 2985       os << 
"\n      Maximum Node Constraint Length= " << max_node_constraint_length
 
 2987          << 
"      Average Node Constraint Length= " << avg_node_constraint_length;
 
 2989 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 
 2993 #endif // LIBMESH_ENABLE_CONSTRAINTS 
 2999 template bool DofMap::is_evaluable<Elem>(
const Elem &, 
unsigned int) 
const;
 
 3000 template bool DofMap::is_evaluable<Node>(
const Node &, 
unsigned int) 
const;