Go to the documentation of this file.
24 #include "libmesh/dof_map.h"
25 #include "libmesh/equation_systems.h"
26 #include "libmesh/int_range.h"
27 #include "libmesh/libmesh_logging.h"
28 #include "libmesh/mesh_base.h"
29 #include "libmesh/numeric_vector.h"
30 #include "libmesh/parameter_vector.h"
31 #include "libmesh/point.h"
32 #include "libmesh/point_locator_base.h"
33 #include "libmesh/qoi_set.h"
34 #include "libmesh/enum_to_string.h"
35 #include "libmesh/system.h"
36 #include "libmesh/system_norm.h"
37 #include "libmesh/utility.h"
38 #include "libmesh/elem.h"
39 #include "libmesh/fe_type.h"
40 #include "libmesh/fe_interface.h"
41 #include "libmesh/fe_compute_data.h"
42 #include "libmesh/auto_ptr.h"
45 #include "libmesh/fe_base.h"
46 #include "libmesh/fe_interface.h"
47 #include "libmesh/parallel.h"
48 #include "libmesh/parallel_algebra.h"
49 #include "libmesh/quadrature.h"
50 #include "libmesh/tensor_value.h"
51 #include "libmesh/vector_value.h"
52 #include "libmesh/tensor_tools.h"
53 #include "libmesh/enum_norm_type.h"
62 const std::string & name_in,
63 const unsigned int number_in) :
66 assemble_before_solve (true),
67 use_fixed_solution (false),
68 extra_quadrature_order (0),
73 _init_system_function (nullptr),
74 _init_system_object (nullptr),
75 _assemble_system_function (nullptr),
76 _assemble_system_object (nullptr),
77 _constrain_system_function (nullptr),
78 _constrain_system_object (nullptr),
79 _qoi_evaluate_function (nullptr),
80 _qoi_evaluate_object (nullptr),
81 _qoi_evaluate_derivative_function (nullptr),
82 _qoi_evaluate_derivative_object (nullptr),
83 _dof_map (libmesh_make_unique<
DofMap>(number_in, es.get_mesh())),
84 _equation_systems (es),
85 _mesh (es.get_mesh()),
87 _sys_number (number_in),
89 _solution_projection (true),
90 _basic_system_only (false),
92 _identify_variable_groups (true),
93 _additional_data_written (false),
94 adjoint_already_solved (false),
105 _equation_systems(other._equation_systems),
107 _sys_number(other._sys_number)
109 libmesh_not_implemented();
116 libmesh_not_implemented();
159 #ifdef LIBMESH_ENABLE_CONSTRAINTS
161 return _dof_map->n_constrained_dofs();
174 #ifdef LIBMESH_ENABLE_CONSTRAINTS
176 return _dof_map->n_local_constrained_dofs();
197 libmesh_assert_less (global_dof_number,
_dof_map->n_dofs());
283 #ifdef LIBMESH_ENABLE_GHOSTED
302 #ifdef LIBMESH_ENABLE_GHOSTED
307 libmesh_error_msg(
"Cannot initialize ghosted vectors when they are not enabled.");
312 pr.second->init (this->
n_dofs(),
false, type);
316 libmesh_assert_equal_to(type,
PARALLEL);
326 #ifdef LIBMESH_ENABLE_AMR
342 #ifdef LIBMESH_ENABLE_GHOSTED
347 libmesh_error_msg(
"Cannot initialize ghosted vectors when they are not enabled.");
355 const std::vector<dof_id_type> & send_list =
_dof_map->get_send_list ();
364 #ifdef LIBMESH_ENABLE_GHOSTED
375 #endif // LIBMESH_ENABLE_AMR
382 #ifdef LIBMESH_ENABLE_AMR
399 #ifdef LIBMESH_ENABLE_CONSTRAINTS
412 const std::vector<dof_id_type> & send_list =
_dof_map->get_send_list ();
418 libmesh_assert_less_equal (send_list.size(),
solution->size());
431 parallel_object_only();
444 libmesh_assert_less_equal (send_list.size(),
solution->size());
456 if (subset !=
nullptr)
457 libmesh_not_implemented();
465 LOG_SCOPE(
"assemble()",
"System");
476 LOG_SCOPE(
"assemble_qoi()",
"System");
485 bool include_liftfunc,
486 bool apply_constraints)
489 LOG_SCOPE(
"assemble_qoi_derivative()",
"System");
503 if (qoi_indices.
size(*
this) > parameters.
size())
515 const Real threshold,
516 const bool verbose)
const
525 libMesh::out <<
" comparing matrices not supported." << std::endl;
530 const int name_result =
_sys_name.compare(other_system.
name());
533 if (name_result == 0)
547 if (solu_result == -1)
548 libMesh::out <<
" identical up to threshold." << std::endl;
550 libMesh::out <<
" first difference occurred at index = "
551 << solu_result <<
"." << std::endl;
557 std::vector<int> ov_result;
563 libMesh::out <<
" Fatal difference. This system handles "
564 << this->
n_vectors() <<
" add'l vectors," << std::endl
565 <<
" while the other system handles "
567 <<
" add'l vectors." << std::endl
568 <<
" Aborting comparison." << std::endl;
584 << pr.first <<
"\" ...";
590 ov_result.push_back(pr.second->compare (other_system_vector,
595 if (ov_result[ov_result.size()-1] == -1)
596 libMesh::out <<
" identical up to threshold." << std::endl;
598 libMesh::out <<
" first difference occurred at" << std::endl
599 <<
" index = " << ov_result[ov_result.size()-1] <<
"." << std::endl;
608 if ((name_result==0) && (solu_result==-1))
610 if (ov_result.size()==0)
611 overall_result =
true;
618 ov_identical = (ov_result[n]==-1);
621 while (ov_identical && n<ov_result.size());
622 overall_result = ov_identical;
626 overall_result =
false;
632 libMesh::out <<
"found no differences." << std::endl << std::endl;
634 libMesh::out <<
"found differences." << std::endl << std::endl;
637 return overall_result;
644 global_soln.resize (
solution->size());
654 global_soln.resize (
solution->size());
656 solution->localize_to_one (global_soln, dest_proc);
662 const bool projections,
671 _vectors.insert (std::make_pair (vec_name, buf));
684 #ifdef LIBMESH_ENABLE_GHOSTED
689 libmesh_error_msg(
"Cannot initialize ghosted vectors when they are not enabled.");
744 unsigned int num = 0;
745 while ((num<vec_num) && (v!=v_end))
761 unsigned int num = 0;
762 while ((num<vec_num) && (v!=v_end))
776 return *(libmesh_map_find(
_vectors, vec_name));
783 return *(libmesh_map_find(
_vectors, vec_name));
792 unsigned int num = 0;
793 while ((num<vec_num) && (v!=v_end))
808 unsigned int num = 0;
809 while ((num<vec_num) && (v!=v_end))
824 unsigned int num = 0;
825 while ((num<vec_num) && (v!=v_end))
839 for (; v != v_end; ++v)
842 if (&vec_reference == v->second)
878 libmesh_assert_greater_equal(qoi_num, -2);
896 std::ostringstream sensitivity_name;
897 sensitivity_name <<
"sensitivity_solution" << i;
899 return this->
add_vector(sensitivity_name.str());
906 std::ostringstream sensitivity_name;
907 sensitivity_name <<
"sensitivity_solution" << i;
909 return this->
get_vector(sensitivity_name.str());
916 std::ostringstream sensitivity_name;
917 sensitivity_name <<
"sensitivity_solution" << i;
919 return this->
get_vector(sensitivity_name.str());
926 return this->
add_vector(
"weighted_sensitivity_solution");
933 return this->
get_vector(
"weighted_sensitivity_solution");
940 return this->
get_vector(
"weighted_sensitivity_solution");
947 std::ostringstream adjoint_name;
948 adjoint_name <<
"adjoint_solution" << i;
959 std::ostringstream adjoint_name;
960 adjoint_name <<
"adjoint_solution" << i;
969 std::ostringstream adjoint_name;
970 adjoint_name <<
"adjoint_solution" << i;
979 std::ostringstream adjoint_name;
980 adjoint_name <<
"weighted_sensitivity_adjoint_solution" << i;
991 std::ostringstream adjoint_name;
992 adjoint_name <<
"weighted_sensitivity_adjoint_solution" << i;
1001 std::ostringstream adjoint_name;
1002 adjoint_name <<
"weighted_sensitivity_adjoint_solution" << i;
1011 std::ostringstream adjoint_rhs_name;
1012 adjoint_rhs_name <<
"adjoint_rhs" << i;
1014 return this->
add_vector(adjoint_rhs_name.str(),
false);
1021 std::ostringstream adjoint_rhs_name;
1022 adjoint_rhs_name <<
"adjoint_rhs" << i;
1024 return this->
get_vector(adjoint_rhs_name.str());
1031 std::ostringstream adjoint_rhs_name;
1032 adjoint_rhs_name <<
"adjoint_rhs" << i;
1034 return this->
get_vector(adjoint_rhs_name.str());
1041 std::ostringstream sensitivity_rhs_name;
1042 sensitivity_rhs_name <<
"sensitivity_rhs" << i;
1044 return this->
add_vector(sensitivity_rhs_name.str(),
false);
1051 std::ostringstream sensitivity_rhs_name;
1052 sensitivity_rhs_name <<
"sensitivity_rhs" << i;
1054 return this->
get_vector(sensitivity_rhs_name.str());
1061 std::ostringstream sensitivity_rhs_name;
1062 sensitivity_rhs_name <<
"sensitivity_rhs" << i;
1064 return this->
get_vector(sensitivity_rhs_name.str());
1071 const std::set<subdomain_id_type> *
const active_subdomains)
1083 libmesh_error_msg(
"ERROR: incompatible variable " << var <<
" has already been added for this system!");
1096 should_be_in_vg =
false;
1103 const std::set<subdomain_id_type> *
const
1104 their_active_subdomains (vg.implicitly_active() ?
1105 nullptr : &vg.active_subdomains());
1108 if (vg.type() != type)
1109 should_be_in_vg =
false;
1112 if (their_active_subdomains &&
1113 (!active_subdomains || (active_subdomains && active_subdomains->empty())))
1114 should_be_in_vg =
false;
1117 if (!their_active_subdomains && (active_subdomains && !active_subdomains->empty()))
1118 should_be_in_vg =
false;
1120 if (their_active_subdomains && active_subdomains)
1122 if (*their_active_subdomains != *active_subdomains)
1123 should_be_in_vg =
false;
1127 if (should_be_in_vg)
1129 const unsigned short curr_n_vars = cast_int<unsigned short>
1134 _variables.push_back(vg(vg.n_variables()-1));
1141 return this->
add_variables (std::vector<std::string>(1, var),
1151 const std::set<subdomain_id_type> *
const active_subdomains)
1162 const std::set<subdomain_id_type> *
const active_subdomains)
1168 for (
auto ovar : vars)
1175 libmesh_error_msg(
"ERROR: incompatible variable " << ovar <<
" has already been added for this system!");
1178 const unsigned short curr_n_vars = cast_int<unsigned short>
1181 const unsigned int next_first_component = this->
n_components();
1186 next_first_component, type) :
1188 next_first_component, type, *active_subdomains));
1200 libmesh_assert_equal_to ((curr_n_vars+vars.size()), this->
n_vars());
1208 return cast_int<unsigned int>(curr_n_vars+vars.size()-1);
1216 const std::set<subdomain_id_type> *
const active_subdomains)
1242 all_variable_numbers.resize(
n_vars());
1245 std::map<std::string, unsigned short int>::const_iterator
1247 std::map<std::string, unsigned short int>::const_iterator
1250 unsigned int count = 0;
1251 for ( ; it != it_end; ++it)
1253 all_variable_numbers[count] = it->second;
1260 std::set<dof_id_type> & var_indices)
const
1263 var_indices.clear();
1265 std::vector<dof_id_type> dof_indices;
1272 for (
const auto & elem : this->
get_mesh().active_local_element_ptr_range())
1278 if (first_local <= dof && dof < end_local)
1279 var_indices.insert(dof);
1288 for (
const auto & node : this->
get_mesh().local_node_ptr_range())
1292 for (
auto dof : dof_indices)
1293 if (first_local <= dof && dof < end_local)
1294 var_indices.insert(dof);
1301 unsigned int var_num)
const
1304 libmesh_assert_less (var_num, this->
n_vars());
1310 const unsigned int sys_num = this->
number();
1313 for (
const auto & node :
mesh.local_node_ptr_range())
1315 unsigned int n_comp = node->n_comp(sys_num,var_num);
1316 for (
unsigned int i=0; i<n_comp; i++)
1318 const dof_id_type index = node->dof_number(sys_num,var_num,i);
1324 for (
const auto & elem :
mesh.active_local_element_ptr_range())
1326 unsigned int n_comp = elem->n_comp(sys_num,var_num);
1327 for (
unsigned int i=0; i<n_comp; i++)
1329 const dof_id_type index = elem->dof_number(sys_num,var_num,i);
1341 std::set<dof_id_type> var_indices;
1351 libmesh_error_msg(
"Invalid norm_type = " << norm_type);
1359 std::set<unsigned int> * skip_dimensions)
const
1368 std::vector<FEMNormType> norms(this->
n_vars(),
L2);
1369 std::vector<Real> weights(this->
n_vars(), 0.0);
1370 norms[var] = norm_type;
1380 std::set<unsigned int> * skip_dimensions)
const
1383 parallel_object_only();
1385 LOG_SCOPE (
"calculate_norm()",
"System");
1390 if (
norm.is_discrete())
1394 unsigned int check_var = 0, check_end = this->
n_vars();
1395 for (; check_var != check_end; ++check_var)
1396 if ((
norm.weight(check_var) != 1.0) || (
norm.type(check_var) != norm_type0))
1400 if (check_var == this->
n_vars())
1409 libmesh_error_msg(
"Invalid norm_type0 = " << norm_type0);
1415 if (norm.
weight(var) == 0.0)
1433 bool using_hilbert_norm =
true,
1434 using_nonhilbert_norm =
true;
1440 Real norm_weight_sq =
norm.weight_sq(var);
1441 if (norm_weight_sq == 0.0)
1443 Real norm_weight =
norm.weight(var);
1447 if ((norm_type==
H1) ||
1453 if (!using_hilbert_norm)
1454 libmesh_not_implemented();
1455 using_nonhilbert_norm =
false;
1457 else if ((norm_type==
L1) ||
1458 (norm_type==
L_INF) ||
1462 if (!using_nonhilbert_norm)
1463 libmesh_not_implemented();
1464 using_hilbert_norm =
false;
1467 libmesh_not_implemented();
1472 std::vector<std::unique_ptr<FEBase>> fe_ptrs(4);
1473 std::vector<std::unique_ptr<QBase>> q_rules(4);
1478 for (
const auto &
dim : elem_dims)
1480 if (skip_dimensions && skip_dimensions->find(
dim) != skip_dimensions->end())
1488 fe_ptrs[
dim]->attach_quadrature_rule (q_rules[
dim].
get());
1491 std::vector<dof_id_type> dof_indices;
1494 for (
const auto & elem : this->
get_mesh().active_local_element_ptr_range())
1496 const unsigned int dim = elem->dim();
1498 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1503 if (elem->infinite() )
1504 libmesh_not_implemented();
1508 if (skip_dimensions && skip_dimensions->find(
dim) != skip_dimensions->end())
1512 QBase * qrule = q_rules[
dim].get();
1516 const std::vector<Real> & JxW = fe->
get_JxW();
1517 const std::vector<std::vector<Real>> * phi =
nullptr;
1518 if (norm_type ==
H1 ||
1525 const std::vector<std::vector<RealGradient>> * dphi =
nullptr;
1526 if (norm_type ==
H1 ||
1531 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1532 const std::vector<std::vector<RealTensor>> * d2phi =
nullptr;
1533 if (norm_type ==
H2 ||
1543 const unsigned int n_qp = qrule->
n_points();
1545 const unsigned int n_sf = cast_int<unsigned int>
1546 (dof_indices.size());
1549 for (
unsigned int qp=0; qp<n_qp; qp++)
1551 if (norm_type ==
L1)
1554 for (
unsigned int i=0; i != n_sf; ++i)
1555 u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]);
1556 v_norm += norm_weight *
1560 if (norm_type ==
L_INF)
1563 for (
unsigned int i=0; i != n_sf; ++i)
1564 u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]);
1565 v_norm = std::max(v_norm, norm_weight *
std::abs(u_h));
1568 if (norm_type ==
H1 ||
1573 for (
unsigned int i=0; i != n_sf; ++i)
1574 u_h += (*phi)[i][qp] * (*local_v)(dof_indices[i]);
1575 v_norm += norm_weight_sq *
1579 if (norm_type ==
H1 ||
1584 for (
unsigned int i=0; i != n_sf; ++i)
1585 grad_u_h.
add_scaled((*dphi)[i][qp], (*local_v)(dof_indices[i]));
1586 v_norm += norm_weight_sq *
1593 for (
unsigned int i=0; i != n_sf; ++i)
1594 grad_u_h.
add_scaled((*dphi)[i][qp], (*local_v)(dof_indices[i]));
1595 v_norm = std::max(v_norm, norm_weight * grad_u_h.
norm());
1598 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1599 if (norm_type ==
H2 ||
1603 for (
unsigned int i=0; i != n_sf; ++i)
1604 hess_u_h.
add_scaled((*d2phi)[i][qp], (*local_v)(dof_indices[i]));
1605 v_norm += norm_weight_sq *
1612 for (
unsigned int i=0; i != n_sf; ++i)
1613 hess_u_h.
add_scaled((*d2phi)[i][qp], (*local_v)(dof_indices[i]));
1614 v_norm = std::max(v_norm, norm_weight * hess_u_h.
norm());
1621 if (using_hilbert_norm)
1623 this->
comm().sum(v_norm);
1628 this->
comm().max(v_norm);
1638 std::ostringstream oss;
1641 const std::string & sys_name = this->
name();
1643 oss <<
" System #" << this->
number() <<
", \"" << sys_name <<
"\"\n"
1651 if (vg_description.
n_variables() > 1) oss <<
"{ ";
1653 oss <<
"\"" << vg_description.
name(vn) <<
"\" ";
1654 if (vg_description.
n_variables() > 1) oss <<
"} ";
1659 oss <<
" Finite Element Types=";
1660 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
1675 oss <<
'\n' <<
" Infinite Element Mapping=";
1684 oss <<
" Approximation Orders=";
1687 #ifndef LIBMESH_ENABLE_INFINITE_ELEMENTS
1702 oss <<
" n_dofs()=" << this->
n_dofs() <<
'\n';
1703 oss <<
" n_local_dofs()=" << this->
n_local_dofs() <<
'\n';
1704 #ifdef LIBMESH_ENABLE_CONSTRAINTS
1709 oss <<
" " <<
"n_vectors()=" << this->
n_vectors() <<
'\n';
1710 oss <<
" " <<
"n_matrices()=" << this->
n_matrices() <<
'\n';
1721 const std::string &
name))
1728 libMesh::out <<
"WARNING: Cannot specify both initialization function and object!"
1744 libMesh::out <<
"WARNING: Cannot specify both initialization object and function!"
1756 const std::string &
name))
1763 libMesh::out <<
"WARNING: Cannot specify both assembly function and object!"
1779 libMesh::out <<
"WARNING: Cannot specify both assembly object and function!"
1791 const std::string &
name))
1798 libMesh::out <<
"WARNING: Cannot specify both constraint function and object!"
1814 libMesh::out <<
"WARNING: Cannot specify both constraint object and function!"
1826 const std::string &,
1834 libMesh::out <<
"WARNING: Cannot specify both QOI function and object!"
1850 libMesh::out <<
"WARNING: Cannot specify both QOI object and function!"
1862 const QoISet &,
bool,
bool))
1869 libMesh::out <<
"WARNING: Cannot specify both QOI derivative function and object!"
1885 libMesh::out <<
"WARNING: Cannot specify both QOI derivative object and function!"
1953 bool include_liftfunc,
1954 bool apply_constraints)
1966 (qoi_indices, include_liftfunc, apply_constraints);
1973 const bool insist_on_success,
1978 parallel_object_only();
1996 std::unique_ptr<PointLocatorBase> locator_ptr =
mesh.sub_point_locator();
1999 if (!insist_on_success || !
mesh.is_serial())
2004 const std::set<subdomain_id_type> & raw_subdomains =
2006 const std::set<subdomain_id_type> * implicit_subdomains =
2007 raw_subdomains.empty() ? nullptr : &raw_subdomains;
2008 const Elem * e = locator(p, implicit_subdomains);
2012 if (e && this->
get_dof_map().is_evaluable(*e, var))
2019 this->
comm().min(lowest_owner);
2025 this->
comm().broadcast(u, lowest_owner);
2052 std::vector<dof_id_type> dof_indices;
2058 const unsigned int num_dofs = cast_int<unsigned int>
2059 (dof_indices.size());
2074 for (
unsigned int l=0; l<num_dofs; l++)
2076 u += fe_data.shape[l] * (*sol)(dof_indices[l]);
2102 const bool insist_on_success,
2107 parallel_object_only();
2125 std::unique_ptr<PointLocatorBase> locator_ptr =
mesh.sub_point_locator();
2128 if (!insist_on_success || !
mesh.is_serial())
2133 const std::set<subdomain_id_type> & raw_subdomains =
2135 const std::set<subdomain_id_type> * implicit_subdomains =
2136 raw_subdomains.empty() ? nullptr : &raw_subdomains;
2137 const Elem * e = locator(p, implicit_subdomains);
2141 if (e && this->
get_dof_map().is_evaluable(*e, var))
2148 this->
comm().min(lowest_owner);
2154 this->
comm().broadcast(grad_u, lowest_owner);
2179 const unsigned int dim = e.
dim();
2185 std::vector<dof_id_type> dof_indices;
2191 const unsigned int num_dofs = cast_int<unsigned int>
2192 (dof_indices.size());
2208 for (
unsigned int l=0; l<num_dofs; l++)
2212 for (std::size_t v=0; v<
dim; v++)
2213 for (std::size_t xyz=0; xyz<LIBMESH_DIM; xyz++)
2218 * (*sol)(dof_indices[l]);
2243 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
2246 const bool insist_on_success,
2251 parallel_object_only();
2269 std::unique_ptr<PointLocatorBase> locator_ptr =
mesh.sub_point_locator();
2272 if (!insist_on_success || !
mesh.is_serial())
2277 const std::set<subdomain_id_type> & raw_subdomains =
2279 const std::set<subdomain_id_type> * implicit_subdomains =
2280 raw_subdomains.empty() ? nullptr : &raw_subdomains;
2281 const Elem * e = locator(p, implicit_subdomains);
2285 if (e && this->
get_dof_map().is_evaluable(*e, var))
2292 this->
comm().min(lowest_owner);
2298 this->
comm().broadcast(hess_u, lowest_owner);
2318 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2320 libmesh_not_implemented();
2330 std::vector<dof_id_type> dof_indices;
2336 const unsigned int num_dofs = cast_int<unsigned int>
2337 (dof_indices.size());
2349 const std::vector<std::vector<RealTensor>> & d2phi = fe->get_d2phi();
2352 fe->reinit (&e, &coor);
2357 for (
unsigned int l=0; l<num_dofs; l++)
2359 hess_u.
add_scaled (d2phi[l][0], (*sol)(dof_indices[l]));
2385 libmesh_error_msg(
"We can only accumulate a hessian with --enable-second");
2392 const NumericVector<Number> *)
const
2394 libmesh_error_msg(
"We can only accumulate a hessian with --enable-second");
2402 libmesh_error_msg(
"We can only accumulate a hessian with --enable-second");
2410 libmesh_error_msg(
"We can only accumulate a hessian with --enable-second");
2416 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
void set_vector_preservation(const std::string &vec_name, bool preserve)
Allows one to set the boolean controlling whether the vector identified by vec_name should be "preser...
The QBase class provides the basic functionality from which various quadrature rules can be derived.
Manages consistently variables, degrees of freedom, and coefficient vectors.
std::map< std::string, bool > _vector_projections
Holds true if a vector by that name should be projected onto a changed grid, false if it should be ze...
unsigned int n_vars() const
bool have_vector(const std::string &vec_name) const
bool is_evaluable(const DofObjectSubclass &obj, unsigned int var_num=libMesh::invalid_uint) const
const EquationSystems & get_equation_systems() const
void project_vector(NumericVector< Number > &new_vector, FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr, int is_adjoint=-1) const
Projects arbitrary functions onto a vector of degree of freedom values for the current system.
virtual void reinit()
Reinitializes degrees of freedom and other required data on the current mesh.
NumericVector< Number > & get_sensitivity_rhs(unsigned int i=0)
bool _is_initialized
true when additional vectors and variables do not require immediate initialization,...
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
void prepare_send_list()
Takes the _send_list vector (which may have duplicate entries) and sorts it.
std::vector< VariableGroup > _variable_groups
The VariableGroup in this System.
std::map< std::string, NumericVector< Number > * >::iterator vectors_iterator
Vector iterator typedefs.
void(* _constrain_system_function)(EquationSystems &es, const std::string &name)
Function to impose constraints.
void attach_QOI_function(void fptr(EquationSystems &es, const std::string &name, const QoISet &qoi_indices))
Register a user function for evaluating the quantities of interest, whose values should be placed in ...
FEFamily family
The type of finite element.
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.
NumericVector< Number > & add_adjoint_solution(unsigned int i=0)
virtual bool compare(const System &other_system, const Real threshold, const bool verbose) const
std::unique_ptr< DofMap > _dof_map
Data structure describing the relationship between nodes, variables, etc...
virtual bool contains_point(const Point &p, Real tol=TOLERANCE) const
unsigned int add_variables(const std::vector< std::string > &vars, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
void init()
Initializes degrees of freedom on the current mesh.
void(* _assemble_system_function)(EquationSystems &es, const std::string &name)
Function that assembles the system.
std::string get_info() const
Gets summary info about the sparsity bandwidth and constraints.
virtual void user_assembly()
Calls user's attached assembly function, or is overridden by the user in derived classes.
unsigned int n_points() const
NumericVector< Number > & get_adjoint_solution(unsigned int i=0)
virtual std::string system_type() const
std::map< std::string, NumericVector< Number > * > _vectors
Some systems need an arbitrary number of vectors.
NumericVector< Number > & get_weighted_sensitivity_adjoint_solution(unsigned int i=0)
std::map< std::string, ParallelType > _vector_types
Holds the type of a vector.
This class implements reference counting.
const std::string & variable_name(const unsigned int i) const
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned short dim() const =0
unsigned int n_vectors() const
NumericVector< Number > & add_sensitivity_solution(unsigned int i=0)
virtual void restrict_vectors()
Restrict vectors after the mesh has coarsened.
This class forms the foundation from which generic finite elements may be derived.
bool identify_variable_groups() const
const std::string & name(unsigned int v) const
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.
virtual bool infinite() const =0
const std::vector< Real > & get_JxW() const
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
const Parallel::Communicator & comm() const
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
NumericVector< Number > & get_adjoint_rhs(unsigned int i=0)
Data structure for holding completed parameter sensitivity calculations.
void process_constraints(MeshBase &)
Postprocesses any constrained degrees of freedom to be constrained only in terms of unconstrained dof...
unsigned int number() const
const FEType & type() const
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO)
After calling this method, any solve will be restricted to the given subdomain.
Real weight(unsigned int var) const
virtual void re_update()
Re-update the local values when the mesh has changed.
dof_id_type first_dof(const processor_id_type proc) const
const std::vector< std::vector< OutputGradient > > & get_dphi() const
const Variable & variable(unsigned int var) const
Return a constant reference to Variable var.
Abstract base class to be used for system initialization.
const std::set< subdomain_id_type > & active_subdomains() const
System & operator=(const System &)
This isn't a copyable object, so let's make sure nobody tries.
virtual ~System()
Destructor.
void local_dof_indices(const unsigned int var, std::set< dof_id_type > &var_indices) const
Fills the std::set with the degrees of freedom on the local processor corresponding the the variable ...
virtual unsigned int n_matrices() const
processor_id_type processor_id() const
virtual void qoi(const QoISet &qoi_indices)=0
Quantity of interest function.
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
void(* _init_system_function)(EquationSystems &es, const std::string &name)
Function that initializes the system.
void attach_constraint_object(Constraint &constrain)
Register a user object for imposing constraints.
void(* _qoi_evaluate_function)(EquationSystems &es, const std::string &name, const QoISet &qoi_indices)
Function to evaluate quantity of interest.
NumericVector< Number > & get_weighted_sensitivity_solution()
const std::string & vector_name(const unsigned int vec_num) const
virtual void assemble()
Prepares matrix and _dof_map for matrix assembly.
virtual numeric_index_type size() const =0
int vector_is_adjoint(const std::string &vec_name) const
virtual void qoi_derivative(const QoISet &qoi_indices, bool include_liftfunc, bool apply_constraints)=0
Quantity of interest derivative function.
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
Number current_solution(const dof_id_type global_dof_number) const
unsigned int n_components() const
bool _basic_system_only
Holds true if the components of more advanced system types (e.g.
void attach_QOI_derivative(void fptr(EquationSystems &es, const std::string &name, const QoISet &qoi_indices, bool include_liftfunc, bool apply_constraints))
Register a user function for evaluating derivatives of a quantity of interest with respect to test fu...
void zero_variable(NumericVector< Number > &v, unsigned int var_num) const
Zeroes all dofs in v that correspond to variable number var_num.
virtual void enable_out_of_mesh_mode()=0
Enables out-of-mesh mode.
NumericVector< Number > & add_adjoint_rhs(unsigned int i=0)
Initialization * _init_system_object
Object that initializes the system.
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
bool has_variable(const std::string &var) const
unsigned int n_variable_groups() const
This is the MeshBase class.
Abstract base class to be used for derivatives of quantities of interest.
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
void attach_init_object(Initialization &init)
Register a user class to use to initialize the system.
virtual void reinit_constraints()
Reinitializes the constraints for this system.
vectors_iterator vectors_end()
End of vectors container.
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
virtual Real subset_l2_norm(const std::set< numeric_index_type > &indices) const
vectors_iterator vectors_begin()
Beginning of vectors container.
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
auto norm_sq() const -> decltype(std::norm(T()))
MeshBase & _mesh
Constant reference to the mesh data structure used for the simulation.
virtual void user_QOI(const QoISet &qoi_indices)
Calls user's attached quantity of interest function, or is overridden by the user in derived classes.
std::vector< Variable > _variables
The Variable in this System.
processor_id_type n_processors() const
dof_id_type n_local_dofs() const
void(* _qoi_evaluate_derivative_function)(EquationSystems &es, const std::string &name, const QoISet &qoi_indices, bool include_liftfunc, bool apply_constraints)
Function to evaluate quantity of interest derivative.
const VariableGroup & variable_group(unsigned int vg) const
Return a constant reference to VariableGroup vg.
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
auto norm_sq() const -> decltype(std::norm(T()))
void attach_QOI_object(QOI &qoi)
Register a user object for evaluating the quantities of interest, whose values should be placed in Sy...
std::string get_info() const
virtual void prolong_vectors()
Prolong vectors after the mesh has refined.
void set_vector_as_adjoint(const std::string &vec_name, int qoi_num)
Allows one to set the QoI index controlling whether the vector identified by vec_name represents a so...
const MeshBase & get_mesh() const
processor_id_type processor_id() const
dof_id_type n_local_constrained_dofs() const
void update_global_solution(std::vector< Number > &global_soln) const
Fill the input vector global_soln so that it contains the global solution on all processors.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
virtual void constrain()=0
Constraint function.
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
NumericVector< Number > & get_sensitivity_solution(unsigned int i=0)
Abstract base class to be used for quantities of interest.
A Point defines a location in LIBMESH_DIM dimensional Real space.
System(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
uint8_t processor_id_type
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
const std::set< unsigned char > & elem_dimensions() const
void add_scaled(const TypeTensor< T2 > &, const T &)
Add a scaled tensor to this tensor without creating a temporary.
const FEType & variable_type(const unsigned int c) const
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
virtual void assemble_qoi(const QoISet &qoi_indices=QoISet())
Calls user qoi function.
This is a base class for classes which represent subsets of the dofs of a System.
virtual Real l1_norm() const =0
const std::vector< dof_id_type > & get_send_list() const
void attach_QOI_derivative_object(QOIDerivative &qoi_derivative)
Register a user object for evaluating derivatives of a quantity of interest with respect to test func...
virtual void qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &sensitivities)
Solves for the derivative of each of the system's quantities of interest q in qoi[qoi_indices] with r...
Abstract base class to be used for system assembly.
virtual void assemble()=0
Assembly function.
This is the EquationSystems class.
Assembly * _assemble_system_object
Object that assembles the system.
void create_dof_constraints(const MeshBase &, Real time=0)
Rebuilds the raw degree of freedom and DofObject constraints.
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Abstract base class to be used for system constraints.
virtual Real linfty_norm() const =0
const FEType & variable_type(const unsigned int i) const
Real discrete_var_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type) const
Finds the discrete norm for the entries in the vector corresponding to Dofs associated with var.
Constraint * _constrain_system_object
Object that constrains the system.
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
const Elem & get(const ElemType type_in)
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Number point_value(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
bool _is_initialized
Flag that tells if init() has been called.
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
QOIDerivative * _qoi_evaluate_derivative_object
Object to compute derivatives of quantities of interest.
This class handles the numbering of degrees of freedom on a mesh.
const std::string & name() const
virtual void init_data()
Initializes the data for the system.
void attach_constraint_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function for imposing constraints.
const NumericVector< Number > * request_vector(const std::string &vec_name) const
dof_id_type n_constrained_dofs() const
void attach_init_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in initializing the system.
This is the base class from which all geometric element types are derived.
unsigned short int variable_number(const std::string &var) const
virtual Real l2_norm() const =0
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
virtual void assemble_qoi_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true)
Calls user qoi derivative function.
EquationSystems & _equation_systems
Constant reference to the EquationSystems object used for the simulation.
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
virtual void clear()
Clear all the data structures associated with the system.
virtual void forward_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &sensitivities)
Solves for parameter sensitivities using the forward method.
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
ParallelType
Defines an enum for parallel data structure types.
const DofMap & get_dof_map() const
dof_id_type n_dofs() const
virtual void user_initialization()
Calls user's attached initialization function, or is overridden by the user in derived classes.
void attach_assemble_object(Assembly &assemble)
Register a user object to use in assembling the system matrix and RHS.
QOI * _qoi_evaluate_object
Object to compute quantities of interest.
const std::vector< std::vector< OutputShape > > & get_phi() const
unsigned int n_variables() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
This is at the core of this class.
An object whose state is distributed along a set of processors.
virtual void user_constrain()
Calls user's attached constraint function, or is overridden by the user in derived classes.
NumericVector< Number > & add_weighted_sensitivity_solution()
void get_all_variable_numbers(std::vector< unsigned int > &all_variable_numbers) const
Fills all_variable_numbers with all the variable numbers for the variables that have been added to th...
bool vector_preservation(const std::string &vec_name) const
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
virtual void initialize()=0
Initialization function.
auto norm() const -> decltype(std::norm(T()))
NumericVector< Number > & add_weighted_sensitivity_adjoint_solution(unsigned int i=0)
void remove_vector(const std::string &vec_name)
Removes the additional vector vec_name from this system.
This class defines a logically grouped set of variables in the system.
virtual void user_QOI_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true)
Calls user's attached quantity of interest derivative function, or is overridden by the user in deriv...
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
This is the base class for point locators.
const NumericVector< Number > & get_vector(const std::string &vec_name) const
std::map< std::string, unsigned short int > _variable_numbers
The variable numbers corresponding to user-specified names, useful for name-based lookups.
auto norm() const -> decltype(std::norm(T()))
virtual void update()
Update the local values to reflect the solution on neighboring processors.
NumericVector< Number > & add_sensitivity_rhs(unsigned int i=0)
bool _solution_projection
Holds true if the solution vector should be projected onto a changed grid, false if it should be zero...
std::size_t size(const System &sys) const
Gradient point_gradient(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Tensor point_hessian(unsigned int var, const Point &p, const bool insist_on_success=true, const NumericVector< Number > *sol=nullptr) const
const std::string _sys_name
A name associated with this system.
std::map< std::string, NumericVector< Number > * >::const_iterator const_vectors_iterator
virtual void adjoint_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &sensitivities)
Solves for parameter sensitivities using the adjoint method.
const VariableGroup & variable_group(const unsigned int c) const
bool closed()
Checks that the library has been closed.
void enable_derivative()
Enable the computation of shape gradients (dshape).
std::vector< std::vector< Real > > local_transform
Storage for local to global mapping at p.
dof_id_type end_dof(const processor_id_type proc) const
static void compute_data(const unsigned int dim, const FEType &fe_t, const Elem *elem, FEComputeData &data)
Lets the appropriate child of FEBase compute the requested data for the input specified in data,...
std::vector< Gradient > dshape
Storage for the computed shape derivative values.
std::map< std::string, int > _vector_is_adjoint
Holds non-negative if a vector by that name should be projected using adjoint constraints/BCs,...