OGS
RichardsMechanicsProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
25
26namespace ProcessLib
27{
28namespace RichardsMechanics
29{
30template <int DisplacementDim>
32 std::string name,
33 MeshLib::Mesh& mesh,
34 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
35 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
36 unsigned const integration_order,
37 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
38 process_variables,
40 SecondaryVariableCollection&& secondary_variables,
41 bool const use_monolithic_scheme)
42 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
43 integration_order, std::move(process_variables),
44 std::move(secondary_variables), use_monolithic_scheme),
45 _process_data(std::move(process_data))
46{
48 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
49
51 mesh, "MassFlowRate", MeshLib::MeshItemType::Node, 1);
52
55 _integration_point_writer, integration_order,
57}
58
59template <int DisplacementDim>
61{
62 return false;
63}
64
65template <int DisplacementDim>
68 const int process_id) const
69{
70 // For the monolithic scheme or the M process (deformation) in the staggered
71 // scheme.
72 if (_use_monolithic_scheme || process_id == 1)
73 {
74 auto const& l = *_local_to_global_index_map;
75 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
76 &l.getGhostIndices(), &this->_sparsity_pattern};
77 }
78
79 // For staggered scheme and H process (pressure).
80 auto const& l = *_local_to_global_index_map_with_base_nodes;
81 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
82 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
83}
84
85template <int DisplacementDim>
87{
88 // Create single component dof in every of the mesh's nodes.
89 _mesh_subset_all_nodes =
90 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
91 // Create single component dof in the mesh's base nodes.
92 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
93 _mesh_subset_base_nodes =
94 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
95
96 // TODO move the two data members somewhere else.
97 // for extrapolation of secondary variables of stress or strain
98 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
99 *_mesh_subset_all_nodes};
100 _local_to_global_index_map_single_component =
101 std::make_unique<NumLib::LocalToGlobalIndexMap>(
102 std::move(all_mesh_subsets_single_component),
103 // by location order is needed for output
105
106 if (_use_monolithic_scheme)
107 {
108 // For pressure, which is the first
109 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
110 *_mesh_subset_base_nodes};
111
112 // For displacement.
113 const int monolithic_process_id = 0;
114 std::generate_n(std::back_inserter(all_mesh_subsets),
115 getProcessVariables(monolithic_process_id)[1]
116 .get()
117 .getNumberOfGlobalComponents(),
118 [&]() { return *_mesh_subset_all_nodes; });
119
120 std::vector<int> const vec_n_components{1, DisplacementDim};
121 _local_to_global_index_map =
122 std::make_unique<NumLib::LocalToGlobalIndexMap>(
123 std::move(all_mesh_subsets), vec_n_components,
125 assert(_local_to_global_index_map);
126 }
127 else
128 {
129 // For displacement equation.
130 const int process_id = 1;
131 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
132 std::generate_n(std::back_inserter(all_mesh_subsets),
133 getProcessVariables(process_id)[0]
134 .get()
135 .getNumberOfGlobalComponents(),
136 [&]() { return *_mesh_subset_all_nodes; });
137
138 std::vector<int> const vec_n_components{DisplacementDim};
139 _local_to_global_index_map =
140 std::make_unique<NumLib::LocalToGlobalIndexMap>(
141 std::move(all_mesh_subsets), vec_n_components,
143
144 // For pressure equation.
145 // Collect the mesh subsets with base nodes in a vector.
146 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
147 *_mesh_subset_base_nodes};
148 _local_to_global_index_map_with_base_nodes =
149 std::make_unique<NumLib::LocalToGlobalIndexMap>(
150 std::move(all_mesh_subsets_base_nodes),
151 // by location order is needed for output
153
154 _sparsity_pattern_with_linear_element = NumLib::computeSparsityPattern(
155 *_local_to_global_index_map_with_base_nodes, _mesh);
156
157 assert(_local_to_global_index_map);
158 assert(_local_to_global_index_map_with_base_nodes);
159 }
160}
161
162template <int DisplacementDim>
164 NumLib::LocalToGlobalIndexMap const& dof_table,
165 MeshLib::Mesh const& mesh,
166 unsigned const integration_order)
167{
170 mesh.getElements(), dof_table, _local_assemblers,
171 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
172 _process_data);
173
175 LocalAssemblerIF::getReflectionDataForOutput(), _secondary_variables,
176 getExtrapolator(), _local_assemblers);
177
178 auto add_secondary_variable = [&](std::string const& name,
179 int const num_components,
180 auto get_ip_values_function)
181 {
182 _secondary_variables.addSecondaryVariable(
183 name,
184 makeExtrapolator(num_components, getExtrapolator(),
185 _local_assemblers,
186 std::move(get_ip_values_function)));
187 };
188
189 //
190 // enable output of internal variables defined by material models
191 //
193 LocalAssemblerIF>(_process_data.solid_materials,
194 add_secondary_variable);
195
198 _process_data.solid_materials, _local_assemblers,
199 _integration_point_writer, integration_order);
200
201 _process_data.element_saturation = MeshLib::getOrCreateMeshProperty<double>(
202 const_cast<MeshLib::Mesh&>(mesh), "saturation_avg",
204
205 _process_data.element_porosity = MeshLib::getOrCreateMeshProperty<double>(
206 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
208
209 _process_data.element_stresses = MeshLib::getOrCreateMeshProperty<double>(
210 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
213 DisplacementDim>::RowsAtCompileTime);
214
215 _process_data.pressure_interpolated =
217 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
219
220 setIPDataInitialConditions(_integration_point_writer, mesh.getProperties(),
221 _local_assemblers);
222
223 // Initialize local assemblers after all variables have been set.
224 GlobalExecutor::executeMemberOnDereferenced(&LocalAssemblerIF::initialize,
225 _local_assemblers,
226 *_local_to_global_index_map);
227}
228
229template <int DisplacementDim>
231 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
232{
233 if (_use_monolithic_scheme)
234 {
235 const int monolithic_process_id = 0;
236 initializeProcessBoundaryConditionsAndSourceTerms(
237 *_local_to_global_index_map, monolithic_process_id, media);
238 return;
239 }
240
241 // Staggered scheme:
242 // for the equations of pressure
243 const int hydraulic_process_id = 0;
244 initializeProcessBoundaryConditionsAndSourceTerms(
245 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id,
246 media);
247
248 // for the equations of deformation.
249 const int mechanical_process_id = 1;
250 initializeProcessBoundaryConditionsAndSourceTerms(
251 *_local_to_global_index_map, mechanical_process_id, media);
252}
253
254template <int DisplacementDim>
256 setInitialConditionsConcreteProcess(std::vector<GlobalVector*>& x,
257 double const t,
258 int const process_id)
259{
260 if (process_id != 0)
261 {
262 return;
263 }
264
265 DBUG("SetInitialConditions RichardsMechanicsProcess.");
266
268 &LocalAssemblerIF::setInitialConditions, _local_assemblers,
269 getActiveElementIDs(), getDOFTables(x.size()), x, t, process_id);
270}
271
272template <int DisplacementDim>
274 const double t, double const dt, std::vector<GlobalVector*> const& x,
275 std::vector<GlobalVector*> const& x_prev, int const process_id,
277{
278 DBUG("Assemble the equations for RichardsMechanics");
279
280 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
281 _local_to_global_index_map.get()};
282
283 // Call global assembler for each local assembly item.
285 _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
286 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
287 &b);
288}
289
290template <int DisplacementDim>
293 const double t, double const dt, std::vector<GlobalVector*> const& x,
294 std::vector<GlobalVector*> const& x_prev, int const process_id,
295 GlobalVector& b, GlobalMatrix& Jac)
296{
297 // For the monolithic scheme
298 if (_use_monolithic_scheme)
299 {
300 DBUG(
301 "Assemble the Jacobian of RichardsMechanics for the monolithic"
302 " scheme.");
303 }
304 else
305 {
306 // For the staggered scheme
307 if (process_id == 0)
308 {
309 DBUG(
310 "Assemble the Jacobian equations of liquid fluid process in "
311 "RichardsMechanics for the staggered scheme.");
312 }
313 else
314 {
315 DBUG(
316 "Assemble the Jacobian equations of mechanical process in "
317 "RichardsMechanics for the staggered scheme.");
318 }
319 }
320
321 auto const dof_tables = getDOFTables(x.size());
324 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
325 process_id, &b, &Jac);
326
327 auto copyRhs = [&](int const variable_id, auto& output_vector)
328 {
329 if (_use_monolithic_scheme)
330 {
331 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
332 output_vector,
333 std::negate<double>());
334 }
335 else
336 {
337 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
338 output_vector,
339 std::negate<double>());
340 }
341 };
342 if (_use_monolithic_scheme || process_id == 0)
343 {
344 copyRhs(0, *_hydraulic_flow);
345 }
346 if (_use_monolithic_scheme || process_id == 1)
347 {
348 copyRhs(1, *_nodal_forces);
349 }
350}
351
352template <int DisplacementDim>
354 std::vector<GlobalVector*> const& x,
355 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
356 const int process_id)
357{
358 if (hasMechanicalProcess(process_id))
359 {
360 DBUG("PostTimestep RichardsMechanicsProcess.");
361
363 &LocalAssemblerIF::postTimestep, _local_assemblers,
364 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
365 process_id);
366 }
367}
368
369template <int DisplacementDim>
371 computeSecondaryVariableConcrete(const double t, const double dt,
372 std::vector<GlobalVector*> const& x,
373 GlobalVector const& x_prev,
374 int const process_id)
375{
376 if (process_id != 0)
377 {
378 return;
379 }
380
381 DBUG("Compute the secondary variables for RichardsMechanicsProcess.");
382
384 &LocalAssemblerIF::computeSecondaryVariable, _local_assemblers,
385 getActiveElementIDs(), getDOFTables(x.size()), t, dt, x, x_prev,
386 process_id);
387}
388
389template <int DisplacementDim>
390std::tuple<NumLib::LocalToGlobalIndexMap*, bool> RichardsMechanicsProcess<
391 DisplacementDim>::getDOFTableForExtrapolatorData() const
392{
393 const bool manage_storage = false;
394 return std::make_tuple(_local_to_global_index_map_single_component.get(),
395 manage_storage);
396}
397
398template <int DisplacementDim>
401 const int process_id) const
402{
403 if (hasMechanicalProcess(process_id))
404 {
405 return *_local_to_global_index_map;
406 }
407
408 // For the equation of pressure
409 return *_local_to_global_index_map_with_base_nodes;
410}
411
412template class RichardsMechanicsProcess<2>;
413template class RichardsMechanicsProcess<3>;
414
415} // namespace RichardsMechanics
416} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Global vector based on Eigen vector.
Definition EigenVector.h:25
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
Properties & getProperties()
Definition Mesh.h:134
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
std::vector< std::unique_ptr< LocalAssemblerIF > > _local_assemblers
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
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
RichardsMechanicsProcess(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, RichardsMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, 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, GlobalVector &b, GlobalMatrix &Jac) override
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, const int process_id) override
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
Handles configuration of several secondary variables from the project file.
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)
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, GlobalVector *b, GlobalMatrix *Jac)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.
void solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, std::vector< std::unique_ptr< LocalAssemblerInterface > > const &local_assemblers, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writer, int const integration_order)
void solidMaterialInternalToSecondaryVariables(std::map< int, std::shared_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
void addReflectedSecondaryVariables(ReflData const &reflection_data, SecondaryVariableCollection &secondary_variables, NumLib::Extrapolator &extrapolator, std::vector< std::unique_ptr< LocAsmIF > > const &local_assemblers)
void addReflectedIntegrationPointWriters(ReflData const &reflection_data, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writers, unsigned const integration_order, std::vector< std::unique_ptr< LocAsmIF > > const &local_assemblers)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
void createLocalAssemblersHM(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
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 executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)