20 #include "libmesh/default_coupling.h" 21 #include "libmesh/dof_map.h" 22 #include "libmesh/eigen_system.h" 23 #include "libmesh/elem.h" 24 #include "libmesh/explicit_system.h" 25 #include "libmesh/fe_interface.h" 26 #include "libmesh/frequency_system.h" 27 #include "libmesh/int_range.h" 28 #include "libmesh/libmesh_logging.h" 29 #include "libmesh/linear_implicit_system.h" 30 #include "libmesh/mesh_base.h" 31 #include "libmesh/mesh_refinement.h" 32 #include "libmesh/newmark_system.h" 33 #include "libmesh/nonlinear_implicit_system.h" 34 #include "libmesh/parallel.h" 35 #include "libmesh/rb_construction.h" 36 #include "libmesh/remote_elem.h" 37 #include "libmesh/transient_rb_construction.h" 38 #include "libmesh/transient_system.h" 47 #include "libmesh/equation_systems.h" 55 _refine_in_reinit(true),
56 _enable_default_ghosting(true)
60 this->
parameters.
set<
unsigned int>(
"linear solver maximum iterations") = 5000;
104 const unsigned int n_sys = this->
n_systems();
106 libmesh_assert_not_equal_to (n_sys, 0);
110 for (
auto & node :
_mesh.node_ptr_range())
111 node->set_n_systems(n_sys);
113 for (
auto & elem :
_mesh.element_ptr_range())
114 elem->set_n_systems(n_sys);
119 #ifdef LIBMESH_ENABLE_AMR 133 parallel_object_only();
135 const unsigned int n_sys = this->
n_systems();
136 libmesh_assert_not_equal_to (n_sys, 0);
139 for (
unsigned int i=0; i != n_sys; ++i)
152 for (
auto & node :
_mesh.node_ptr_range())
153 node->set_n_systems(n_sys);
156 for (
auto & elem :
_mesh.element_ptr_range())
157 elem->set_n_systems(n_sys);
161 for (
unsigned int i=0; i != n_sys; ++i)
164 #ifdef LIBMESH_ENABLE_AMR 166 bool mesh_changed =
false;
172 for (
unsigned int i=0; i != n_sys; ++i)
243 #endif // #ifdef LIBMESH_ENABLE_AMR 264 const unsigned int n_sys = this->
n_systems();
266 libmesh_assert_not_equal_to (n_sys, 0);
273 for (
auto & node :
_mesh.node_ptr_range())
274 node->set_n_systems(n_sys);
276 for (
auto & elem :
_mesh.element_ptr_range())
277 elem->set_n_systems(n_sys);
302 mesh.add_ghosting_functor(
mesh.default_ghosting());
304 mesh.remove_ghosting_functor(
mesh.default_ghosting());
320 LOG_SCOPE(
"update()",
"EquationSystems");
330 std::string_view
name)
341 else if (sys_type ==
"Basic")
342 this->add_system<System> (
name);
345 else if (sys_type ==
"Newmark")
346 this->add_system<NewmarkSystem> (
name);
349 else if ((sys_type ==
"Explicit"))
350 this->add_system<ExplicitSystem> (
name);
353 else if ((sys_type ==
"Implicit") ||
354 (sys_type ==
"Steady" ))
355 this->add_system<ImplicitSystem> (
name);
358 else if ((sys_type ==
"Transient") ||
359 (sys_type ==
"TransientImplicit") ||
360 (sys_type ==
"TransientLinearImplicit"))
361 this->add_system<TransientLinearImplicitSystem> (
name);
364 else if (sys_type ==
"TransientNonlinearImplicit")
365 this->add_system<TransientNonlinearImplicitSystem> (
name);
368 else if (sys_type ==
"TransientExplicit")
369 this->add_system<TransientExplicitSystem> (
name);
372 else if (sys_type ==
"LinearImplicit")
373 this->add_system<LinearImplicitSystem> (
name);
376 else if (sys_type ==
"NonlinearImplicit")
377 this->add_system<NonlinearImplicitSystem> (
name);
380 else if (sys_type ==
"RBConstruction")
381 this->add_system<RBConstruction> (
name);
384 else if (sys_type ==
"TransientRBConstruction")
385 this->add_system<TransientRBConstruction> (
name);
387 #ifdef LIBMESH_HAVE_SLEPC 389 else if (sys_type ==
"Eigen")
390 this->add_system<EigenSystem> (
name);
391 else if (sys_type ==
"TransientEigenSystem")
392 this->add_system<TransientEigenSystem> (
name);
395 #if defined(LIBMESH_USE_COMPLEX_NUMBERS) 397 else if (sys_type ==
"Frequency")
398 this->add_system<FrequencySystem> (
name);
402 libmesh_error_msg(
"ERROR: Unknown system type: " << sys_type);
426 this->
get_system(i).sensitivity_solve(parameters_in);
435 for (
unsigned int i=this->
n_systems(); i != 0; --i)
436 this->
get_system(i-1).adjoint_solve(qoi_indices);
443 const std::set<std::string> * system_names)
const 446 unsigned int var_num = var_names.size();
450 std::unordered_multiset<std::string> seen_names;
456 unsigned int n_scalar_vars = 0;
457 unsigned int n_vector_vars = 0;
459 for (
const auto & [sys_name, sys_ptr] :
_systems)
462 bool use_current_system = (system_names ==
nullptr);
463 if (!use_current_system)
464 use_current_system = system_names->count(sys_name);
465 if (!use_current_system || sys_ptr->hide_output())
468 seen_names.insert(sys_ptr->variable_name(vn));
474 seen_names.insert(sys_ptr->variable_name(vn));
485 unsigned int nv = n_scalar_vars +
dim*n_vector_vars;
489 libmesh_assert_less_equal ( nv, (
dim > 0 ?
dim : 1)*this->
n_vars() );
494 var_names.resize( nv );
497 for (
const auto & [sys_name, sys_ptr] :
_systems)
500 bool use_current_system = (system_names ==
nullptr);
501 if (!use_current_system)
502 use_current_system = system_names->count(sys_name);
503 if (!use_current_system || sys_ptr->hide_output())
508 const std::string & var_name = sys_ptr->variable_name(vn);
509 const FEType & fe_type = sys_ptr->variable_type(vn);
514 if (type ==
nullptr || (type && *type == fe_type))
522 var_names[var_num++] = var_name;
523 libmesh_error_msg_if(seen_names.count(var_name) > 1,
524 "Duplicate variable name "+var_name);
527 var_names[var_num++] = var_name+
"_x";
528 var_names[var_num++] = var_name+
"_y";
529 libmesh_error_msg_if(seen_names.count(var_name+
"_x"),
530 "Duplicate variable name "+var_name+
"_x");
531 libmesh_error_msg_if(seen_names.count(var_name+
"_y"),
532 "Duplicate variable name "+var_name+
"_y");
535 var_names[var_num++] = var_name+
"_x";
536 var_names[var_num++] = var_name+
"_y";
537 var_names[var_num++] = var_name+
"_z";
538 libmesh_error_msg_if(seen_names.count(var_name+
"_x"),
539 "Duplicate variable name "+var_name+
"_x");
540 libmesh_error_msg_if(seen_names.count(var_name+
"_y"),
541 "Duplicate variable name "+var_name+
"_y");
542 libmesh_error_msg_if(seen_names.count(var_name+
"_z"),
543 "Duplicate variable name "+var_name+
"_z");
546 libmesh_error_msg(
"Invalid dim in build_variable_names");
550 var_names[var_num++] = var_name;
555 var_names.resize(var_num);
562 std::string_view)
const 565 libmesh_not_implemented();
571 std::unique_ptr<NumericVector<Number>>
573 bool add_sides)
const 575 LOG_SCOPE(
"build_parallel_solution_vector()",
"EquationSystems");
578 parallel_object_only();
598 unsigned int n_scalar_vars = 0;
599 unsigned int n_vector_vars = 0;
600 for (
const auto & [sys_name, sys_ptr] :
_systems)
603 bool use_current_system = (system_names ==
nullptr);
604 if (!use_current_system)
605 use_current_system = system_names->count(sys_name);
606 if (!use_current_system || sys_ptr->hide_output())
619 nv = n_scalar_vars +
dim*n_vector_vars;
625 _mesh.local_nodes_end()));
636 const dof_id_type gaps_per_processor = n_gaps / n_proc;
637 const dof_id_type remainder_gaps = n_gaps % n_proc;
639 n_local_nodes = n_local_nodes +
641 (my_pid < remainder_gaps);
646 added_side_nodes = 0;
649 std::vector<dof_id_type> others_added_side_nodes;
653 std::map<std::tuple<dof_id_type, unsigned short, unsigned short>,
661 std::vector<dof_id_type> true_node_offsets;
663 std::vector<dof_id_type> added_node_offsets;
665 auto node_id_to_vec_id =
666 [&true_node_offsets, &added_node_offsets]
669 if (true_node_offsets.empty())
673 const auto lb = std::upper_bound(true_node_offsets.begin(),
674 true_node_offsets.end(), node_id);
678 return node_id + added_node_offsets[p];
683 true_node_offsets.resize(n_proc);
684 added_node_offsets.resize(n_proc);
687 for (
const auto & elem :
_mesh.active_element_ptr_range())
689 for (
auto s : elem->side_index_range())
694 const std::vector<unsigned int> side_nodes =
695 elem->nodes_on_side(s);
698 local_added_side_nodes += side_nodes.size();
702 others_added_side_nodes.resize(n_proc);
704 others_added_side_nodes);
706 added_side_nodes = std::accumulate(others_added_side_nodes.begin(),
707 others_added_side_nodes.end(), 0,
712 true_node_offsets[p+1] += true_node_offsets[p];
717 added_node_offsets[0] = 0;
719 added_node_offsets[p+1] =
720 added_node_offsets[p] + others_added_side_nodes[p];
725 dof_id_type node_counter = true_node_offsets[my_pid];
727 node_counter += others_added_side_nodes[p];
730 for (
const auto & elem :
_mesh.active_local_element_ptr_range())
732 for (
auto s : elem->side_index_range())
737 const std::vector<unsigned int> side_nodes =
738 elem->nodes_on_side(s);
741 discontinuous_node_indices
742 [std::make_tuple(elem->id(),s,n)] = node_counter++;
748 n_global_vals = (max_nn + added_side_nodes) * nv,
749 n_local_vals = (n_local_nodes + local_added_side_nodes) * nv;
754 parallel_soln.
init(n_global_vals, n_local_vals,
false,
PARALLEL);
760 repeat_count.
init(n_global_vals, n_local_vals,
false,
PARALLEL);
762 repeat_count.
close();
764 unsigned int var_num=0;
771 for (
const auto & [sys_name, sys_ptr] :
_systems)
774 bool use_current_system = (system_names ==
nullptr);
775 if (!use_current_system)
776 use_current_system = system_names->count(sys_name);
777 if (!use_current_system || sys_ptr->hide_output())
780 const System & system = *sys_ptr;
781 const unsigned int nv_sys = system.
n_vars();
782 const unsigned int sys_num = system.
number();
785 unsigned int n_scalar_vars = 0;
786 unsigned int n_vector_vars = 0;
797 unsigned int nv_sys_split = n_scalar_vars +
dim*n_vector_vars;
808 if (!non_const_sys.
solution->closed())
817 std::vector<Number> elem_soln;
818 std::vector<Number> nodal_soln;
819 std::vector<dof_id_type> dof_indices;
821 unsigned var_inc = 0;
822 for (
unsigned int var=0; var<nv_sys; var++)
830 for (
const auto & elem :
_mesh.active_local_element_ptr_range())
835 sys_soln.
get(dof_indices, elem_soln);
846 if (!elem->infinite())
848 libmesh_assert_equal_to (nodal_soln.size(), n_vec_dim*elem->n_nodes());
850 for (
auto n : elem->node_index_range())
852 const Node & node = elem->node_ref(n);
855 nv * node_id_to_vec_id(node.
id());
857 for (
unsigned int d=0; d < n_vec_dim; d++)
861 parallel_soln.
add(node_idx + (var_inc+d + var_num), nodal_soln[n_vec_dim*n+d]);
864 repeat_count.
add(node_idx + (var_inc+d + var_num), 1);
870 for (
auto s : elem->side_index_range())
878 (fe_type, elem, s, elem_soln,
879 nodal_soln, add_p_level, n_vec_dim);
882 const std::vector<unsigned int> side_nodes =
883 elem->nodes_on_side(s);
885 libmesh_assert_equal_to
893 std::size_t node_index =
894 nv * libmesh_map_find(discontinuous_node_indices,
895 std::make_tuple(elem->id(), s, n));
897 for (
unsigned int d=0; d < n_vec_dim; d++)
899 parallel_soln.
add(node_index + (var_inc+d + var_num), nodal_soln[n_vec_dim*n+d]);
900 repeat_count.
add(node_index + (var_inc+d + var_num), 1);
908 for (
auto n : elem->node_index_range())
910 const Node & node = elem->node_ref(n);
914 if (!node.
n_dofs(sys_num, var))
917 nv * node_id_to_vec_id(node.
id());
919 for (
unsigned int d=0; d < n_vec_dim; d++)
920 repeat_count.
add(node_idx + (var_inc+d + var_num), 1);
925 var_inc += n_vec_dim;
928 var_num += nv_sys_split;
932 parallel_soln.
close();
933 repeat_count.
close();
947 if (std::abs(repeat_count(i)) <
TOLERANCE)
948 repeat_count.
set(i, 1.);
953 repeat_count.
close();
957 parallel_soln /= repeat_count;
959 return parallel_soln_ptr;
965 const std::set<std::string> * system_names,
966 bool add_sides)
const 968 LOG_SCOPE(
"build_solution_vector()",
"EquationSystems");
971 std::unique_ptr<NumericVector<Number>> parallel_soln =
975 parallel_soln->localize_to_one(soln);
981 std::vector<std::set<subdomain_id_type>> & vars_active_subdomains)
const 983 vars_active_subdomains.clear();
984 vars_active_subdomains.resize(names.size());
988 const auto & sys_ptr = pr.second;
991 const std::string & var_name = sys_ptr->variable_name(vn);
993 auto names_it = std::find(names.begin(), names.end(), var_name);
994 if(names_it != names.end())
996 const Variable & variable = sys_ptr->variable(vn);
997 const std::set<subdomain_id_type> & active_subdomains = variable.
active_subdomains();
998 vars_active_subdomains[
std::distance(names.begin(), names_it)] = active_subdomains;
1006 #ifdef LIBMESH_ENABLE_DEPRECATED 1008 std::vector<std::string> & names)
const 1010 libmesh_deprecated();
1013 #endif // LIBMESH_ENABLE_DEPRECATED 1019 std::vector<std::string> & names)
const 1022 std::unique_ptr<NumericVector<Number>> parallel_soln =
1032 parallel_soln->localize_to_one(soln);
1035 std::vector<std::pair<unsigned int, unsigned int>>
1037 (std::vector<std::string> & names,
const FEType * type,
const std::vector<FEType> * types)
const 1040 parallel_object_only();
1045 libmesh_assert_msg(!type || !types,
1046 "Input 'type', 'types', or neither in find_variable_numbers, but not both.");
1048 std::vector<FEType> type_filter;
1050 type_filter.push_back(*type);
1052 type_filter = *types;
1058 std::vector<std::string> name_filter = names;
1059 bool is_names_empty = name_filter.empty();
1066 const std::vector<std::string> component_suffix = {
"_x",
"_y",
"_z"};
1067 unsigned int dim = _mesh.spatial_dimension();
1068 libmesh_error_msg_if(
dim > 3,
"Invalid dim in find_variable_numbers");
1072 std::vector<std::pair<unsigned int, unsigned int>> var_nums;
1074 for (
const auto & pr : _systems)
1076 const System & system = *(pr.second);
1082 if (type_filter.size() &&
1083 std::find(type_filter.begin(), type_filter.end(), var_type) == type_filter.end())
1089 std::vector<std::string> component_names;
1090 for (
unsigned int comp = 0; comp <
dim; ++comp)
1093 if (is_names_empty ||
1094 (std::find(name_filter.begin(), name_filter.end(),
name) != name_filter.end()))
1095 component_names.push_back(
name);
1098 if (! component_names.empty())
1099 names.insert(names.end(), component_names.begin(), component_names.end());
1106 if (is_names_empty ||
1107 (std::find(name_filter.begin(), name_filter.end(),
name) != name_filter.end()))
1108 names.push_back(
name);
1114 var_nums.emplace_back(system.
number(), var);
1119 std::vector<unsigned int> sort_index(var_nums.size());
1120 std::iota(sort_index.begin(), sort_index.end(), 0);
1121 std::sort(sort_index.begin(), sort_index.end(),
1122 [&](
const unsigned int & lhs,
const unsigned int & rhs)
1123 {
return this->get_system(var_nums[lhs].first).variable_name(var_nums[lhs].second) <
1124 this->get_system(var_nums[rhs].first).variable_name(var_nums[rhs].second);});
1126 std::vector<std::pair<unsigned int, unsigned int>> var_nums_sorted(var_nums.size());
1129 var_nums_sorted[i].first = var_nums[sort_index[i]].first;
1130 var_nums_sorted[i].second = var_nums[sort_index[i]].second;
1134 std::sort(names.begin(), names.end());
1137 return var_nums_sorted;
1141 std::unique_ptr<NumericVector<Number>>
1148 std::vector<std::pair<unsigned int, unsigned int>> var_nums =
1151 const std::size_t nv = names.size();
1157 return std::unique_ptr<NumericVector<Number>>(
nullptr);
1171 parallel_soln_local_size = div+1;
1176 parallel_soln.
init(parallel_soln_global_size,
1177 parallel_soln_local_size,
1181 unsigned int sys_ctr = 0;
1182 unsigned int var_ctr = 0;
1185 std::pair<unsigned int, unsigned int> var_num = var_nums[i];
1189 if (sys_ctr != var_num.first)
1198 if (!non_const_sys.
solution->closed())
1201 sys_ctr = var_num.first;
1207 std::vector<dof_id_type> dof_indices;
1209 const unsigned int var = var_num.second;
1218 const unsigned int n_comps =
1222 for (
const auto & elem :
_mesh.active_local_element_ptr_range())
1229 libmesh_assert_equal_to(dof_indices.size(), n_comps);
1231 for (
unsigned int comp = 0; comp < n_comps; comp++)
1232 parallel_soln.
set(ne * (var_ctr + comp) + elem->id(), sys_soln(dof_indices[comp]));
1243 libmesh_assert_equal_to(names.size(), var_ctr);
1245 parallel_soln.
close();
1246 return parallel_soln_ptr;
1253 (std::vector<Number> & soln,
1254 const std::set<std::string> * system_names,
1255 const std::vector<std::string> * var_names,
1257 bool add_sides)
const 1259 LOG_SCOPE(
"build_discontinuous_solution_vector()",
"EquationSystems");
1265 unsigned int nv = 0;
1267 for (
const auto & [sys_name, sys_ptr] : _systems)
1270 bool use_current_system = (system_names ==
nullptr);
1271 if (!use_current_system)
1272 use_current_system = system_names->count(sys_name);
1273 if (!use_current_system || sys_ptr->hide_output())
1278 for (
auto var_id :
make_range(sys_ptr->n_vars()))
1280 bool use_current_var = (var_names ==
nullptr);
1281 if (!use_current_var)
1282 use_current_var = std::count(var_names->begin(),
1284 sys_ptr->variable_name(var_id));
1288 if (use_current_var)
1296 for (
const auto & elem : _mesh.active_element_ptr_range())
1298 tw += vertices_only ? elem->n_vertices() : elem->n_nodes();
1302 for (
auto s : elem->side_index_range())
1304 if (redundant_added_side(*elem,s))
1307 const std::vector<unsigned int> side_nodes =
1308 elem->nodes_on_side(s);
1311 tw += side_nodes.size();
1314 if (elem->is_vertex(side_nodes[n]))
1322 if (_mesh.processor_id() == 0)
1325 std::vector<Number> sys_soln;
1329 unsigned int var_offset = 0;
1336 for (
const auto & [sys_name, system] : _systems)
1339 bool use_current_system = (system_names ==
nullptr);
1340 if (!use_current_system)
1341 use_current_system = system_names->count(sys_name);
1342 if (!use_current_system || system->hide_output())
1345 const unsigned int nv_sys = system->n_vars();
1346 const auto & dof_map = system->get_dof_map();
1348 system->update_global_solution (sys_soln, 0);
1351 unsigned int n_vars_written_current_system = 0;
1353 if (_mesh.processor_id() == 0)
1355 std::vector<Number> soln_coeffs;
1356 std::vector<Number> nodal_soln;
1357 std::vector<dof_id_type> dof_indices;
1364 for (
unsigned int var=0; var<nv_sys; var++)
1366 bool use_current_var = (var_names ==
nullptr);
1367 if (!use_current_var)
1368 use_current_var = std::count(var_names->begin(),
1370 system->variable_name(var));
1374 if (!use_current_var)
1377 const FEType & fe_type = system->variable_type(var);
1378 const Variable & var_description = system->variable(var);
1379 const auto vg = dof_map.var_group_from_var_number(var);
1380 const bool add_p_level = dof_map.should_p_refine(vg);
1384 for (
auto & elem : _mesh.active_element_ptr_range())
1388 system->get_dof_map().dof_indices (elem, dof_indices, var);
1390 soln_coeffs.resize(dof_indices.size());
1393 soln_coeffs[i] = sys_soln[dof_indices[i]];
1407 if (!elem->infinite())
1409 libmesh_assert_equal_to (nodal_soln.size(), elem->n_nodes());
1411 const unsigned int n_vals =
1412 vertices_only ? elem->n_vertices() : elem->n_nodes();
1414 for (
unsigned int n=0; n<n_vals; n++)
1418 nv * (nn++) + (n_vars_written_current_system + var_offset);
1420 soln[index] += nodal_soln[n];
1425 nn += vertices_only ? elem->n_vertices() : elem->n_nodes();
1436 std::vector<std::vector<const Elem *>>
1437 elems_by_pid(_mesh.n_processors());
1439 for (
const auto & elem : _mesh.active_element_ptr_range())
1440 elems_by_pid[elem->processor_id()].push_back(elem);
1443 for (
const Elem * elem : elems_by_pid[p])
1447 system->get_dof_map().dof_indices (elem, dof_indices, var);
1449 soln_coeffs.resize(dof_indices.size());
1452 soln_coeffs[i] = sys_soln[dof_indices[i]];
1454 for (
auto s : elem->side_index_range())
1456 if (redundant_added_side(*elem,s))
1459 const std::vector<unsigned int> side_nodes =
1460 elem->nodes_on_side(s);
1467 (fe_type, elem, s, soln_coeffs,
1468 nodal_soln, add_p_level,
1471 libmesh_assert_equal_to
1487 neigh->
level() == elem->level() &&
1490 std::vector<dof_id_type> neigh_indices;
1491 system->get_dof_map().dof_indices (neigh, neigh_indices, var);
1492 std::vector<Number> neigh_coeffs(neigh_indices.size());
1495 neigh_coeffs[i] = sys_soln[neigh_indices[i]];
1497 const unsigned int s_neigh =
1499 std::vector<Number> neigh_soln;
1501 (fe_type, neigh, s_neigh,
1502 neigh_coeffs, neigh_soln, add_p_level,
1505 const std::vector<unsigned int> neigh_nodes =
1509 if (neigh->
node_ptr(neigh_nodes[neigh_n])
1510 == elem->node_ptr(side_nodes[n]))
1512 nodal_soln[n] += neigh_soln[neigh_n];
1519 if (vertices_only &&
1520 !elem->is_vertex(n))
1525 nv * (nn++) + (n_vars_written_current_system + var_offset);
1527 soln[index] += nodal_soln[n];
1533 nn += vertices_only ? elem->n_vertices() : elem->n_nodes();
1535 for (
auto s : elem->side_index_range())
1537 if (redundant_added_side(*elem,s))
1540 const std::vector<unsigned int> side_nodes =
1541 elem->nodes_on_side(s);
1545 if (vertices_only &&
1546 !elem->is_vertex(n))
1556 n_vars_written_current_system++;
1562 var_offset += n_vars_written_current_system;
1583 if (!neigh->active())
1588 return (neigh->id() < elem.
id());
1594 const Real threshold,
1595 const bool verbose)
const 1599 std::vector<bool> os_result;
1605 libMesh::out <<
" Fatal difference. This system handles " 1606 << this->
n_systems() <<
" systems," << std::endl
1607 <<
" while the other system handles " 1609 <<
" systems." << std::endl
1610 <<
" Aborting comparison." << std::endl;
1617 for (
const auto & [sys_name, sys_ptr] :
_systems)
1622 os_result.push_back (sys_ptr->compare (other_system, threshold, verbose));
1630 if (os_result.size()==0)
1638 os_identical = os_result[n];
1641 while (os_identical && n<os_result.size());
1642 return os_identical;
1650 std::ostringstream oss;
1652 oss <<
" EquationSystems\n" 1653 <<
" n_systems()=" << this->
n_systems() <<
'\n';
1657 oss << pr.second->get_info();
1703 tot += pr.second->n_vars();
1715 tot += pr.second->n_dofs();
1728 tot += pr.second->n_active_dofs();
1737 for (
auto & node :
_mesh.node_ptr_range())
1741 for (
auto & elem :
_mesh.element_ptr_range())
1747 this->
get_system(sys_num).get_dof_map().remove_default_ghosting();
EquationSystems(MeshBase &mesh)
Constructor.
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
void allgather(const T &send_data, std::vector< T, A > &recv_data) const
This is the EquationSystems class.
unsigned int n_vars() const
void build_variable_names(std::vector< std::string > &var_names, const FEType *type=nullptr, const std::set< std::string > *system_names=nullptr) const
Fill the input vector var_names with the names of the variables for each system.
A Node is like a Point, but with more information.
unsigned int n_systems() const
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
bool _enable_default_ghosting
Flag for whether to enable default ghosting on newly added Systems.
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
Access multiple components at once.
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
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.
std::unique_ptr< NumericVector< Number > > build_parallel_elemental_solution_vector(std::vector< std::string > &names) const
Builds a parallel vector of CONSTANT MONOMIAL and/or components of CONSTANT MONOMIAL_VEC solution val...
std::size_t distribute_dofs(MeshBase &)
Distribute dofs on the current mesh.
static constexpr Real TOLERANCE
std::size_t n_dofs() const
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
void build_discontinuous_solution_vector(std::vector< Number > &soln, const std::set< std::string > *system_names=nullptr, const std::vector< std::string > *var_names=nullptr, bool vertices_only=false, bool add_sides=false) const
Fill the input vector soln with solution values.
void add_default_ghosting()
Add the default functor(s) for coupling and algebraic ghosting.
bool refine_elements()
Only refines the user-requested elements.
bool reinit_solutions()
Handle any mesh changes and project any solutions onto the updated mesh.
virtual void clear()
Restores the data structure to a pristine state.
virtual void allgather()
Gathers all elements and nodes of the mesh onto every processor.
virtual ~EquationSystems()
Destructor.
This is the base class from which all geometric element types are derived.
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
static FEFieldType field_type(const FEType &fe_type)
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
processor_id_type rank() const
const Parallel::Communicator & comm() const
virtual void enable_default_ghosting(bool enable)
Enable or disable default ghosting functors on the Mesh and on all Systems.
void clean_refinement_flags()
Sets the refinement flag to Elem::DO_NOTHING for each element in the mesh.
void should_p_refine(unsigned int g, bool p_refine)
Describe whether the given variable group should be p-refined.
virtual void adjoint_solve(const QoISet &qoi_indices=QoISet())
Call adjoint_solve on all the individual equation systems.
const Parallel::Communicator & _communicator
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(std::string_view name) const
Real distance(const Point &p)
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
Change the dimension of the vector to n.
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
static void side_nodal_soln(const FEType &fe_t, const Elem *elem, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, const bool add_p_level=true, const unsigned int vdim=1)
Build the nodal soln on one side from the (full) element soln.
uint8_t processor_id_type
This is the MeshBase class.
virtual std::string get_info() const
virtual void prolong_vectors()
Prolong vectors after the mesh has refined.
MeshBase & _mesh
The mesh data structure.
bool coarsen_elements()
Only coarsens the user-requested elements.
unsigned int var_group_from_var_number(unsigned int var_num) const
signed char & underrefined_boundary_limit()
If underrefined_boundary_limit is set to a nonnegative value, then refinement and coarsening will pro...
unsigned int n_dofs(const unsigned int s, const unsigned int var=libMesh::invalid_uint) const
This class handles the numbering of degrees of freedom on a mesh.
void allgather()
Serializes a distributed mesh and its associated degree of freedom numbering for all systems...
void get_vars_active_subdomains(const std::vector< std::string > &names, std::vector< std::set< subdomain_id_type >> &vars_active_subdomains) const
Retrieve vars_active_subdomains, which indicates the active subdomains for each variable in names...
processor_id_type size() const
virtual void sensitivity_solve(const ParameterVector ¶meters)
Call sensitivity_solve on all the individual equation systems.
virtual void reinit_systems()
Reinitialize all systems on the current mesh.
processor_id_type n_processors() const
virtual bool is_serial() const
unsigned int number() const
This class defines the notion of a variable in the system.
const std::set< subdomain_id_type > & active_subdomains() const
Implements (adaptive) mesh refinement algorithms for a MeshBase.
dof_id_type numeric_index_type
unsigned int which_neighbor_am_i(const Elem *e) const
This function tells you which neighbor e is.
Manages consistently variables, degrees of freedom, and coefficient vectors.
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
virtual dof_id_type max_elem_id() const =0
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
const std::string & variable_name(const unsigned int i) const
void build_solution_vector(std::vector< Number > &soln, std::string_view system_name, std::string_view variable_name="all_vars") const
Fill the input vector soln with the solution values for the system named name.
virtual void restrict_vectors()
Restrict vectors after the mesh has coarsened.
bool active_on_subdomain(subdomain_id_type sid) const
bool _refine_in_reinit
Flag for whether to call coarsen/refine in reinit().
virtual void reinit_mesh()
Handle the association of a completely new mesh with the EquationSystem and all the Systems assigned ...
unsigned char & face_level_mismatch_limit()
If face_level_mismatch_limit is set to a nonzero value, then refinement and coarsening will produce m...
An object whose state is distributed along a set of processors.
static unsigned int n_vec_dim(const MeshBase &mesh, const FEType &fe_type)
virtual void reinit_constraints()
Reinitializes the constraints for this system.
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
static FEContinuity get_continuity(const FEType &fe_type)
Returns the input FEType's FEContinuity based on the underlying FEFamily and potentially the Order...
virtual void solve()
Call solve on all the individual equation systems.
virtual numeric_index_type first_local_index() const =0
const Elem * neighbor_ptr(unsigned int i) const
void get_solution(std::vector< Number > &soln, std::vector< std::string > &names) const
Retrieve the solution data for CONSTANT MONOMIALs and/or components of CONSTANT MONOMIAL_VECs.
virtual void update()
Update the local values to reflect the solution on neighboring processors.
unsigned int level() const
const FEType & variable_type(const unsigned int i) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool contract()=0
Delete subactive (i.e.
subdomain_id_type subdomain_id() const
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
const Node * node_ptr(const unsigned int i) const
unsigned int spatial_dimension() const
T & set(const std::string &)
std::unique_ptr< NumericVector< Number > > build_parallel_solution_vector(const std::set< std::string > *system_names=nullptr, bool add_sides=false) const
A version of build_solution_vector which is appropriate for "parallel" output formats like Nemesis...
virtual void clear()
Clears internal data structures & frees any allocated memory.
signed char & overrefined_boundary_limit()
If overrefined_boundary_limit is set to a nonnegative value, then refinement and coarsening will prod...
void remove_default_ghosting()
Remove any default ghosting functor(s).
void _remove_default_ghosting(unsigned int sys_num)
This just calls DofMap::remove_default_ghosting() but using a shim lets us forward-declare DofMap...
const MeshBase & get_mesh() const
std::size_t n_active_dofs() const
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC)
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
std::vector< std::pair< unsigned int, unsigned int > > find_variable_numbers(std::vector< std::string > &names, const FEType *type=nullptr, const std::vector< FEType > *types=nullptr) const
Finds system and variable numbers for any variables of 'type' or of 'types' corresponding to the entr...
Parameters parameters
Data structure holding arbitrary parameters.
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
void prepare_send_list()
Takes the _send_list vector (which may have duplicate entries) and sorts it.
virtual void init()
Initialize all the systems.
virtual bool compare(const EquationSystems &other_es, const Real threshold, const bool verbose) const
unsigned int n_vars() const
virtual dof_id_type max_node_id() const =0
virtual dof_id_type n_elem() const =0
virtual void add(const numeric_index_type i, const T value)=0
Adds value to the vector entry specified by i.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
processor_id_type processor_id() const
static void nodal_soln(const unsigned int dim, const FEType &fe_t, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, const bool add_p_level=true, const unsigned int vdim=1)
Build the nodal soln from the element soln.
void _add_system_to_nodes_and_elems()
This function is used in the implementation of add_system, it loops over the nodes and elements of th...
const DofMap & get_dof_map() const
virtual numeric_index_type last_local_index() const =0
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
void update()
Updates local values for all the systems.
virtual std::vector< unsigned int > nodes_on_side(const unsigned int) const =0
void build_elemental_solution_vector(std::vector< Number > &soln, std::vector< std::string > &names) const
Retrieve the solution data for CONSTANT MONOMIALs and/or components of CONSTANT MONOMIAL_VECs.
std::map< std::string, std::unique_ptr< System >, std::less<> > _systems
Data structure holding the systems.
static bool redundant_added_side(const Elem &elem, unsigned int side)
const RemoteElem * remote_elem