29template <
int DisplacementDim>
33 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
34 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
35 unsigned const integration_order,
36 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
40 bool const use_monolithic_scheme)
41 :
Process(std::move(
name), mesh, std::move(jacobian_assembler), parameters,
42 integration_order, std::move(process_variables),
43 std::move(secondary_variables), use_monolithic_scheme),
46 _process_data(std::move(process_data))
58 _integration_point_writer.emplace_back(
59 std::make_unique<MeshLib::IntegrationPointWriter>(
61 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) ,
62 integration_order, local_assemblers_,
65 _integration_point_writer.emplace_back(
66 std::make_unique<MeshLib::IntegrationPointWriter>(
68 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) ,
69 integration_order, local_assemblers_,
72 _integration_point_writer.emplace_back(
73 std::make_unique<MeshLib::IntegrationPointWriter>(
75 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) ,
76 integration_order, local_assemblers_,
80template <
int DisplacementDim>
86template <
int DisplacementDim>
89 const int process_id)
const
96 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
102 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
106template <
int DisplacementDim>
111 std::make_unique<MeshLib::MeshSubset>(
_mesh,
_mesh.getNodes());
119 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
122 std::make_unique<NumLib::LocalToGlobalIndexMap>(
123 std::move(all_mesh_subsets_single_component),
130 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
137 const int monolithic_process_id = 0;
138 std::generate_n(std::back_inserter(all_mesh_subsets),
141 .getNumberOfGlobalComponents(),
144 std::vector<int>
const vec_n_components{1, 1, DisplacementDim};
146 std::make_unique<NumLib::LocalToGlobalIndexMap>(
147 std::move(all_mesh_subsets), vec_n_components,
154 const int process_id = 2;
155 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
156 std::generate_n(std::back_inserter(all_mesh_subsets),
159 .getNumberOfGlobalComponents(),
162 std::vector<int>
const vec_n_components{DisplacementDim};
164 std::make_unique<NumLib::LocalToGlobalIndexMap>(
165 std::move(all_mesh_subsets), vec_n_components,
170 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
173 std::make_unique<NumLib::LocalToGlobalIndexMap>(
174 std::move(all_mesh_subsets_base_nodes),
186template <
int DisplacementDim>
190 unsigned const integration_order)
198 auto add_secondary_variable = [&](std::string
const&
name,
199 int const num_components,
200 auto get_ip_values_function)
206 std::move(get_ip_values_function)));
209 add_secondary_variable(
212 DisplacementDim>::RowsAtCompileTime,
215 add_secondary_variable(
218 DisplacementDim>::RowsAtCompileTime,
221 add_secondary_variable(
224 DisplacementDim>::RowsAtCompileTime,
227 add_secondary_variable(
230 DisplacementDim>::RowsAtCompileTime,
233 add_secondary_variable(
236 DisplacementDim>::RowsAtCompileTime,
239 add_secondary_variable(
240 "ice_volume_fraction", 1,
243 add_secondary_variable(
247 add_secondary_variable(
251 add_secondary_variable(
272 DisplacementDim>::RowsAtCompileTime);
281 const_cast<MeshLib::Mesh&
>(mesh),
"temperature_interpolated",
289 add_secondary_variable);
305template <
int DisplacementDim>
307 std::map<
int, std::shared_ptr<MaterialPropertyLib::Medium>>
const& media)
311 const int process_id_of_thermohydromechancs = 0;
320 const int thermal_process_id = 0;
325 const int hydraulic_process_id = 1;
331 const int mechanical_process_id = 2;
336template <
int DisplacementDim>
340 int const process_id)
342 DBUG(
"SetInitialConditions ThermoHydroMechanicsProcess.");
349template <
int DisplacementDim>
351 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
352 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
355 DBUG(
"Assemble the equations for ThermoHydroMechanics");
358 t, dt, x, x_prev, process_id, M, K, b);
361template <
int DisplacementDim>
364 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
365 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
372 "Assemble the Jacobian of ThermoHydroMechanics for the monolithic "
381 "Assemble the Jacobian equations of heat transport process in "
382 "ThermoHydroMechanics for the staggered scheme.");
384 else if (process_id == 1)
387 "Assemble the Jacobian equations of liquid fluid process in "
388 "ThermoHydroMechanics for the staggered scheme.");
393 "Assemble the Jacobian equations of mechanical process in "
394 "ThermoHydroMechanics for the staggered scheme.");
401 assembleWithJacobian(t, dt, x, x_prev, process_id, b, Jac);
403 auto copyRhs = [&](
int const variable_id,
auto& output_vector)
407 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
409 std::negate<double>());
413 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
415 std::negate<double>());
432template <
int DisplacementDim>
434 std::vector<GlobalVector*>
const& x,
double const t,
double const dt,
435 const int process_id)
437 DBUG(
"PreTimestep ThermoHydroMechanicsProcess.");
447 updateActiveElements();
451template <
int DisplacementDim>
453 std::vector<GlobalVector*>
const& x,
454 std::vector<GlobalVector*>
const& x_prev,
double const t,
double const dt,
455 const int process_id)
462 DBUG(
"PostTimestep ThermoHydroMechanicsProcess.");
467 x_prev, t, dt, process_id);
470template <
int DisplacementDim>
471std::vector<std::vector<std::string>>
473 std::vector<std::reference_wrapper<MeshLib::Mesh>>
const& meshes)
475 INFO(
"ThermoHydroMechanicsProcess process initializeSubmeshOutput().");
476 std::vector<std::vector<std::string>> per_process_residuum_names;
479 per_process_residuum_names = {
480 {
"HeatFlowRate",
"VolumetricFlowRate",
"NodalForces"}};
484 per_process_residuum_names = {
485 {
"HeatFlowRate"}, {
"VolumetricFlowRate"}, {
"NodalForces"}};
489 initializeAssemblyOnSubmeshes(meshes, per_process_residuum_names);
491 return per_process_residuum_names;
494template <
int DisplacementDim>
497 std::vector<GlobalVector*>
const& x,
499 const int process_id)
506 DBUG(
"Compute the secondary variables for ThermoHydroMechanicsProcess.");
511 x, x_prev, process_id);
514template <
int DisplacementDim>
518 const bool manage_storage =
false;
523template <
int DisplacementDim>
526 const int process_id)
const
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
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()
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 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_prev, 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)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void setInitialConditions(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
void assemble(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) final
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< 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, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
std::vector< std::size_t > const & getActiveElementIDs() const
SecondaryVariableCollection _secondary_variables
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
NumLib::Extrapolator & getExtrapolator() const
GlobalSparsityPattern _sparsity_pattern
const bool _use_monolithic_scheme
Handles configuration of several secondary variables from the project file.
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
bool hasMechanicalProcess(int const process_id) const
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
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().
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
MeshLib::PropertyVector< double > * _hydraulic_flow
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
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)
void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
MeshLib::PropertyVector< double > * _heat_flux
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
MeshLib::PropertyVector< double > * _nodal_forces
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
std::vector< MeshLib::Node * > _base_nodes
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > local_assemblers_
bool isLinear() const override
ThermoHydroMechanicsProcessData< DisplacementDim > _process_data
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
GlobalSparsityPattern _sparsity_pattern_with_linear_element
void constructDofTable() override
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)
@ 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.
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 executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtFluidDensity(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 & getIntPtEpsilon0(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
virtual std::vector< double > getEpsilon() 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 > getEpsilonM() const =0
virtual std::vector< double > const & getIntPtSigmaIce(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 & getIntPtIceVolume(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 & getIntPtViscosity(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 & getIntPtEpsilonM(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