OGS
HydroMechanicsProcess.cpp
Go to the documentation of this file.
1
12
20#include "MeshLib/Mesh.h"
21#include "MeshLib/Properties.h"
30
31namespace ProcessLib
32{
33namespace LIE
34{
35namespace HydroMechanics
36{
37template <int GlobalDim>
39 std::string name,
40 MeshLib::Mesh& mesh,
41 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
42 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
43 unsigned const integration_order,
44 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
45 process_variables,
47 SecondaryVariableCollection&& secondary_variables,
48 bool const use_monolithic_scheme)
49 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
50 integration_order, std::move(process_variables),
51 std::move(secondary_variables), use_monolithic_scheme),
52 _process_data(std::move(process_data))
53{
54 INFO("[LIE/HM] looking for fracture elements in the given mesh");
55 std::vector<int> vec_fracture_mat_IDs;
56 std::vector<std::vector<MeshLib::Element*>> vec_vec_fracture_elements;
57 std::vector<std::vector<MeshLib::Element*>>
58 vec_vec_fracture_matrix_elements;
59 std::vector<std::vector<MeshLib::Node*>> vec_vec_fracture_nodes;
60 std::vector<std::pair<std::size_t, std::vector<int>>>
61 vec_branch_nodeID_matIDs;
62 std::vector<std::pair<std::size_t, std::vector<int>>>
63 vec_junction_nodeID_matIDs;
65 mesh, _vec_matrix_elements, vec_fracture_mat_IDs,
66 vec_vec_fracture_elements, vec_vec_fracture_matrix_elements,
67 vec_vec_fracture_nodes, vec_branch_nodeID_matIDs,
68 vec_junction_nodeID_matIDs);
70 vec_vec_fracture_elements[0].begin(),
71 vec_vec_fracture_elements[0].end());
74 vec_vec_fracture_matrix_elements[0].begin(),
75 vec_vec_fracture_matrix_elements[0].end());
77 vec_vec_fracture_nodes[0].begin(),
78 vec_vec_fracture_nodes[0].end());
79
80 if (!_vec_fracture_elements.empty())
81 {
82 // set fracture property assuming a fracture forms a straight line
83 setFractureProperty(GlobalDim,
85 *_process_data.fracture_property);
86 }
87
88 //
89 // If Neumann BCs for the displacement_jump variable are required they need
90 // special treatment because of the levelset function. The implementation
91 // exists in the version 6.1.0 (e54815cc07ee89c81f953a4955b1c788595dd725)
92 // and was removed due to lack of applications.
93 //
94
95 if (!_process_data.deactivate_matrix_in_flow)
96 {
97 _process_data.p_element_status =
98 std::make_unique<MeshLib::ElementStatus>(&mesh);
99 }
100 else
101 {
102 auto const range =
104 if (!range)
105 {
106 OGS_FATAL(
107 "Could not get minimum/maximum ranges values for the "
108 "MaterialIDs property in the mesh '{:s}'.",
109 mesh.getName());
110 }
111
112 std::vector<int> vec_p_inactive_matIDs;
113 for (int matID = range->first; matID <= range->second; matID++)
114 {
115 if (std::find(vec_fracture_mat_IDs.begin(),
116 vec_fracture_mat_IDs.end(),
117 matID) == vec_fracture_mat_IDs.end())
118 {
119 vec_p_inactive_matIDs.push_back(matID);
120 }
121 }
122 _process_data.p_element_status =
123 std::make_unique<MeshLib::ElementStatus>(&mesh,
124 vec_p_inactive_matIDs);
125
126 const int monolithic_process_id = 0;
127 ProcessVariable const& pv_p =
128 getProcessVariables(monolithic_process_id)[0];
130 }
131}
132
133template <int GlobalDim>
135{
136 //------------------------------------------------------------
137 // prepare mesh subsets to define DoFs
138 //------------------------------------------------------------
139 // for extrapolation
140 _mesh_subset_all_nodes =
141 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
142 // pressure
143 _mesh_nodes_p = MeshLib::getBaseNodes(
144 _process_data.p_element_status->getActiveElements());
145 _mesh_subset_nodes_p =
146 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh_nodes_p);
147 // regular u
148 _mesh_subset_matrix_nodes =
149 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
150 if (!_vec_fracture_nodes.empty())
151 {
152 // u jump
153 _mesh_subset_fracture_nodes =
154 std::make_unique<MeshLib::MeshSubset>(_mesh, _vec_fracture_nodes);
155 }
156
157 // Collect the mesh subsets in a vector.
158 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
159 std::vector<int> vec_n_components;
160 std::vector<std::vector<MeshLib::Element*> const*> vec_var_elements;
161 // pressure
162 vec_n_components.push_back(1);
163 all_mesh_subsets.emplace_back(*_mesh_subset_nodes_p);
164 if (!_process_data.deactivate_matrix_in_flow)
165 {
166 vec_var_elements.push_back(&_mesh.getElements());
167 }
168 else
169 {
170 // TODO set elements including active nodes for pressure.
171 // cannot use ElementStatus
172 vec_var_elements.push_back(&_vec_fracture_matrix_elements);
173 }
174 // regular displacement
175 vec_n_components.push_back(GlobalDim);
176 std::generate_n(std::back_inserter(all_mesh_subsets), GlobalDim,
177 [&]() { return *_mesh_subset_matrix_nodes; });
178 vec_var_elements.push_back(&_vec_matrix_elements);
179 if (!_vec_fracture_nodes.empty())
180 {
181 // displacement jump
182 vec_n_components.push_back(GlobalDim);
183 std::generate_n(std::back_inserter(all_mesh_subsets), GlobalDim,
184 [&]() { return *_mesh_subset_fracture_nodes; });
185 vec_var_elements.push_back(&_vec_fracture_matrix_elements);
186 }
187
188 INFO("[LIE/HM] creating a DoF table");
189 _local_to_global_index_map =
190 std::make_unique<NumLib::LocalToGlobalIndexMap>(
191 std::move(all_mesh_subsets),
192 vec_n_components,
193 vec_var_elements,
195
196 DBUG("created {:d} DoF", _local_to_global_index_map->size());
197}
198
199template <int GlobalDim>
201 NumLib::LocalToGlobalIndexMap const& dof_table,
202 MeshLib::Mesh const& mesh,
203 unsigned const integration_order)
204{
205 assert(mesh.getDimension() == GlobalDim);
206 INFO("[LIE/HM] creating local assemblers");
211 mesh.getElements(), dof_table, _local_assemblers,
212 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
213 _process_data);
214
215 auto add_secondary_variable = [&](std::string const& name,
216 int const num_components,
217 auto get_ip_values_function)
218 {
219 _secondary_variables.addSecondaryVariable(
220 name,
221 makeExtrapolator(num_components, getExtrapolator(),
222 _local_assemblers,
223 std::move(get_ip_values_function)));
224 };
225
226 add_secondary_variable(
227 "sigma",
229 &LocalAssemblerInterface::getIntPtSigma);
230
231 add_secondary_variable(
232 "epsilon",
234 &LocalAssemblerInterface::getIntPtEpsilon);
235
236 add_secondary_variable("velocity", GlobalDim,
237 &LocalAssemblerInterface::getIntPtDarcyVelocity);
238
239 add_secondary_variable("fracture_velocity", GlobalDim,
240 &LocalAssemblerInterface::getIntPtFractureVelocity);
241
242 add_secondary_variable("fracture_stress", GlobalDim,
243 &LocalAssemblerInterface::getIntPtFractureStress);
244
245 add_secondary_variable("fracture_aperture", 1,
246 &LocalAssemblerInterface::getIntPtFractureAperture);
247
248 add_secondary_variable(
249 "fracture_permeability", 1,
250 &LocalAssemblerInterface::getIntPtFracturePermeability);
251
252 _process_data.element_stresses = MeshLib::getOrCreateMeshProperty<double>(
253 const_cast<MeshLib::Mesh&>(mesh), "sigma_avg",
256
257 _process_data.element_velocities = MeshLib::getOrCreateMeshProperty<double>(
258 const_cast<MeshLib::Mesh&>(mesh), "velocity_avg",
259 MeshLib::MeshItemType::Cell, GlobalDim);
260
261 if (!_vec_fracture_elements.empty())
262 {
263 auto mesh_prop_levelset = MeshLib::getOrCreateMeshProperty<double>(
264 const_cast<MeshLib::Mesh&>(mesh), "levelset1",
266 mesh_prop_levelset->resize(mesh.getNumberOfElements());
267 for (MeshLib::Element const* e : _mesh.getElements())
268 {
269 if (e->getDimension() < GlobalDim)
270 {
271 continue;
272 }
273
274 std::vector<FractureProperty*> fracture_props(
275 {_process_data.fracture_property.get()});
276 std::vector<JunctionProperty*> junction_props;
277 std::unordered_map<int, int> fracID_to_local({{0, 0}});
278 std::vector<double> levelsets = uGlobalEnrichments(
279 fracture_props, junction_props, fracID_to_local,
280 Eigen::Vector3d(MeshLib::getCenterOfGravity(*e).data()));
281 (*mesh_prop_levelset)[e->getID()] = levelsets[0];
282 }
283
284 _process_data.element_local_jumps =
285 MeshLib::getOrCreateMeshProperty<double>(
286 const_cast<MeshLib::Mesh&>(mesh), "local_jump_w_avg",
287 MeshLib::MeshItemType::Cell, GlobalDim);
288
289 _process_data.element_fracture_stresses =
290 MeshLib::getOrCreateMeshProperty<double>(
291 const_cast<MeshLib::Mesh&>(mesh), "fracture_stress_avg",
292 MeshLib::MeshItemType::Cell, GlobalDim);
293
294 _process_data.element_fracture_velocities =
295 MeshLib::getOrCreateMeshProperty<double>(
296 const_cast<MeshLib::Mesh&>(mesh), "fracture_velocity_avg",
297 MeshLib::MeshItemType::Cell, GlobalDim);
298
299 auto mesh_prop_b = MeshLib::getOrCreateMeshProperty<double>(
300 const_cast<MeshLib::Mesh&>(mesh), "fracture_aperture_avg",
302 mesh_prop_b->resize(mesh.getNumberOfElements());
303
304 auto const mesh_prop_matid = materialIDs(mesh);
305 if (!mesh_prop_matid)
306 {
307 OGS_FATAL("Could not access MaterialIDs property from mesh.");
308 }
309 auto const& frac = _process_data.fracture_property;
310 for (MeshLib::Element const* e : _mesh.getElements())
311 {
312 if (e->getDimension() == GlobalDim)
313 {
314 continue;
315 }
316 if ((*mesh_prop_matid)[e->getID()] != frac->mat_id)
317 {
318 continue;
319 }
320 // Mean value for the element. This allows usage of node based
321 // properties for aperture.
322 (*mesh_prop_b)[e->getID()] =
323 frac->aperture0
324 .getNodalValuesOnElement(*e, /*time independent*/ 0)
325 .mean();
326 }
327 _process_data.mesh_prop_b = mesh_prop_b;
328
329 auto mesh_prop_k_f = MeshLib::getOrCreateMeshProperty<double>(
330 const_cast<MeshLib::Mesh&>(mesh), "fracture_permeability_avg",
332 mesh_prop_k_f->resize(mesh.getNumberOfElements());
333 _process_data.mesh_prop_k_f = mesh_prop_k_f;
334
335 auto mesh_prop_fracture_shear_failure =
336 MeshLib::getOrCreateMeshProperty<double>(
337 const_cast<MeshLib::Mesh&>(mesh), "f_shear_failure",
339 mesh_prop_fracture_shear_failure->resize(mesh.getNumberOfElements());
340 _process_data.mesh_prop_fracture_shear_failure =
341 mesh_prop_fracture_shear_failure;
342
343 auto mesh_prop_nodal_p = MeshLib::getOrCreateMeshProperty<double>(
344 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
346 mesh_prop_nodal_p->resize(mesh.getNumberOfNodes());
347 _process_data.mesh_prop_nodal_p = mesh_prop_nodal_p;
348
349 _process_data.mesh_prop_nodal_forces =
350 MeshLib::getOrCreateMeshProperty<double>(
351 const_cast<MeshLib::Mesh&>(mesh), "NodalForces",
352 MeshLib::MeshItemType::Node, GlobalDim);
353 assert(_process_data.mesh_prop_nodal_forces->size() ==
354 GlobalDim * mesh.getNumberOfNodes());
355
356 _process_data.mesh_prop_nodal_forces_jump =
357 MeshLib::getOrCreateMeshProperty<double>(
358 const_cast<MeshLib::Mesh&>(mesh), "NodalForcesJump",
359 MeshLib::MeshItemType::Node, GlobalDim);
360 assert(_process_data.mesh_prop_nodal_forces_jump->size() ==
361 GlobalDim * mesh.getNumberOfNodes());
362
363 _process_data.mesh_prop_hydraulic_flow =
364 MeshLib::getOrCreateMeshProperty<double>(
365 const_cast<MeshLib::Mesh&>(mesh), "MassFlowRate",
367 assert(_process_data.mesh_prop_hydraulic_flow->size() ==
368 mesh.getNumberOfNodes());
369 }
370}
371
372template <int GlobalDim>
374 std::vector<GlobalVector*> const& x,
375 std::vector<GlobalVector*> const& x_prev, const double t, double const dt,
376 int const process_id)
377{
378 if (process_id == 0)
379 {
380 DBUG("PostTimestep HydroMechanicsProcess.");
381
383 getProcessVariables(process_id)[0];
386 _local_assemblers, pv.getActiveElementIDs(), getDOFTables(x.size()),
387 x, x_prev, t, dt, process_id);
388 }
389
390 DBUG("Compute the secondary variables for HydroMechanicsProcess.");
391
392 const auto& dof_table = getDOFTable(process_id);
393
394 // Copy displacement jumps in a solution vector to mesh property
395 // Remark: the copy is required because mesh properties for primary
396 // variables are set during output and are not ready yet when this function
397 // is called.
398 int g_variable_id = 0;
399 {
400 const int monolithic_process_id = 0;
401 auto const& pvs = getProcessVariables(monolithic_process_id);
402 auto const it =
403 std::find_if(pvs.begin(), pvs.end(),
404 [](ProcessVariable const& pv)
405 { return pv.getName() == "displacement_jump1"; });
406 if (it == pvs.end())
407 {
408 OGS_FATAL(
409 "Didn't find expected 'displacement_jump1' process variable.");
410 }
411 g_variable_id = static_cast<int>(std::distance(pvs.begin(), it));
412 }
413
415
416 const int monolithic_process_id = 0;
417 ProcessVariable& pv_g =
418 this->getProcessVariables(monolithic_process_id)[g_variable_id];
419 auto const num_comp = pv_g.getNumberOfGlobalComponents();
420 auto& mesh_prop_g = *MeshLib::getOrCreateMeshProperty<double>(
421 _mesh, pv_g.getName(), MeshLib::MeshItemType::Node, num_comp);
422 for (int component_id = 0; component_id < num_comp; ++component_id)
423 {
424 auto const& mesh_subset =
425 dof_table.getMeshSubset(g_variable_id, component_id);
426 auto const mesh_id = mesh_subset.getMeshID();
427 for (auto const* node : mesh_subset.getNodes())
428 {
430 node->getID());
431
432 auto const global_index =
433 dof_table.getGlobalIndex(l, g_variable_id, component_id);
434 mesh_prop_g[node->getID() * num_comp + component_id] =
435 (*x[process_id])[global_index];
436 }
437 }
438}
439
440template <int GlobalDim>
442{
443 return false;
444}
445
446template <int GlobalDim>
448 const double t, double const dt, std::vector<GlobalVector*> const& x,
449 std::vector<GlobalVector*> const& x_prev, int const process_id,
451{
452 DBUG("Assemble HydroMechanicsProcess.");
453
454 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
455 _local_to_global_index_map.get()};
456 // Call global assembler for each local assembly item.
458 _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
459 dof_table, t, dt, x, x_prev, process_id, M, K, b);
460}
461
462template <int GlobalDim>
464 const double t, double const dt, std::vector<GlobalVector*> const& x,
465 std::vector<GlobalVector*> const& x_prev, int const process_id,
467{
468 DBUG("AssembleWithJacobian HydroMechanicsProcess.");
469
470 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
471
472 // Call global assembler for each local assembly item.
473 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
474 _local_to_global_index_map.get()};
477 _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x,
478 x_prev, process_id, M, K, b, Jac);
479
480 auto copyRhs = [&](int const variable_id, auto& output_vector)
481 {
482 transformVariableFromGlobalVector(b, variable_id,
483 *_local_to_global_index_map,
484 output_vector, std::negate<double>());
485 };
486 copyRhs(0, *_process_data.mesh_prop_hydraulic_flow);
487 copyRhs(1, *_process_data.mesh_prop_nodal_forces);
488 copyRhs(2, *_process_data.mesh_prop_nodal_forces_jump);
489}
490
491template <int GlobalDim>
493 std::vector<GlobalVector*> const& x, double const t, double const dt,
494 int const process_id)
495{
496 DBUG("PreTimestep HydroMechanicsProcess.");
497
498 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
499
502 pv.getActiveElementIDs(), *_local_to_global_index_map, *x[process_id],
503 t, dt);
504}
505
506// ------------------------------------------------------------------------------------
507// template instantiation
508// ------------------------------------------------------------------------------------
509template class HydroMechanicsProcess<2>;
510template class HydroMechanicsProcess<3>;
511
512} // namespace HydroMechanics
513} // namespace LIE
514} // namespace ProcessLib
Definition of the ElementStatus class.
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of the MeshInformation class.
Definition of the class Properties that implements a container of properties.
Definition of the Mesh class.
Global vector based on Eigen vector.
Definition EigenVector.h:25
const double * data() const
Definition Point3d.h:59
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
const std::string getName() const
Get name of the mesh.
Definition Mesh.h:103
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:97
static std::optional< std::pair< T, T > > const getValueBounds(MeshLib::PropertyVector< T > const &property)
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, int const process_id) override
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
std::vector< MeshLib::Element * > _vec_fracture_matrix_elements
HydroMechanicsProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, HydroMechanicsProcessData< GlobalDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) override
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
std::vector< std::size_t > const & getActiveElementIDs() const
ParameterLib::Parameter< double > const & getInitialCondition() const
std::string const & getName() const
int getNumberOfGlobalComponents() const
Returns the number of components of the process variable.
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
Handles configuration of several secondary variables from the project file.
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
MathLib::Point3d getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition Element.cpp:124
@ BY_COMPONENT
Ordering data by component type.
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, NumLib::IntegrationOrder const integration_order, ExtraCtorArgs &&... extra_ctor_args)
void setFractureProperty(int const dim, MeshLib::Element const &e, FractureProperty &frac_prop)
void getFractureMatrixDataInMesh(MeshLib::Mesh const &mesh, std::vector< MeshLib::Element * > &vec_matrix_elements, std::vector< int > &vec_fracture_mat_IDs, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_elements, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_matrix_elements, std::vector< std::vector< MeshLib::Node * > > &vec_fracture_nodes, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_branch_nodeID_matIDs, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_junction_nodeID_matIDs)
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberDereferenced(Object &object, Method method, Container const &container, Args &&... args)