26template <
int DisplacementDim>
30 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
31 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
32 unsigned const integration_order,
33 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
37 bool const use_monolithic_scheme)
38 :
Process(std::move(
name), mesh, std::move(jacobian_assembler), parameters,
39 integration_order, std::move(process_variables),
40 std::move(secondary_variables), use_monolithic_scheme),
41 _process_data(std::move(process_data))
53 std::make_unique<IntegrationPointWriter>(
60 std::make_unique<IntegrationPointWriter>(
67 std::make_unique<IntegrationPointWriter>(
68 "saturation_ip", 1 , integration_order,
72 std::make_unique<IntegrationPointWriter>(
79template <
int DisplacementDim>
85template <
int DisplacementDim>
88 const int process_id)
const
92 if (_use_monolithic_scheme || process_id == deformation_process_id)
94 auto const& l = *_local_to_global_index_map;
95 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
96 &l.getGhostIndices(), &this->_sparsity_pattern};
100 auto const& l = *_local_to_global_index_map_with_base_nodes;
101 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
102 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
105template <
int DisplacementDim>
109 _mesh_subset_all_nodes = std::make_unique<MeshLib::MeshSubset>(
110 _mesh, _mesh.getNodes(), _process_data.use_TaylorHood_elements);
113 _mesh_subset_base_nodes = std::make_unique<MeshLib::MeshSubset>(
114 _mesh, _base_nodes, _process_data.use_TaylorHood_elements);
118 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
119 *_mesh_subset_all_nodes};
120 _local_to_global_index_map_single_component =
121 std::make_unique<NumLib::LocalToGlobalIndexMap>(
122 std::move(all_mesh_subsets_single_component),
126 if (_use_monolithic_scheme)
129 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
130 *_mesh_subset_base_nodes};
133 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
136 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
140 std::back_inserter(all_mesh_subsets),
141 getProcessVariables(monolithic_process_id)[deformation_process_id]
143 .getNumberOfGlobalComponents(),
144 [&]() {
return *_mesh_subset_all_nodes; });
146 std::vector<int>
const vec_n_components{
147 n_gas_pressure_components, n_capillary_pressure_components,
148 n_temperature_components, n_displacement_components};
150 _local_to_global_index_map =
151 std::make_unique<NumLib::LocalToGlobalIndexMap>(
152 std::move(all_mesh_subsets), vec_n_components,
154 assert(_local_to_global_index_map);
158 OGS_FATAL(
"A Staggered version of TH2M is not implemented.");
162template <
int DisplacementDim>
166 unsigned const integration_order)
168 ProcessLib::createLocalAssemblersHM<DisplacementDim, TH2MLocalAssembler>(
172 auto add_secondary_variable = [&](std::string
const&
name,
173 int const num_components,
174 auto get_ip_values_function)
176 _secondary_variables.addSecondaryVariable(
180 std::move(get_ip_values_function)));
183 add_secondary_variable(
"sigma",
185 DisplacementDim>::RowsAtCompileTime,
187 add_secondary_variable(
"swelling_stress",
189 DisplacementDim>::RowsAtCompileTime,
191 add_secondary_variable(
"epsilon",
193 DisplacementDim>::RowsAtCompileTime,
195 add_secondary_variable(
"velocity_gas", mesh.
getDimension(),
197 add_secondary_variable(
200 add_secondary_variable(
"saturation", 1,
202 add_secondary_variable(
"vapour_pressure", 1,
204 add_secondary_variable(
"porosity", 1,
206 add_secondary_variable(
"gas_density", 1,
208 add_secondary_variable(
"solid_density", 1,
210 add_secondary_variable(
"liquid_density", 1,
212 add_secondary_variable(
"mole_fraction_gas", 1,
214 add_secondary_variable(
"mass_fraction_gas", 1,
216 add_secondary_variable(
217 "mass_fraction_liquid", 1,
220 add_secondary_variable(
221 "relative_permeability_gas", 1,
223 add_secondary_variable(
224 "relative_permeability_liquid", 1,
227 add_secondary_variable(
"enthalpy_gas", 1,
230 add_secondary_variable(
"enthalpy_liquid", 1,
233 add_secondary_variable(
"enthalpy_solid", 1,
236 _process_data.element_saturation = MeshLib::getOrCreateMeshProperty<double>(
240 _process_data.gas_pressure_interpolated =
241 MeshLib::getOrCreateMeshProperty<double>(
242 const_cast<MeshLib::Mesh&
>(mesh),
"gas_pressure_interpolated",
245 _process_data.capillary_pressure_interpolated =
246 MeshLib::getOrCreateMeshProperty<double>(
247 const_cast<MeshLib::Mesh&
>(mesh),
"capillary_pressure_interpolated",
250 _process_data.liquid_pressure_interpolated =
251 MeshLib::getOrCreateMeshProperty<double>(
252 const_cast<MeshLib::Mesh&
>(mesh),
"liquid_pressure_interpolated",
255 _process_data.temperature_interpolated =
256 MeshLib::getOrCreateMeshProperty<double>(
257 const_cast<MeshLib::Mesh&
>(mesh),
"temperature_interpolated",
261 for (
auto const& ip_writer : _integration_point_writer)
264 auto const&
name = ip_writer->name();
269 auto const& _meshproperty =
273 if (_meshproperty.getMeshItemType() !=
279 auto const ip_meta_data =
283 if (ip_meta_data.n_components !=
284 _meshproperty.getNumberOfGlobalComponents())
287 "Different number of components in meta data ({:d}) than in "
288 "the integration point field data for '{:s}': {:d}.",
289 ip_meta_data.n_components,
name,
290 _meshproperty.getNumberOfGlobalComponents());
295 std::size_t position = 0;
296 for (
auto& local_asm : _local_assemblers)
298 std::size_t
const integration_points_read =
299 local_asm->setIPDataInitialConditions(
300 name, &_meshproperty[position],
301 ip_meta_data.integration_order);
302 if (integration_points_read == 0)
305 "No integration points read in the integration point "
306 "initial conditions set function.");
308 position += integration_points_read * ip_meta_data.n_components;
315 *_local_to_global_index_map);
318template <
int DisplacementDim>
321 if (_use_monolithic_scheme)
323 initializeProcessBoundaryConditionsAndSourceTerms(
324 *_local_to_global_index_map, monolithic_process_id);
329 OGS_FATAL(
"A Staggered version of TH2M is not implemented.");
332template <
int DisplacementDim>
334 std::vector<GlobalVector*>& x,
double const t,
int const process_id)
341 DBUG(
"Set initial conditions of TH2MProcess.");
345 getDOFTable(process_id), *x[process_id],
t, _use_monolithic_scheme,
349template <
int DisplacementDim>
351 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
352 std::vector<GlobalVector*>
const& xdot,
int const process_id,
355 DBUG(
"Assemble the equations for TH2M");
357 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
358 dof_table = {std::ref(*_local_to_global_index_map)};
368template <
int DisplacementDim>
370 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
371 std::vector<GlobalVector*>
const& xdot,
int const process_id,
374 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
377 if (_use_monolithic_scheme)
379 DBUG(
"Assemble the Jacobian of TH2M for the monolithic scheme.");
380 dof_tables.emplace_back(*_local_to_global_index_map);
385 OGS_FATAL(
"A Staggered version of TH2M is not implemented.");
393 process_id, M, K, b, Jac);
395 auto copyRhs = [&](
int const variable_id,
auto& output_vector)
397 if (_use_monolithic_scheme)
401 std::negate<double>());
407 std::negate<double>());
410 if (_use_monolithic_scheme || process_id == 1)
412 copyRhs(0, *_hydraulic_flow);
414 if (_use_monolithic_scheme || process_id == 2)
416 copyRhs(1, *_nodal_forces);
420template <
int DisplacementDim>
422 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
423 const int process_id)
425 DBUG(
"PreTimestep TH2MProcess.");
427 if (hasMechanicalProcess(process_id))
430 getProcessVariables(process_id)[0];
435 *x[process_id],
t, dt);
439template <
int DisplacementDim>
441 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
442 const int process_id)
444 DBUG(
"PostTimestep TH2MProcess.");
445 auto const dof_tables = getDOFTables(x.size());
452template <
int DisplacementDim>
454 double const t,
double const dt, std::vector<GlobalVector*>
const& x,
462 DBUG(
"Compute the secondary variables for TH2MProcess.");
463 auto const dof_tables = getDOFTables(x.size());
471template <
int DisplacementDim>
472std::tuple<NumLib::LocalToGlobalIndexMap*, bool>
475 const bool manage_storage =
false;
476 return std::make_tuple(_local_to_global_index_map_single_component.get(),
480template <
int DisplacementDim>
482 const int process_id)
const
484 if (hasMechanicalProcess(process_id))
486 return *_local_to_global_index_map;
490 return *_local_to_global_index_map_with_base_nodes;
493template <
int DisplacementDim>
494std::vector<NumLib::LocalToGlobalIndexMap const*>
497 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
498 dof_tables.reserve(number_of_processes);
499 std::generate_n(std::back_inserter(dof_tables), number_of_processes,
500 [&]() {
return &getDOFTable(dof_tables.size()); });
void DBUG(char const *fmt, Args const &... args)
Global vector based on Eigen vector.
bool isAxiallySymmetric() const
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Properties & getProperties()
bool existsPropertyVector(std::string const &name) const
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
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)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
void setInitialConditions(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Handles configuration of several secondary variables from the project file.
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
bool isLinear() const override
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
MeshLib::PropertyVector< double > * _nodal_forces
TH2MProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const ¶meters, 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)
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double dt, int const) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
MeshLib::PropertyVector< double > * _hydraulic_flow
void constructDofTable() override
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, const int process_id) 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< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
void initializeBoundaryConditions() override
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
@ BY_LOCATION
Ordering data by spatial location.
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor mapFunction)
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Properties const &properties, std::string const &name)
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 executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtLiquidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > getEpsilon() const =0
virtual std::vector< double > const & getIntPtRelativePermeabilityLiquid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEnthalpyGas(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtVapourPressure(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtRelativePermeabilityGas(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSolidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtMoleFractionGas(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtMassFractionGas(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSwellingStress(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEnthalpyLiquid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > getSaturation() const =0
virtual std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEnthalpySolid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocityGas(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocityLiquid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > getSwellingStress() const =0
virtual std::vector< double > const & getIntPtGasDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtMassFractionLiquid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > getSigma() const =0