25 namespace ThermoHydroMechanics
27 template <
int DisplacementDim>
31 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
32 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
33 unsigned const integration_order,
34 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
38 bool const use_monolithic_scheme)
39 :
Process(std::move(
name), mesh, std::move(jacobian_assembler), parameters,
40 integration_order, std::move(process_variables),
41 std::move(secondary_variables), use_monolithic_scheme),
42 _process_data(std::move(process_data))
50 _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
54 std::make_unique<IntegrationPointWriter>(
61 std::make_unique<IntegrationPointWriter>(
68 template <
int DisplacementDim>
74 template <
int DisplacementDim>
77 const int process_id)
const
81 if (_use_monolithic_scheme || process_id == 2)
83 auto const& l = *_local_to_global_index_map;
84 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
85 &l.getGhostIndices(), &this->_sparsity_pattern};
89 auto const& l = *_local_to_global_index_map_with_base_nodes;
90 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
91 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
94 template <
int DisplacementDim>
98 _mesh_subset_all_nodes =
99 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
102 _mesh_subset_base_nodes =
103 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
107 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
108 *_mesh_subset_all_nodes};
109 _local_to_global_index_map_single_component =
110 std::make_unique<NumLib::LocalToGlobalIndexMap>(
111 std::move(all_mesh_subsets_single_component),
115 if (_use_monolithic_scheme)
118 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
119 *_mesh_subset_base_nodes};
122 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
125 const int monolithic_process_id = 0;
126 std::generate_n(std::back_inserter(all_mesh_subsets),
127 getProcessVariables(monolithic_process_id)[2]
129 .getNumberOfGlobalComponents(),
130 [&]() {
return *_mesh_subset_all_nodes; });
132 std::vector<int>
const vec_n_components{1, 1, DisplacementDim};
133 _local_to_global_index_map =
134 std::make_unique<NumLib::LocalToGlobalIndexMap>(
135 std::move(all_mesh_subsets), vec_n_components,
137 assert(_local_to_global_index_map);
142 const int process_id = 2;
143 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
144 std::generate_n(std::back_inserter(all_mesh_subsets),
145 getProcessVariables(process_id)[0]
147 .getNumberOfGlobalComponents(),
148 [&]() {
return *_mesh_subset_all_nodes; });
150 std::vector<int>
const vec_n_components{DisplacementDim};
151 _local_to_global_index_map =
152 std::make_unique<NumLib::LocalToGlobalIndexMap>(
153 std::move(all_mesh_subsets), vec_n_components,
158 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
159 *_mesh_subset_base_nodes};
160 _local_to_global_index_map_with_base_nodes =
161 std::make_unique<NumLib::LocalToGlobalIndexMap>(
162 std::move(all_mesh_subsets_base_nodes),
167 *_local_to_global_index_map_with_base_nodes, _mesh);
169 assert(_local_to_global_index_map);
170 assert(_local_to_global_index_map_with_base_nodes);
174 template <
int DisplacementDim>
178 unsigned const integration_order)
186 _secondary_variables.addSecondaryVariable(
189 DisplacementDim>::RowsAtCompileTime,
190 getExtrapolator(), _local_assemblers,
193 _secondary_variables.addSecondaryVariable(
196 DisplacementDim>::RowsAtCompileTime,
197 getExtrapolator(), _local_assemblers,
200 _secondary_variables.addSecondaryVariable(
206 _process_data.pressure_interpolated =
207 MeshLib::getOrCreateMeshProperty<double>(
211 _process_data.temperature_interpolated =
212 MeshLib::getOrCreateMeshProperty<double>(
213 const_cast<MeshLib::Mesh&
>(mesh),
"temperature_interpolated",
222 *_local_to_global_index_map);
225 template <
int DisplacementDim>
227 DisplacementDim>::initializeBoundaryConditions()
229 if (_use_monolithic_scheme)
231 const int process_id_of_thermohydromechancs = 0;
232 initializeProcessBoundaryConditionsAndSourceTerms(
233 *_local_to_global_index_map, process_id_of_thermohydromechancs);
239 const int thermal_process_id = 0;
240 initializeProcessBoundaryConditionsAndSourceTerms(
241 *_local_to_global_index_map_with_base_nodes, thermal_process_id);
244 const int hydraulic_process_id = 1;
245 initializeProcessBoundaryConditionsAndSourceTerms(
246 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id);
249 const int mechanical_process_id = 2;
250 initializeProcessBoundaryConditionsAndSourceTerms(
251 *_local_to_global_index_map, mechanical_process_id);
254 template <
int DisplacementDim>
256 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
257 std::vector<GlobalVector*>
const& xdot,
int const process_id,
260 DBUG(
"Assemble the equations for ThermoHydroMechanics");
262 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
263 dof_table = {std::ref(*_local_to_global_index_map)};
273 template <
int DisplacementDim>
276 std::vector<GlobalVector*>
const& x,
277 std::vector<GlobalVector*>
const& xdot,
282 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
285 if (_use_monolithic_scheme)
288 "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
290 dof_tables.emplace_back(*_local_to_global_index_map);
298 "Assemble the Jacobian equations of heat transport process in "
299 "ThermoHydroMechanics for the staggered scheme.");
301 else if (process_id == 1)
304 "Assemble the Jacobian equations of liquid fluid process in "
305 "ThermoHydroMechanics for the staggered scheme.");
310 "Assemble the Jacobian equations of mechanical process in "
311 "ThermoHydroMechanics for the staggered scheme.");
313 dof_tables.emplace_back(*_local_to_global_index_map_with_base_nodes);
314 dof_tables.emplace_back(*_local_to_global_index_map_with_base_nodes);
315 dof_tables.emplace_back(*_local_to_global_index_map);
323 process_id, M, K, b, Jac);
325 auto copyRhs = [&](
int const variable_id,
auto& output_vector)
327 if (_use_monolithic_scheme)
331 std::negate<double>());
337 std::negate<double>());
340 if (_use_monolithic_scheme || process_id == 0)
342 copyRhs(0, *_heat_flux);
344 if (_use_monolithic_scheme || process_id == 1)
346 copyRhs(1, *_hydraulic_flow);
348 if (_use_monolithic_scheme || process_id == 2)
350 copyRhs(2, *_nodal_forces);
354 template <
int DisplacementDim>
356 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
357 const int process_id)
359 DBUG(
"PreTimestep ThermoHydroMechanicsProcess.");
361 if (hasMechanicalProcess(process_id))
364 getProcessVariables(process_id)[0];
369 *x[process_id],
t, dt);
373 template <
int DisplacementDim>
375 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
376 const int process_id)
383 DBUG(
"PostTimestep ThermoHydroMechanicsProcess.");
384 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
385 auto const n_processes = x.size();
386 dof_tables.reserve(n_processes);
387 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
389 dof_tables.push_back(&getDOFTable(process_id));
399 template <
int DisplacementDim>
403 double const dt,
const int process_id)
405 if (!hasMechanicalProcess(process_id))
410 DBUG(
"PostNonLinearSolver ThermoHydroMechanicsProcess.");
418 _use_monolithic_scheme, process_id);
421 template <
int DisplacementDim>
424 std::vector<GlobalVector*>
const& x,
426 const int process_id)
433 DBUG(
"Compute the secondary variables for ThermoHydroMechanicsProcess.");
434 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
435 auto const n_processes = x.size();
436 dof_tables.reserve(n_processes);
437 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
439 dof_tables.push_back(&getDOFTable(process_id));
448 template <
int DisplacementDim>
450 DisplacementDim>::getDOFTableForExtrapolatorData()
const
452 const bool manage_storage =
false;
453 return std::make_tuple(_local_to_global_index_map_single_component.get(),
457 template <
int DisplacementDim>
460 const int process_id)
const
462 if (hasMechanicalProcess(process_id))
464 return *_local_to_global_index_map;
468 return *_local_to_global_index_map_with_base_nodes;
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Global vector based on Eigen vector.
bool isAxiallySymmetric() const
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Properties & getProperties()
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)
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)
void postNonLinearSolver(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, 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.
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
void postNonLinearSolverConcreteProcess(GlobalVector const &x, GlobalVector const &xdot, const double t, double const dt, int const process_id) override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
MeshLib::PropertyVector< double > * _hydraulic_flow
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 postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
MeshLib::PropertyVector< double > * _heat_flux
MeshLib::PropertyVector< double > * _nodal_forces
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, const int process_id) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
bool isLinear() const override
ThermoHydroMechanicsProcess(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, ThermoHydroMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
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
void constructDofTable() override
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)
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)
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.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.
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)
void setIPDataInitialConditions(std::vector< std::unique_ptr< 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, NumLib::IntegrationOrder const integration_order, ExtraCtorArgs &&... extra_ctor_args)
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 > getSigma() const =0
virtual std::vector< double > getEpsilon() 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 & getIntPtDarcyVelocity(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