OGS
TH2MProcess.cpp
Go to the documentation of this file.
1
11#include "TH2MProcess.h"
12
13#include <cassert>
14
16#include "MaterialLib/SolidModels/MechanicsBase.h" // for the instantiation of process data
22#include "ProcessLib/Process.h"
27#include "TH2MProcessData.h"
28
29namespace ProcessLib
30{
31namespace TH2M
32{
33template <int DisplacementDim>
35 std::string name,
36 MeshLib::Mesh& mesh,
37 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
38 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
39 unsigned const integration_order,
40 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
41 process_variables,
43 SecondaryVariableCollection&& secondary_variables,
44 bool const use_monolithic_scheme)
45 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
46 integration_order, std::move(process_variables),
47 std::move(secondary_variables), use_monolithic_scheme),
49 _process_data(std::move(process_data))
50{
52 DisplacementDim>(
54 _integration_point_writer, integration_order, local_assemblers_);
55}
56
57template <int DisplacementDim>
59{
60 return false;
61}
62
63template <int DisplacementDim>
66 const int process_id) const
67{
68 // For the monolithic scheme or the M process (deformation) in the staggered
69 // scheme.
70 if (_use_monolithic_scheme || process_id == deformation_process_id)
71 {
72 auto const& l = *_local_to_global_index_map;
73 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
74 &l.getGhostIndices(), &this->_sparsity_pattern};
75 }
76
77 // For staggered scheme and T or H process (pressure).
78 auto const& l = *_local_to_global_index_map_with_base_nodes;
79 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
80 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
81}
82
83template <int DisplacementDim>
85{
86 // Create single component dof in every of the mesh's nodes.
87 _mesh_subset_all_nodes = std::make_unique<MeshLib::MeshSubset>(
88 _mesh, _mesh.getNodes(), _process_data.use_TaylorHood_elements);
89 // Create single component dof in the mesh's base nodes.
90 _base_nodes = MeshLib::getBaseNodes(_mesh.getElements());
91 _mesh_subset_base_nodes = std::make_unique<MeshLib::MeshSubset>(
92 _mesh, _base_nodes, _process_data.use_TaylorHood_elements);
93
94 // TODO move the two data members somewhere else.
95 // for extrapolation of secondary variables of stress or strain
96 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
97 *_mesh_subset_all_nodes};
98 _local_to_global_index_map_single_component =
99 std::make_unique<NumLib::LocalToGlobalIndexMap>(
100 std::move(all_mesh_subsets_single_component),
101 // by location order is needed for output
103
104 if (_use_monolithic_scheme)
105 {
106 // For gas pressure, which is the first
107 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
108 *_mesh_subset_base_nodes};
109
110 // For capillary pressure, which is the second
111 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
112
113 // For temperature, which is the third
114 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
115
116 // For displacement.
117 std::generate_n(
118 std::back_inserter(all_mesh_subsets),
119 getProcessVariables(monolithic_process_id)[deformation_process_id]
120 .get()
121 .getNumberOfGlobalComponents(),
122 [&]() { return *_mesh_subset_all_nodes; });
123
124 std::vector<int> const vec_n_components{
125 n_gas_pressure_components, n_capillary_pressure_components,
126 n_temperature_components, n_displacement_components};
127
128 _local_to_global_index_map =
129 std::make_unique<NumLib::LocalToGlobalIndexMap>(
130 std::move(all_mesh_subsets), vec_n_components,
132 assert(_local_to_global_index_map);
133 }
134 else
135 {
136 OGS_FATAL("A Staggered version of TH2M is not implemented.");
137 }
138}
139
140template <int DisplacementDim>
142 NumLib::LocalToGlobalIndexMap const& dof_table,
143 MeshLib::Mesh const& mesh,
144 unsigned const integration_order)
145{
146 createLocalAssemblers<DisplacementDim>(
147 mesh.getElements(), dof_table, local_assemblers_,
148 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
149 _process_data);
150
151 auto add_secondary_variable = [&](std::string const& name,
152 int const num_components,
153 auto get_ip_values_function)
154 {
155 _secondary_variables.addSecondaryVariable(
156 name,
157 makeExtrapolator(num_components, getExtrapolator(),
158 local_assemblers_,
159 std::move(get_ip_values_function)));
160 };
161
162 ProcessLib::Reflection::addReflectedSecondaryVariables<DisplacementDim>(
164 _secondary_variables, getExtrapolator(), local_assemblers_);
165
166 add_secondary_variable(
167 "velocity_gas", mesh.getDimension(),
169 add_secondary_variable(
170 "velocity_liquid", mesh.getDimension(),
172 add_secondary_variable(
173 "diffusion_velocity_vapour_gas", mesh.getDimension(),
175 DisplacementDim>::getIntPtDiffusionVelocityVapourGas);
176 add_secondary_variable(
177 "diffusion_velocity_gas_gas", mesh.getDimension(),
179 DisplacementDim>::getIntPtDiffusionVelocityGasGas);
180 add_secondary_variable(
181 "diffusion_velocity_solute_liquid", mesh.getDimension(),
183 DisplacementDim>::getIntPtDiffusionVelocitySoluteLiquid);
184 add_secondary_variable(
185 "diffusion_velocity_liquid_liquid", mesh.getDimension(),
187 DisplacementDim>::getIntPtDiffusionVelocityLiquidLiquid);
188
189 add_secondary_variable(
190 "enthalpy_solid", 1,
192
193 //
194 // enable output of internal variables defined by material models
195 //
197 LocalAssemblerInterface<DisplacementDim>>(_process_data.solid_materials,
198 add_secondary_variable);
199
202 _process_data.solid_materials, local_assemblers_,
203 _integration_point_writer, integration_order);
204
205 _process_data.element_saturation = MeshLib::getOrCreateMeshProperty<double>(
206 const_cast<MeshLib::Mesh&>(mesh), "saturation_avg",
208
209 _process_data.gas_pressure_interpolated =
210 MeshLib::getOrCreateMeshProperty<double>(
211 const_cast<MeshLib::Mesh&>(mesh), "gas_pressure_interpolated",
213
214 _process_data.capillary_pressure_interpolated =
215 MeshLib::getOrCreateMeshProperty<double>(
216 const_cast<MeshLib::Mesh&>(mesh), "capillary_pressure_interpolated",
218
219 _process_data.liquid_pressure_interpolated =
220 MeshLib::getOrCreateMeshProperty<double>(
221 const_cast<MeshLib::Mesh&>(mesh), "liquid_pressure_interpolated",
223
224 _process_data.temperature_interpolated =
225 MeshLib::getOrCreateMeshProperty<double>(
226 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
228
229 setIPDataInitialConditions(_integration_point_writer, mesh.getProperties(),
230 local_assemblers_);
231
232 // Initialize local assemblers after all variables have been set.
235 local_assemblers_, *_local_to_global_index_map);
236}
237
238template <int DisplacementDim>
240 std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
241{
242 if (_use_monolithic_scheme)
243 {
244 initializeProcessBoundaryConditionsAndSourceTerms(
245 *_local_to_global_index_map, monolithic_process_id, media);
246 return;
247 }
248
249 // Staggered scheme:
250 OGS_FATAL("A Staggered version of TH2M is not implemented.");
251}
252
253template <int DisplacementDim>
255 std::vector<GlobalVector*>& x, double const t, int const process_id)
256{
257 if (process_id != 0)
258 {
259 return;
260 }
261
262 DBUG("Set initial conditions of TH2MProcess.");
263
266 local_assemblers_, getDOFTables(x.size()), x, t, process_id);
267}
268
269template <int DisplacementDim>
271 const double t, double const dt, std::vector<GlobalVector*> const& x,
272 std::vector<GlobalVector*> const& x_prev, int const process_id,
274{
275 DBUG("Assemble the equations for TH2M");
276
277 AssemblyMixin<TH2MProcess<DisplacementDim>>::assemble(t, dt, x, x_prev,
278 process_id, M, K, b);
279}
280
281template <int DisplacementDim>
283 const double t, double const dt, std::vector<GlobalVector*> const& x,
284 std::vector<GlobalVector*> const& x_prev, int const process_id,
286{
287 if (!_use_monolithic_scheme)
288 {
289 OGS_FATAL("A Staggered version of TH2M is not implemented.");
290 }
291
292 AssemblyMixin<TH2MProcess<DisplacementDim>>::assembleWithJacobian(
293 t, dt, x, x_prev, process_id, M, K, b, Jac);
294}
295
296template <int DisplacementDim>
298 std::vector<GlobalVector*> const& x, double const t, double const dt,
299 const int process_id)
300{
301 DBUG("PreTimestep TH2MProcess.");
302
303 if (hasMechanicalProcess(process_id))
304 {
306 getProcessVariables(process_id)[0];
307
310 local_assemblers_, pv.getActiveElementIDs(),
311 *_local_to_global_index_map, *x[process_id], t, dt);
312 }
313
314 AssemblyMixin<TH2MProcess<DisplacementDim>>::updateActiveElements(
315 process_id);
316}
317
318template <int DisplacementDim>
320 std::vector<GlobalVector*> const& x,
321 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
322 const int process_id)
323{
324 DBUG("PostTimestep TH2MProcess.");
325
326 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
329 local_assemblers_, pv.getActiveElementIDs(), getDOFTables(x.size()), x,
330 x_prev, t, dt, process_id);
331}
332
333template <int DisplacementDim>
335 double const t, double const dt, std::vector<GlobalVector*> const& x,
336 GlobalVector const& x_prev, const int process_id)
337{
338 if (process_id != 0)
339 {
340 return;
341 }
342
343 DBUG("Compute the secondary variables for TH2MProcess.");
344
345 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
348 local_assemblers_, pv.getActiveElementIDs(), getDOFTables(x.size()), t,
349 dt, x, x_prev, process_id);
350}
351
352template <int DisplacementDim>
353std::vector<std::string>
355 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
356{
357 INFO("TH2M process initializeSubmeshOutput().");
358 const int process_id = 0;
359 std::vector<std::string> residuum_names{
360 "GasMassFlowRate", "LiquidMassFlowRate", "HeatFlowRate", "NodalForces"};
361
362 AssemblyMixin<TH2MProcess<DisplacementDim>>::initializeAssemblyOnSubmeshes(
363 process_id, meshes, residuum_names);
364
365 return residuum_names;
366}
367
368template <int DisplacementDim>
369std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
371{
372 const bool manage_storage = false;
373 return std::make_tuple(_local_to_global_index_map_single_component.get(),
374 manage_storage);
375}
376
377template <int DisplacementDim>
379 const int process_id) const
380{
381 if (hasMechanicalProcess(process_id))
382 {
383 return *_local_to_global_index_map;
384 }
385
386 // For the equation of pressure
387 return *_local_to_global_index_map_with_base_nodes;
388}
389
390template class TH2MProcess<2>;
391template class TH2MProcess<3>;
392
393} // namespace TH2M
394} // namespace ProcessLib
#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
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
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
Properties & getProperties()
Definition Mesh.h:134
std::vector< std::size_t > const & getActiveElementIDs() const
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:366
Handles configuration of several secondary variables from the project file.
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const) override
bool isLinear() const override
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
TH2MProcess(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, TH2MProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
std::vector< std::string > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) override
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 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
void constructDofTable() override
NumLib::LocalToGlobalIndexMap const & getDOFTable(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().
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ BY_LOCATION
Ordering data by spatial location.
void solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::unique_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::unique_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
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)
static void executeSelectedMemberOnDereferenced(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)