50 #include <unordered_map> 57 std::unique_ptr<SparsityPattern::Build>
60 libmesh_assert (
mesh.is_prepared());
62 LOG_SCOPE(
"build_sparsity()",
"DofMap");
85 std::unique_ptr<SparsityPattern::Build> sp
90 implicit_neighbor_dofs,
94 mesh.active_local_elements_end()), *sp);
103 libmesh_assert_equal_to (sp->sparsity_pattern.size(), n_dofs_on_proc);
111 libMesh::out <<
"WARNING: You have specified both an extra sparsity function and object.\n" 112 <<
" Are you sure this is what you meant to do??" 117 (sp->sparsity_pattern, sp->n_nz,
123 (sp->sparsity_pattern, sp->n_nz, sp->n_oz);
125 return std::unique_ptr<SparsityPattern::Build>(sp.release());
133 _dof_coupling(nullptr),
134 _error_on_cyclic_constraint(false),
137 _variable_group_numbers(),
145 _augment_sparsity_pattern(nullptr),
146 _extra_sparsity_function(nullptr),
147 _extra_sparsity_context(nullptr),
148 _augment_send_list(nullptr),
149 _extra_send_list_function(nullptr),
150 _extra_send_list_context(nullptr),
153 need_full_sparsity_pattern(false),
158 #ifdef LIBMESH_ENABLE_AMR
162 _first_old_scalar_df()
164 #ifdef LIBMESH_ENABLE_CONSTRAINTS
166 , _stashed_dof_constraints()
167 , _primal_constraint_values()
168 , _adjoint_constraint_values()
170 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
171 , _node_constraints()
173 #ifdef LIBMESH_ENABLE_PERIODIC
176 #ifdef LIBMESH_ENABLE_DIRICHLET
178 , _adjoint_dirichlet_boundaries()
180 , _implicit_neighbor_dofs_initialized(false),
181 _implicit_neighbor_dofs(false)
189 #ifdef LIBMESH_ENABLE_PERIODIC 210 #ifdef LIBMESH_ENABLE_DIRICHLET 217 #ifdef LIBMESH_ENABLE_PERIODIC 253 for (
unsigned int var=0; var<new_var_group.
n_variables(); var++)
264 parallel_object_only();
277 bool computed_sparsity_already =
280 this->
comm().
max(computed_sparsity_already);
281 if (computed_sparsity_already &&
287 libmesh_assert(
_sp.get());
308 return mesh.node_ptr(i);
315 return mesh.elem_ptr(i);
320 template <
typename iterator_type>
322 iterator_type objects_end,
324 dofobject_accessor objects)
327 parallel_object_only();
331 std::unordered_map<processor_id_type, dof_id_type> ghost_objects_from_proc;
333 iterator_type it = objects_begin;
335 for (; it != objects_end; ++it)
344 ghost_objects_from_proc[obj_procid]++;
349 std::map<processor_id_type, std::vector<dof_id_type>>
354 auto ghost_end = ghost_objects_from_proc.end();
358 auto ghost_it = ghost_objects_from_proc.find(p);
359 if (ghost_it != ghost_end)
360 requested_ids[p].reserve(ghost_it->second);
363 for (it = objects_begin; it != objects_end; ++it)
372 if (ghost_objects_from_proc.count(p))
373 libmesh_assert_equal_to (requested_ids[p].size(), ghost_objects_from_proc[p]);
375 libmesh_assert(!requested_ids.count(p));
379 typedef std::vector<dof_id_type> datum;
381 auto gather_functor =
382 [
this, &
mesh, &objects]
384 const std::vector<dof_id_type> & ids,
385 std::vector<datum> &
data)
392 const std::size_t query_size = ids.size();
394 data.resize(query_size);
395 for (
auto & d :
data)
396 d.resize(2 * n_var_groups);
398 for (std::size_t i=0; i != query_size; ++i)
401 libmesh_assert(requested);
403 libmesh_assert_equal_to (requested->
n_var_groups(sys_num), n_var_groups);
404 for (
unsigned int vg=0; vg != n_var_groups; ++vg)
406 unsigned int n_comp_g =
408 data[i][vg] = n_comp_g;
412 data[i][n_var_groups+vg] = my_first_dof;
417 auto action_functor =
418 [
this, &
mesh, &objects]
420 const std::vector<dof_id_type> & ids,
421 const std::vector<datum> &
data)
431 libmesh_assert(requested);
432 libmesh_assert_equal_to (requested->
processor_id(), pid);
433 for (
unsigned int vg=0; vg != n_var_groups; ++vg)
435 unsigned int n_comp_g =
436 cast_int<unsigned int>(
data[i][vg]);
443 (sys_num, vg, my_first_dof);
449 datum * ex =
nullptr;
451 (this->
comm(), requested_ids, gather_functor, action_functor, ex);
455 for (it = objects_begin; it != objects_end; ++it)
458 libmesh_assert (obj);
460 for (
unsigned int v=0; v != num_variables; ++v)
462 unsigned int n_comp =
476 libmesh_assert (
mesh.is_prepared());
478 LOG_SCOPE(
"reinit()",
"DofMap");
493 unsigned int standard_n_levels =
509 std::vector<unsigned int> n_vars_per_group; n_vars_per_group.reserve (n_var_groups);
511 for (
unsigned int vg=0; vg<n_var_groups; vg++)
512 n_vars_per_group.push_back (this->variable_group(vg).n_variables());
514 #ifdef LIBMESH_ENABLE_AMR 519 for (
auto & node :
mesh.node_ptr_range())
521 node->clear_old_dof_object();
522 libmesh_assert (!node->old_dof_object);
525 for (
auto & elem :
mesh.element_ptr_range())
527 elem->clear_old_dof_object();
528 libmesh_assert (!elem->old_dof_object);
536 for (
auto & elem :
mesh.element_ptr_range())
542 for (
unsigned int n=0; n<elem->n_nodes(); n++)
544 Node & node = elem->node_ref(n);
551 libmesh_assert (!elem->old_dof_object);
553 if (elem->has_dofs(sys_num))
554 elem->set_old_dof_object();
557 #endif // #ifdef LIBMESH_ENABLE_AMR 566 for (
auto & node :
mesh.node_ptr_range())
567 node->set_n_vars_per_group(sys_num, n_vars_per_group);
570 for (
auto & elem :
mesh.element_ptr_range())
571 elem->set_n_vars_per_group(sys_num, n_vars_per_group);
578 for (
unsigned int vg=0; vg<n_var_groups; vg++)
582 const unsigned int n_var_in_group = vg_description.
n_variables();
583 const FEType & base_fe_type = vg_description.
type();
594 const bool extra_hanging_dofs =
598 for (
auto & elem :
mesh.active_element_ptr_range())
600 libmesh_assert(elem);
608 const unsigned int dim = elem->dim();
610 FEType fe_type = base_fe_type;
612 #ifdef LIBMESH_ENABLE_AMR 615 if (elem->p_level() + base_fe_type.
order >
618 libmesh_assert_less_msg(static_cast<unsigned int>(base_fe_type.
order.
get_order()),
620 "ERROR: Finite element " 622 <<
" on geometric element " 624 <<
"\nonly supports FEInterface::max_order = " 626 <<
", not fe_type.order = " 627 << base_fe_type.
order);
632 <<
" on geometric element " 634 <<
"could not be p refined past FEInterface::max_order = " 639 - base_fe_type.
order);
647 for (
unsigned int n=0; n<elem->n_nodes(); n++)
649 Node & node = elem->node_ref(n);
651 if (elem->is_vertex(n))
653 const unsigned int old_node_dofs =
656 const unsigned int vertex_dofs =
662 if (vertex_dofs > old_node_dofs)
687 for (
auto & elem :
mesh.active_element_ptr_range())
689 libmesh_assert(elem);
697 const unsigned int dim = elem->dim();
699 FEType fe_type = base_fe_type;
704 for (
unsigned int n=0; n<elem->n_nodes(); n++)
706 Node & node = elem->node_ref(n);
708 const unsigned int old_node_dofs =
711 const unsigned int vertex_dofs = old_node_dofs?
712 cast_int<unsigned int>(node.
vg_dof_base (sys_num,vg)):0;
714 const unsigned int new_node_dofs =
718 if (elem->is_vertex(n))
720 libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
730 libmesh_assert_greater_equal (vertex_dofs, new_node_dofs);
750 else if (vertex_dofs == 0)
752 if (new_node_dofs > old_node_dofs)
764 else if (extra_hanging_dofs)
766 if (new_node_dofs > old_node_dofs - vertex_dofs)
769 vertex_dofs + new_node_dofs);
779 libmesh_assert_greater_equal (old_node_dofs, vertex_dofs);
780 if (new_node_dofs > old_node_dofs)
792 const unsigned int dofs_per_elem =
796 elem->set_n_comp_group(sys_num, vg, dofs_per_elem);
815 const unsigned int sys_num = this->
sys_number();
818 for (
auto & node :
mesh.node_ptr_range())
819 node->invalidate_dofs(sys_num);
822 for (
auto & elem :
mesh.active_element_ptr_range())
823 elem->invalidate_dofs(sys_num);
846 this->_coupling_functors.clear();
863 this->_algebraic_ghosting_functors.clear();
882 #ifdef LIBMESH_ENABLE_AMR 905 parallel_object_only();
908 LOG_SCOPE(
"distribute_dofs()",
"DofMap");
910 libmesh_assert (
mesh.is_prepared());
916 libmesh_assert_less (proc_id, n_proc);
940 std::vector<dof_id_type> dofs_on_proc(n_proc, 0);
944 #ifdef LIBMESH_ENABLE_AMR 970 libmesh_assert_equal_to (next_free_dof,
_end_df[proc_id]);
999 for (
auto & node :
mesh.node_ptr_range())
1004 for (
unsigned int v=0; v != dofobj->
n_vars(sys_num); ++v)
1005 for (
unsigned int c=0; c != dofobj->
n_comp(sys_num,v); ++c)
1008 libmesh_assert_greater_equal (dofid, this->
first_dof(obj_proc_id));
1009 libmesh_assert_less (dofid, this->
end_dof(obj_proc_id));
1013 for (
auto & elem :
mesh.element_ptr_range())
1018 for (
unsigned int v=0; v != dofobj->
n_vars(sys_num); ++v)
1019 for (
unsigned int c=0; c != dofobj->
n_comp(sys_num,v); ++c)
1022 libmesh_assert_greater_equal (dofid, this->
first_dof(obj_proc_id));
1023 libmesh_assert_less (dofid, this->
end_dof(obj_proc_id));
1031 #ifdef LIBMESH_ENABLE_AMR 1044 for (
unsigned int v=0; v<this->
n_variables(); v++)
1055 gf->dofmap_reinit();
1061 gf->dofmap_reinit();
1078 unsigned int var_num)
const 1084 const unsigned int sys_num = this->
sys_number();
1092 for (
auto & elem :
mesh.active_local_element_ptr_range())
1099 const unsigned int n_nodes = elem->n_nodes();
1102 for (
unsigned int n=0; n<
n_nodes; n++)
1104 Node & node = elem->node_ref(n);
1109 const unsigned int n_comp = node.
n_comp(sys_num, var_num);
1110 for(
unsigned int i=0; i<n_comp; i++)
1115 if (
idx.empty() || index >
idx.back())
1116 idx.push_back(index);
1121 const unsigned int n_comp = elem->n_comp(sys_num, var_num);
1122 for (
unsigned int i=0; i<n_comp; i++)
1124 const dof_id_type index = elem->dof_number(sys_num,var_num,i);
1125 if (
idx.empty() || index >
idx.back())
1126 idx.push_back(index);
1139 for (
const auto & node :
mesh.local_node_ptr_range())
1141 libmesh_assert(node);
1143 const unsigned int n_comp = node->n_comp(sys_num, var_num);
1144 for (
unsigned int i=0; i<n_comp; i++)
1146 const dof_id_type index = node->dof_number(sys_num,var_num,i);
1147 if (
idx.empty() || index >
idx.back())
1148 idx.push_back(index);
1156 std::vector<dof_id_type> di_scalar;
1158 idx.insert(
idx.end(), di_scalar.begin(), di_scalar.end());
1166 const unsigned int sys_num = this->
sys_number();
1174 for (
auto & elem :
mesh.active_local_element_ptr_range())
1178 const unsigned int n_nodes = elem->n_nodes();
1181 for (
unsigned int n=0; n<
n_nodes; n++)
1183 Node & node = elem->node_ref(n);
1185 for (
unsigned vg=0; vg<n_var_groups; vg++)
1210 for (
unsigned vg=0; vg<n_var_groups; vg++)
1216 if (elem->n_comp_group(sys_num,vg) > 0)
1218 libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
1221 elem->set_vg_dof_base(sys_num,
1226 elem->n_comp(sys_num,vg));
1240 for (
auto & node :
mesh.local_node_ptr_range())
1241 for (
unsigned vg=0; vg<n_var_groups; vg++)
1245 if (node->n_comp_group(sys_num,vg))
1248 node->set_vg_dof_base (sys_num,
1253 node->n_comp(sys_num,vg));
1259 for (
unsigned vg=0; vg<n_var_groups; vg++)
1282 MeshTools::libmesh_assert_valid_procids<Node>(
mesh);
1284 for (
auto & node :
mesh.local_node_ptr_range())
1286 unsigned int n_var_g = node->n_var_groups(this->
sys_number());
1287 for (
unsigned int vg=0; vg != n_var_g; ++vg)
1289 unsigned int n_comp_g =
1292 node->vg_dof_base(this->
sys_number(), vg) : 0;
1305 const unsigned int sys_num = this->
sys_number();
1313 for (
unsigned vg=0; vg<n_var_groups; vg++)
1317 const unsigned int n_vars_in_group = vg_description.
n_variables();
1323 for (
auto & elem :
mesh.active_local_element_ptr_range())
1331 const unsigned int n_nodes = elem->n_nodes();
1334 for (
unsigned int n=0; n<
n_nodes; n++)
1336 Node & node = elem->node_ref(n);
1347 next_free_dof += (n_vars_in_group*
1353 if (elem->n_comp_group(sys_num,vg) > 0)
1355 libmesh_assert_equal_to (elem->vg_dof_base(sys_num,vg),
1358 elem->set_vg_dof_base(sys_num,
1362 next_free_dof += (n_vars_in_group*
1363 elem->n_comp_group(sys_num,vg));
1375 for (
auto & node :
mesh.local_node_ptr_range())
1376 if (node->n_comp_group(sys_num,vg))
1379 node->set_vg_dof_base (sys_num,
1383 next_free_dof += (n_vars_in_group*
1384 node->n_comp_group(sys_num,vg));
1390 for (
unsigned vg=0; vg<n_var_groups; vg++)
1410 MeshTools::libmesh_assert_valid_procids<Node>(
mesh);
1412 for (
auto & node :
mesh.local_node_ptr_range())
1414 unsigned int n_var_g = node->n_var_groups(this->
sys_number());
1415 for (
unsigned int vg=0; vg != n_var_g; ++vg)
1417 unsigned int n_comp_g =
1420 node->vg_dof_base(this->
sys_number(), vg) : 0;
1433 std::set<CouplingMatrix *> & temporary_coupling_matrices,
1434 const std::set<GhostingFunctor *>::iterator & gf_begin,
1435 const std::set<GhostingFunctor *>::iterator & gf_end,
1440 for (
const auto & gf :
as_range(gf_begin, gf_end))
1445 (*gf)(elems_begin, elems_end, p, more_elements_to_ghost);
1447 for (
const auto & pr : more_elements_to_ghost)
1449 GhostingFunctor::map_type::iterator existing_it =
1450 elements_to_ghost.find (pr.first);
1451 if (existing_it == elements_to_ghost.end())
1452 elements_to_ghost.insert(pr);
1455 if (existing_it->second)
1462 if (temporary_coupling_matrices.empty() ||
1463 temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second)) == temporary_coupling_matrices.end())
1466 temporary_coupling_matrices.insert(cm);
1467 existing_it->second = cm;
1469 const_cast<CouplingMatrix &
>(*existing_it->second) &= *pr.second;
1481 std::set<CouplingMatrix *>::iterator temp_it =
1482 temporary_coupling_matrices.find(const_cast<CouplingMatrix *>(existing_it->second));
1483 if (temp_it != temporary_coupling_matrices.end())
1484 temporary_coupling_matrices.erase(temp_it);
1486 existing_it->second =
nullptr;
1501 LOG_SCOPE(
"add_neighbors_to_send_list()",
"DofMap");
1506 =
mesh.active_local_elements_begin();
1508 =
mesh.active_local_elements_end();
1513 std::set<CouplingMatrix *> temporary_coupling_matrices;
1519 temporary_coupling_matrices,
1522 local_elem_it, local_elem_end,
mesh.processor_id());
1525 temporary_coupling_matrices,
1528 local_elem_it, local_elem_end,
mesh.processor_id());
1533 std::map<const CouplingMatrix *, std::vector<unsigned int>>
1534 column_variable_lists;
1536 for (
auto & pr : elements_to_send)
1538 const Elem *
const partner = pr.first;
1541 libmesh_assert_not_equal_to
1551 libmesh_assert_equal_to (ghost_coupling->size(), n_var);
1554 std::map<const CouplingMatrix *, std::vector<unsigned int>>::const_iterator
1555 column_variable_list = column_variable_lists.find(ghost_coupling);
1558 if (column_variable_list == column_variable_lists.end())
1560 std::pair<std::map<const CouplingMatrix *, std::vector<unsigned int>>::iterator,
bool>
1561 inserted_variable_list_pair = column_variable_lists.insert(std::make_pair(ghost_coupling,
1562 std::vector<unsigned int>()));
1563 column_variable_list = inserted_variable_list_pair.first;
1565 std::vector<unsigned int> & new_variable_list =
1566 inserted_variable_list_pair.first->second;
1568 std::vector<unsigned char> has_variable(n_var,
false);
1570 for (
unsigned int vi = 0; vi != n_var; ++vi)
1574 for (
const auto & vj : ccr)
1575 has_variable[vj] =
true;
1577 for (
unsigned int vj = 0; vj != n_var; ++vj)
1579 if (has_variable[vj])
1580 new_variable_list.push_back(vj);
1584 const std::vector<unsigned int> & variable_list =
1585 column_variable_list->second;
1587 for (
const auto & vj : variable_list)
1589 std::vector<dof_id_type> di;
1600 std::vector<dof_id_type> di;
1604 for (
const auto & dof : di)
1612 for (
auto & mat : temporary_coupling_matrices)
1623 for ( ; local_elem_it != local_elem_end; ++local_elem_it)
1625 const Elem * elem = *local_elem_it;
1627 std::vector<dof_id_type> di;
1631 for (
const auto & dof : di)
1641 LOG_SCOPE(
"prepare_send_list()",
"DofMap");
1649 libMesh::out <<
"WARNING: You have specified both an extra send list function and object.\n" 1650 <<
" Are you sure this is what you meant to do??" 1666 std::vector<dof_id_type>::iterator new_end =
1686 bool implicit_neighbor_dofs =
1693 if (implicit_neighbor_dofs)
1715 if (!implicit_neighbor_dofs)
1723 bool all_discontinuous_dofs =
true;
1725 for (
unsigned int var=0; var<this->
n_variables(); var++)
1728 all_discontinuous_dofs =
false;
1730 if (all_discontinuous_dofs)
1731 implicit_neighbor_dofs =
true;
1734 return implicit_neighbor_dofs;
1752 mat->update_sparsity_pattern (
_sp->sparsity_pattern);
1759 _n_nz =
new std::vector<dof_id_type>();
1762 _n_oz =
new std::vector<dof_id_type>();
1775 libmesh_assert(
_sp.get());
1782 libmesh_assert(!
_sp.get());
1849 const std::vector<dof_id_type> & dof_indices_in,
1852 const unsigned int n_original_dofs = dof_indices_in.size();
1854 #ifdef LIBMESH_ENABLE_AMR 1857 libmesh_assert_equal_to (dof_indices_in.size(), Ue.
size());
1858 bool has_constrained_dofs =
false;
1860 for (
unsigned int il=0; il != n_original_dofs; ++il)
1866 libmesh_assert_less (ig, Ug.
size());
1874 if (has_constrained_dofs)
1877 std::vector<dof_id_type> constrained_dof_indices(dof_indices_in);
1884 libmesh_assert_equal_to (dof_indices_in.size(), C.
m());
1885 libmesh_assert_equal_to (constrained_dof_indices.size(), C.
n());
1891 for (
unsigned int i=0; i != n_original_dofs; i++)
1895 const unsigned int n_constrained =
1896 cast_int<unsigned int>(constrained_dof_indices.size());
1897 for (
unsigned int j=0; j<n_constrained; j++)
1899 const dof_id_type jg = constrained_dof_indices[j];
1907 Ue.
el(i) += C(i,j)*Ug(jg);
1916 libmesh_assert_equal_to (n_original_dofs, Ue.
size());
1918 for (
unsigned int il=0; il<n_original_dofs; il++)
1931 std::vector<dof_id_type> & di)
const 1939 libmesh_assert(!elem || elem->
active());
1941 LOG_SCOPE(
"dof_indices()",
"DofMap");
1950 std::size_t tot_size = 0;
1962 std::vector<const Node *> elem_nodes;
1966 for (
unsigned int vg=0; vg<n_var_groups; vg++)
1969 const unsigned int vars_in_group = var.
n_variables();
1974 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
1979 std::vector<dof_id_type> di_new;
1981 di.insert( di.end(), di_new.begin(), di_new.end());
1985 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
1989 cast_int<unsigned int>(elem_nodes.size())
1991 , var.
number(vig), tot_size
2003 for (
unsigned int vg=0; vg<n_var_groups; vg++)
2006 const unsigned int vars_in_group = var.
n_variables();
2012 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
2017 std::vector<dof_id_type> di_new;
2019 di.insert( di.end(), di_new.begin(), di_new.end());
2023 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
2028 , var.
number(vig), tot_size
2035 libmesh_assert_equal_to (tot_size, di.size());
2041 std::vector<dof_id_type> & di,
2042 const unsigned int vn,
2048 LOG_SCOPE(
"dof_indices()",
"DofMap");
2054 if (p_level == -12345)
2055 p_level = elem ? elem->
p_level() : 0;
2059 const unsigned int vig = vn - var.
number();
2063 std::size_t tot_size = 0;
2075 std::vector<const Node *> elem_nodes;
2078 _dof_indices(*elem, p_level, di, vg, vig, elem_nodes.data(),
2079 cast_int<unsigned int>(elem_nodes.size())
2097 std::vector<dof_id_type> di_new;
2099 di.insert( di.end(), di_new.begin(), di_new.end());
2110 libmesh_assert_equal_to (tot_size, di.size());
2116 std::vector<dof_id_type> & di)
const 2121 LOG_SCOPE(
"dof_indices(Node)",
"DofMap");
2127 const unsigned int sys_num = this->
sys_number();
2130 for (
unsigned int vg=0; vg<n_var_groups; vg++)
2133 const unsigned int vars_in_group = var.
n_variables();
2137 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
2139 std::vector<dof_id_type> di_new;
2141 di.insert( di.end(), di_new.begin(), di_new.end());
2147 for (
unsigned int vig=0; vig != vars_in_group; ++vig)
2149 for (
int i=0; i != n_comp; ++i)
2152 node->
dof_number(sys_num, vg, vig, i, n_comp);
2153 libmesh_assert_not_equal_to
2164 std::vector<dof_id_type> & di,
2165 const unsigned int vn)
const 2176 LOG_SCOPE(
"dof_indices(Node)",
"DofMap");
2181 const unsigned int sys_num = this->
sys_number();
2189 std::vector<dof_id_type> di_new;
2191 di.insert( di.end(), di_new.begin(), di_new.end());
2195 const unsigned int vig = vn - var.
number();
2197 for (
int i=0; i != n_comp; ++i)
2200 node->
dof_number(sys_num, vg, vig, i, n_comp);
2201 libmesh_assert_not_equal_to
2211 std::vector<dof_id_type> & di,
2212 const unsigned int vg,
2213 const unsigned int vig,
2214 const Node *
const * nodes,
2218 const unsigned int v,
2219 std::size_t & tot_size
2228 const unsigned int sys_num = this->
sys_number();
2229 const unsigned int dim = elem.
dim();
2230 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 2231 const bool is_inf = elem.
infinite();
2238 const bool extra_hanging_dofs =
2253 for (
unsigned int n=0; n !=
n_nodes; n++)
2255 const Node * node = nodes[n];
2260 const std::pair<unsigned int, unsigned int>
2262 libmesh_assert_equal_to (vg, vg_and_offset.first);
2263 libmesh_assert_equal_to (vig, vg_and_offset.second);
2265 const unsigned int n_comp = node->
n_comp_group(sys_num,vg);
2271 const unsigned int nc =
2272 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 2276 ndan (type, fe_type.
order, n);
2282 if (extra_hanging_dofs && !elem.
is_vertex(n))
2284 const int dof_offset = n_comp - nc;
2291 libmesh_assert(!elem.
active());
2295 for (
int i=n_comp-1; i>=dof_offset; i--)
2298 node->
dof_number(sys_num, vg, vig, i, n_comp);
2307 for (
unsigned int i=0; i<nc; i++)
2310 node->
dof_number(sys_num, vg, vig, i, n_comp);
2325 const unsigned int n_comp = elem.
n_comp_group(sys_num,vg);
2326 if (elem.
n_systems() > sys_num && nc <= n_comp)
2328 for (
unsigned int i=0; i<nc; i++)
2331 elem.
dof_number(sys_num, vg, vig, i, n_comp);
2349 const unsigned int vn,
2350 #ifdef LIBMESH_ENABLE_AMR
2357 LOG_SCOPE(
"SCALAR_dof_indices()",
"DofMap");
2361 #ifdef LIBMESH_ENABLE_AMR 2377 di.resize(n_dofs_vn);
2378 for (
int i = 0; i != n_dofs_vn; ++i)
2403 for (
const auto & di : dof_indices_in)
2412 template <
typename DofObjectSub
class>
2414 unsigned int var_num)
const 2420 std::vector<dof_id_type> di;
2432 #ifdef LIBMESH_ENABLE_AMR 2435 std::vector<dof_id_type> & di,
2436 const unsigned int vn)
const 2438 LOG_SCOPE(
"old_dof_indices()",
"DofMap");
2440 libmesh_assert(elem);
2443 const unsigned int sys_num = this->
sys_number();
2445 const unsigned int dim = elem->
dim();
2446 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 2447 const bool is_inf = elem->
infinite();
2453 libmesh_assert(!elem->
has_dofs(sys_num) ||
2461 std::vector<const Node *> elem_nodes;
2462 const Node *
const * nodes_ptr;
2469 nodes_ptr = elem_nodes.data();
2470 n_nodes = cast_int<unsigned int>(elem_nodes.size());
2480 for (
unsigned int vg=0; vg<n_var_groups; vg++)
2483 const unsigned int vars_in_group = var.
n_variables();
2485 for (
unsigned int vig=0; vig<vars_in_group; vig++)
2487 const unsigned int v = var.
number(vig);
2495 std::vector<dof_id_type> di_new;
2497 di.insert( di.end(), di_new.begin(), di_new.end());
2507 int p_adjustment = 0;
2510 libmesh_assert_greater (elem->
p_level(), 0);
2522 const bool extra_hanging_dofs =
2529 for (
unsigned int n=0; n<
n_nodes; n++)
2531 const Node * node = nodes_ptr[n];
2533 libmesh_assert(old_dof_obj);
2539 const unsigned int nc =
2540 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 2544 ndan (type, fe_type.
order, n);
2546 const int n_comp = old_dof_obj->
n_comp_group(sys_num,vg);
2552 if (extra_hanging_dofs && !elem->
is_vertex(n))
2554 const int dof_offset = n_comp - nc;
2566 for (
int i=n_comp-1; i>=dof_offset; i--)
2569 old_dof_obj->
dof_number(sys_num, vg, vig, i, n_comp);
2578 for (
unsigned int i=0; i<nc; i++)
2581 old_dof_obj->
dof_number(sys_num, vg, vig, i, n_comp);
2598 libmesh_assert(old_dof_obj);
2600 const unsigned int n_comp =
2603 if (old_dof_obj->
n_systems() > sys_num &&
2607 for (
unsigned int i=0; i<nc; i++)
2610 old_dof_obj->
dof_number(sys_num, vg, vig, i, n_comp);
2628 #endif // LIBMESH_ENABLE_AMR 2631 #ifdef LIBMESH_ENABLE_CONSTRAINTS 2635 typedef std::set<dof_id_type> RCSet;
2638 RCSet dof_set (elem_dofs.begin(), elem_dofs.end());
2647 for (
const auto & dof : elem_dofs)
2651 DofConstraints::const_iterator
2665 for (
const auto & pr : constraint_row)
2666 if (!dof_set.count (pr.first))
2668 dof_set.insert (pr.first);
2680 elem_dofs.insert (elem_dofs.end(),
2681 dof_set.begin(), dof_set.end());
2692 #endif // LIBMESH_ENABLE_CONSTRAINTS 2705 std::ostringstream os;
2710 const char * may_equal =
" <= ";
2715 if (mat->need_full_sparsity_pattern())
2719 long double avg_n_nz = 0, avg_n_oz = 0;
2723 for (
const auto & val : *
_n_nz)
2725 max_n_nz =
std::max(max_n_nz, val);
2729 std::size_t n_nz_size =
_n_nz->size();
2735 avg_n_nz /=
std::max(n_nz_size,std::size_t(1));
2737 libmesh_assert(
_n_oz);
2739 for (
const auto & val : *
_n_oz)
2741 max_n_oz =
std::max(max_n_oz, val);
2745 std::size_t n_oz_size =
_n_oz->size();
2751 avg_n_oz /=
std::max(n_oz_size,std::size_t(1));
2754 os <<
" DofMap Sparsity\n Average On-Processor Bandwidth" 2755 << may_equal << avg_n_nz <<
'\n';
2757 os <<
" Average Off-Processor Bandwidth" 2758 << may_equal << avg_n_oz <<
'\n';
2760 os <<
" Maximum On-Processor Bandwidth" 2761 << may_equal << max_n_nz <<
'\n';
2763 os <<
" Maximum Off-Processor Bandwidth" 2764 << may_equal << max_n_oz << std::endl;
2766 #ifdef LIBMESH_ENABLE_CONSTRAINTS 2768 std::size_t n_constraints = 0, max_constraint_length = 0,
2770 long double avg_constraint_length = 0.;
2780 std::size_t rowsize = row.size();
2782 max_constraint_length =
std::max(max_constraint_length,
2784 avg_constraint_length += rowsize;
2791 this->
comm().
sum(n_constraints);
2793 this->
comm().
sum(avg_constraint_length);
2794 this->
comm().
max(max_constraint_length);
2796 os <<
" DofMap Constraints\n Number of DoF Constraints = " 2800 <<
" Number of Heterogenous Constraints= " << n_rhss;
2803 avg_constraint_length /= n_constraints;
2806 <<
" Average DoF Constraint Length= " << avg_constraint_length;
2809 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 2810 std::size_t n_node_constraints = 0, max_node_constraint_length = 0,
2812 long double avg_node_constraint_length = 0.;
2817 const Node * node = pr.first;
2822 std::size_t rowsize = row.
size();
2824 max_node_constraint_length =
std::max(max_node_constraint_length,
2826 avg_node_constraint_length += rowsize;
2827 n_node_constraints++;
2829 if (pr.second.second !=
Point(0))
2833 this->
comm().
sum(n_node_constraints);
2834 this->
comm().
sum(n_node_rhss);
2835 this->
comm().
sum(avg_node_constraint_length);
2836 this->
comm().
max(max_node_constraint_length);
2838 os <<
"\n Number of Node Constraints = " << n_node_constraints;
2841 <<
" Number of Heterogenous Node Constraints= " << n_node_rhss;
2842 if (n_node_constraints)
2844 avg_node_constraint_length /= n_node_constraints;
2845 os <<
"\n Maximum Node Constraint Length= " << max_node_constraint_length
2847 <<
" Average Node Constraint Length= " << avg_node_constraint_length;
2849 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 2853 #endif // LIBMESH_ENABLE_CONSTRAINTS 2859 template bool DofMap::is_evaluable<Elem>(
const Elem &,
unsigned int)
const;
2860 template bool DofMap::is_evaluable<Node>(
const Node &,
unsigned int)
const;
std::vector< VariableGroup > _variable_groups
void find_connected_dofs(std::vector< dof_id_type > &elem_dofs) const
std::unique_ptr< SparsityPattern::Build > _sp
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
Manages the family, order, etc. parameters for a given FE.
RefinementState refinement_flag() const
bool _error_on_cyclic_constraint
dof_id_type vg_dof_base(const unsigned int s, const unsigned int vg) const
bool _implicit_neighbor_dofs_initialized
void distribute_local_dofs_node_major(dof_id_type &next_free_dof, MeshBase &mesh)
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
void local_variable_indices(std::vector< dof_id_type > &idx, const MeshBase &mesh, unsigned int var_num) const
unsigned int n_variable_groups() const
bool _implicit_neighbor_dofs
dof_id_type dof_number(const unsigned int s, const unsigned int var, const unsigned int comp) const
DefaultCoupling & default_coupling()
A geometric point in (x,y,z) space associated with a DOF.
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
void set_old_dof_object()
dof_id_type n_SCALAR_dofs() const
std::set< GhostingFunctor * >::const_iterator coupling_functors_begin() const
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
const unsigned int invalid_uint
dof_id_type n_dofs_on_processor(const processor_id_type proc) const
unsigned int n_comp(const unsigned int s, const unsigned int var) const
unsigned int n_var_groups(const unsigned int s) const
void * _extra_sparsity_context
void extract_local_vector(const NumericVector< Number > &Ug, const std::vector< dof_id_type > &dof_indices, DenseVectorBase< Number > &Ue) const
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
void set_implicit_neighbor_dofs(bool implicit_neighbor_dofs)
void add_default_ghosting()
Maps between boundary ids and PeriodicBoundaryBase objects.
bool is_periodic_boundary(const boundary_id_type boundaryid) const
void allgather(const T &send, std::vector< T, A > &recv) const
virtual numeric_index_type size() const =0
RefinementState p_refinement_flag() const
std::vector< dof_id_type > _send_list
const FEType & variable_type(const unsigned int c) const
bool is_attached(SparseMatrix< Number > &matrix)
void attach_matrix(SparseMatrix< Number > &matrix)
void remove_ghosting_functor(GhostingFunctor &ghosting_functor)
The base class for all geometric element types.
IntRange< std::size_t > index_range(const std::vector< T > &vec)
uint8_t processor_id_type
unsigned int n_variables() const
DefaultCoupling & default_algebraic_ghosting()
AugmentSparsityPattern * _augment_sparsity_pattern
const Parallel::Communicator & comm() const
std::set< GhostingFunctor * > _algebraic_ghosting_functors
std::vector< dof_id_type > _end_old_df
unsigned int p_level() const
bool use_coupled_neighbor_dofs(const MeshBase &mesh) const
Utility class for defining generic ranges for threading.
void set_vg_dof_base(const unsigned int s, const unsigned int vg, const dof_id_type db)
virtual void augment_send_list(std::vector< dof_id_type > &send_list)=0
dof_id_type n_old_dofs() const
bool has_dofs(const unsigned int s=libMesh::invalid_uint) const
long double max(long double a, double b)
void add_coupling_functor(GhostingFunctor &coupling_functor, bool to_mesh=true)
std::vector< dof_id_type > _first_old_df
void old_dof_indices(const Elem *const elem, std::vector< dof_id_type > &di, const unsigned int vn=libMesh::invalid_uint) const
std::string get_info() const
static unsigned int max_order(const FEType &fe_t, const ElemType &el_t)
unsigned int sys_number() const
void SCALAR_dof_indices(std::vector< dof_id_type > &di, const unsigned int vn, const bool old_dofs=false) const
std::vector< dof_id_type > _first_scalar_df
AdjointDofConstraintValues _adjoint_constraint_values
std::vector< dof_id_type > * _n_nz
AugmentSendList * _augment_send_list
std::vector< dof_id_type > * _n_oz
bool need_full_sparsity_pattern
virtual node_iterator nodes_begin()=0
virtual element_iterator elements_begin()=0
processor_id_type n_processors() const
const dof_id_type n_nodes
dof_id_type n_dofs() const
A variable which is solved for in a System of equations.
void add_neighbors_to_send_list(MeshBase &mesh)
std::vector< dof_id_type > _first_old_scalar_df
const Variable & variable(const unsigned int c) const
void(* _extra_sparsity_function)(SparsityPattern::Graph &, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz, void *)
dof_id_type _n_SCALAR_dofs
void set_nonlocal_dof_objects(iterator_type objects_begin, iterator_type objects_end, MeshBase &mesh, dofobject_accessor objects)
unsigned int n_variables() const
static bool extra_hanging_dofs(const FEType &fe_t)
void pull_parallel_vector_data(const Communicator &comm, const MapToVectors &queries, RequestContainer &reqs, GatherFunctor &gather_data, ActionFunctor &act_on_data, const datum *example)
virtual unsigned int n_nodes() const =0
virtual unsigned int size() const =0
static const processor_id_type invalid_processor_id
unsigned int n_vars(const unsigned int s, const unsigned int vg) const
unsigned int n_systems() const
std::unique_ptr< SparsityPattern::Build > build_sparsity(const MeshBase &mesh) const
const Node *const * get_nodes() const
SimpleRange< I > as_range(const std::pair< I, I > &p)
bool is_constrained_dof(const dof_id_type dof) const
std::unordered_map< const Elem *, const CouplingMatrix * > map_type
A surface shell element used in mechanics calculations.
static const dof_id_type invalid_id
void distribute_dofs(MeshBase &)
const VariableGroup & variable_group(const unsigned int c) const
DofConstraints _dof_constraints
void reinit(MeshBase &mesh)
DofMap(const unsigned int sys_number, MeshBase &mesh)
CouplingMatrix * _dof_coupling
OStreamProxy err(std::cerr)
std::vector< DirichletBoundaries * > _adjoint_dirichlet_boundaries
bool active_on_subdomain(subdomain_id_type sid) const
T command_line_next(std::string name, T value)
void print_info(std::ostream &os=libMesh::out) const
void * _extra_send_list_context
An object whose state is distributed along a set of processors.
void set_n_comp_group(const unsigned int s, const unsigned int vg, const unsigned int ncomp)
bool is_evaluable(const DofObjectSubclass &obj, unsigned int var_num=libMesh::invalid_uint) const
dof_id_type end_dof() const
static n_dofs_at_node_ptr n_dofs_at_node_function(const unsigned int dim, const FEType &fe_t)
std::string enum_to_string(const T e)
std::vector< SparseMatrix< Number > *> _matrices
DofObject * elem_ptr(MeshBase &mesh, dof_id_type i) const
std::set< GhostingFunctor * >::const_iterator algebraic_ghosting_functors_end() const
void attach_dof_map(const DofMap &dof_map)
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
virtual numeric_index_type first_local_index() const =0
virtual void augment_sparsity_pattern(SparsityPattern::Graph &sparsity, std::vector< dof_id_type > &n_nz, std::vector< dof_id_type > &n_oz)=0
std::unique_ptr< PeriodicBoundaries > _periodic_boundaries
virtual bool need_full_sparsity_pattern() const
void remove_coupling_functor(GhostingFunctor &coupling_functor)
unsigned int(* n_dofs_at_node_ptr)(const ElemType, const Order, const unsigned int)
virtual void update_sparsity_pattern(const SparsityPattern::Graph &)
DofObject * old_dof_object
subdomain_id_type subdomain_id() const
virtual unsigned short dim() const =0
void invalidate_dofs(MeshBase &mesh) const
virtual bool is_vertex(const unsigned int i) const =0
std::unique_ptr< DefaultCoupling > _default_coupling
void swap(Iterator &lhs, Iterator &rhs)
DofConstraintValueMap _primal_constraint_values
std::vector< Variable > _variables
DofConstraints _stashed_dof_constraints
void remove_default_ghosting()
dof_id_type first_dof() const
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
std::map< dof_id_type, Real, std::less< dof_id_type >, Threads::scalable_allocator< std::pair< const dof_id_type, Real > > > DofConstraintRow
void parallel_reduce(const Range &range, Body &body)
std::set< GhostingFunctor * > _coupling_functors
std::vector< dof_id_type > _end_df
unsigned int number(unsigned int v) const
std::vector< dof_id_type > _first_df
std::pair< unsigned int, unsigned int > var_to_vg_and_offset(const unsigned int s, const unsigned int var) const
void(* _extra_send_list_function)(std::vector< dof_id_type > &, void *)
DofObject * node_ptr(MeshBase &mesh, dof_id_type i) const
std::set< GhostingFunctor * >::const_iterator coupling_functors_end() const
std::map< const Node *, Real, std::less< const Node * >, Threads::scalable_allocator< std::pair< const Node *const, Real > > > NodeConstraintRow
bool all_semilocal_indices(const std::vector< dof_id_type > &dof_indices) const
bool on_command_line(std::string arg)
virtual bool infinite() const =0
std::unique_ptr< DefaultCoupling > _default_evaluating
std::vector< unsigned int > _variable_group_numbers
void remove_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor)
std::set< GhostingFunctor * >::const_iterator algebraic_ghosting_functors_begin() const
void compute_sparsity(const MeshBase &)
void add_algebraic_ghosting_functor(GhostingFunctor &evaluable_functor, bool to_mesh=true)
unsigned int n_comp_group(const unsigned int s, const unsigned int vg) const
processor_id_type processor_id() const
virtual T el(const unsigned int i) const =0
OStreamProxy out(std::cout)
void set_error_on_cyclic_constraint(bool error_on_cyclic_constraint)
processor_id_type processor_id() const
virtual ElemType type() const =0
void add_variable_group(const VariableGroup &var_group)
A geometric point in (x,y,z) space.
virtual numeric_index_type last_local_index() const =0
void distribute_local_dofs_var_major(dof_id_type &next_free_dof, MeshBase &mesh)
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)
bool semilocal_index(dof_id_type dof_index) const
void add_ghosting_functor(GhostingFunctor &ghosting_functor)
bool local_index(dof_id_type dof_index) const
const FEType & type() const
NodeConstraints _node_constraints
Defines the coupling between variables of a System.