OGS
ProcessLib::TH2M::TH2MProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::TH2M::TH2MProcess< DisplacementDim >

Thermally induced deformation process in linear kinematics poro-mechanical/triphasic model.

Definition at line 28 of file TH2MProcess.h.

#include <TH2MProcess.h>

Inheritance diagram for ProcessLib::TH2M::TH2MProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::TH2M::TH2MProcess< DisplacementDim >:
[legend]

Public Member Functions

 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)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< 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, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
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
 
void assembleWithJacobian (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) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
MeshLib::MeshgetMesh () const
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
bool requiresNormalization () const override
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual ~SubmeshAssemblySupport ()=default
 

Private Member Functions

void constructDofTable () override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
void initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
 
void setInitialConditionsConcreteProcess (std::vector< GlobalVector * > &x, double const t, int const 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, 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, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const) override
 
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override
 
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, const int process_id) override
 
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const override
 
bool hasMechanicalProcess (int const process_id) const
 
- Private Member Functions inherited from ProcessLib::AssemblyMixin< TH2MProcess< DisplacementDim > >
void initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
 
void updateActiveElements ()
 
void assemble (double const t, double const dt, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const process_id, GlobalMatrix &, GlobalMatrix &, GlobalVector &)
 
void assembleWithJacobian (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac)
 

Private Attributes

std::vector< MeshLib::Node * > _base_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
 
TH2MProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > local_assemblers_
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_with_base_nodes
 
GlobalSparsityPattern _sparsity_pattern_with_linear_element
 

Static Private Attributes

static constexpr int monolithic_process_id = 0
 
static constexpr int deformation_process_id = 3
 
static constexpr int n_gas_pressure_components = 1
 
static constexpr int n_capillary_pressure_components = 1
 
static constexpr int n_temperature_components = 1
 
static constexpr int n_displacement_components = DisplacementDim
 

Friends

class AssemblyMixin< TH2MProcess< DisplacementDim > >
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
CellAverageData cell_average_data_
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Constructor & Destructor Documentation

◆ TH2MProcess()

template<int DisplacementDim>
ProcessLib::TH2M::TH2MProcess< DisplacementDim >::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 )

Definition at line 35 of file TH2MProcess.cpp.

46 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
47 integration_order, std::move(process_variables),
48 std::move(secondary_variables), use_monolithic_scheme),
49 AssemblyMixin<TH2MProcess<DisplacementDim>>{*_jacobian_assembler},
50 _process_data(std::move(process_data))
51{
53 DisplacementDim>(
56}
std::string const name
Definition Process.h:362
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:376
TH2MProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface< DisplacementDim > > > 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)

References ProcessLib::Process::_jacobian_assembler.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 240 of file TH2MProcess.cpp.

244{
245 DBUG("Assemble the equations for TH2M");
246
247 AssemblyMixin<TH2MProcess<DisplacementDim>>::assemble(t, dt, x, x_prev,
248 process_id, M, K, b);
249}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
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
Definition Process.cpp:266

References DBUG().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 252 of file TH2MProcess.cpp.

256{
258 {
259 OGS_FATAL("A Staggered version of TH2M is not implemented.");
260 }
261
262 AssemblyMixin<TH2MProcess<DisplacementDim>>::assembleWithJacobian(
263 t, dt, x, x_prev, process_id, b, Jac);
264}
#define OGS_FATAL(...)
Definition Error.h:26
void assembleWithJacobian(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) final
Definition Process.cpp:287
const bool _use_monolithic_scheme
Definition Process.h:379

References OGS_FATAL.

◆ computeSecondaryVariableConcrete()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::computeSecondaryVariableConcrete ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 299 of file TH2MProcess.cpp.

302{
303 if (process_id != 0)
304 {
305 return;
306 }
307
308 DBUG("Compute the secondary variables for TH2MProcess.");
309
313 x, x_prev, process_id);
314
315 computeCellAverages<DisplacementDim>(cell_average_data_, local_assemblers_);
316}
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)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
CellAverageData cell_average_data_
Definition Process.h:372
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::constructDofTable ( )
overrideprivatevirtual

This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map. Other cases can be considered by overloading this member function in the derived class.

Reimplemented from ProcessLib::Process.

Definition at line 85 of file TH2MProcess.cpp.

86{
87 // Create single component dof in every of the mesh's nodes.
88 _mesh_subset_all_nodes = std::make_unique<MeshLib::MeshSubset>(
89 _mesh, _mesh.getNodes(), _process_data.use_TaylorHood_elements);
90 // Create single component dof in the mesh's base nodes.
92 _mesh_subset_base_nodes = std::make_unique<MeshLib::MeshSubset>(
93 _mesh, _base_nodes, _process_data.use_TaylorHood_elements);
94
95 // TODO move the two data members somewhere else.
96 // for extrapolation of secondary variables of stress or strain
97 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
100 std::make_unique<NumLib::LocalToGlobalIndexMap>(
101 std::move(all_mesh_subsets_single_component),
102 // by location order is needed for output
104
106 {
107 // For gas pressure, which is the first
108 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
110
111 // For capillary pressure, which is the second
112 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
113
114 // For temperature, which is the third
115 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
116
117 // For displacement.
118 std::generate_n(
119 std::back_inserter(all_mesh_subsets),
121 .get()
122 .getNumberOfGlobalComponents(),
123 [&]() { return *_mesh_subset_all_nodes; });
124
125 std::vector<int> const vec_n_components{
128
130 std::make_unique<NumLib::LocalToGlobalIndexMap>(
131 std::move(all_mesh_subsets), vec_n_components,
134 }
135 else
136 {
137 OGS_FATAL("A Staggered version of TH2M is not implemented.");
138 }
139}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:366
MeshLib::Mesh & _mesh
Definition Process.h:365
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
static constexpr int monolithic_process_id
static constexpr int n_temperature_components
static constexpr int n_displacement_components
std::vector< MeshLib::Node * > _base_nodes
static constexpr int n_gas_pressure_components
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
static constexpr int n_capillary_pressure_components
static constexpr int deformation_process_id
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ BY_LOCATION
Ordering data by spatial location.
auto & get(Tuples &... ts)
Definition Get.h:59

References NumLib::BY_LOCATION, MeshLib::getBaseNodes(), and OGS_FATAL.

◆ getDOFTable()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::TH2M::TH2MProcess< DisplacementDim >::getDOFTable ( const int process_id) const
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 344 of file TH2MProcess.cpp.

346{
347 if (hasMechanicalProcess(process_id))
348 {
350 }
351
352 // For the equation of pressure
354}
bool hasMechanicalProcess(int const process_id) const
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes

◆ getDOFTableForExtrapolatorData()

template<int DisplacementDim>
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::TH2M::TH2MProcess< DisplacementDim >::getDOFTableForExtrapolatorData ( ) const
overrideprivatevirtual

Get the address of a LocalToGlobalIndexMap, and the status of its memory. If the LocalToGlobalIndexMap is created as new in this function, the function also returns a true boolean value to let Extrapolator manage the memory by the address returned by this function.

Returns
Address of a LocalToGlobalIndexMap and its memory status.

Reimplemented from ProcessLib::Process.

Definition at line 336 of file TH2MProcess.cpp.

337{
338 const bool manage_storage = false;
339 return std::make_tuple(_local_to_global_index_map_single_component.get(),
340 manage_storage);
341}

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::TH2M::TH2MProcess< DisplacementDim >::getMatrixSpecifications ( const int process_id) const
override

Get the size and the sparse pattern of the global matrix in order to create the global matrices and vectors for the system equations of this process.

Parameters
process_idProcess ID. If the monolithic scheme is applied, process_id = 0. For the staggered scheme, process_id = 0 represents the hydraulic (H) process, while process_id = 1 represents the mechanical (M) process.
Returns
Matrix specifications including size and sparse pattern.

Definition at line 66 of file TH2MProcess.cpp.

68{
69 // For the monolithic scheme or the M process (deformation) in the staggered
70 // scheme.
72 {
73 auto const& l = *_local_to_global_index_map;
74 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
75 &l.getGhostIndices(), &this->_sparsity_pattern};
76 }
77
78 // For staggered scheme and T or H process (pressure).
80 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
81 &l.getGhostIndices(), &_sparsity_pattern_with_linear_element};
82}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:392
GlobalSparsityPattern _sparsity_pattern_with_linear_element

◆ hasMechanicalProcess()

template<int DisplacementDim>
bool ProcessLib::TH2M::TH2MProcess< DisplacementDim >::hasMechanicalProcess ( int const process_id) const
inlineprivate

Check whether the process represented by process_id is/has mechanical process. In the present implementation, the mechanical process has process_id == 1 in the staggered scheme.

Definition at line 145 of file TH2MProcess.h.

146 {
147 return _use_monolithic_scheme || process_id == deformation_process_id;
148 }

References ProcessLib::Process::_use_monolithic_scheme, and ProcessLib::TH2M::TH2MProcess< DisplacementDim >::deformation_process_id.

◆ initializeAssemblyOnSubmeshes()

template<int DisplacementDim>
std::vector< std::vector< std::string > > ProcessLib::TH2M::TH2MProcess< DisplacementDim >::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & meshes)
overrideprivatevirtual

Initializes the assembly on submeshes

Parameters
meshesthe submeshes on whom the assembly shall proceed.
Attention
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!
Returns
The names of the residuum vectors that will be assembled for each process: outer vector of size 1 for monolithic schemes and greater for staggered schemes.

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 320 of file TH2MProcess.cpp.

322{
323 INFO("TH2M process initializeSubmeshOutput().");
324 std::vector<std::vector<std::string>> residuum_names{
325 {"GasMassFlowRate", "LiquidMassFlowRate", "HeatFlowRate",
326 "NodalForces"}};
327
328 AssemblyMixin<TH2MProcess<DisplacementDim>>::initializeAssemblyOnSubmeshes(
329 meshes, residuum_names);
330
331 return residuum_names;
332}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override

References INFO().

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::initializeBoundaryConditions ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
overrideprivatevirtual

Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().

Reimplemented from ProcessLib::Process.

Definition at line 209 of file TH2MProcess.cpp.

211{
213 {
216 return;
217 }
218
219 // Staggered scheme:
220 OGS_FATAL("A Staggered version of TH2M is not implemented.");
221}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:90

References OGS_FATAL.

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 142 of file TH2MProcess.cpp.

146{
147 createLocalAssemblers<DisplacementDim>(
148 mesh.getElements(), dof_table, local_assemblers_,
149 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
151
152 auto add_secondary_variable = [&](std::string const& name,
153 int const num_components,
154 auto get_ip_values_function)
155 {
157 name,
158 makeExtrapolator(num_components, getExtrapolator(),
160 std::move(get_ip_values_function)));
161 };
162
163 ProcessLib::Reflection::addReflectedSecondaryVariables<DisplacementDim>(
166
167 //
168 // enable output of internal variables defined by material models
169 //
171 LocalAssemblerInterface<DisplacementDim>>(_process_data.solid_materials,
172 add_secondary_variable);
173
176 _process_data.solid_materials, local_assemblers_,
177 _integration_point_writer, integration_order);
178
179 _process_data.gas_pressure_interpolated =
180 MeshLib::getOrCreateMeshProperty<double>(
181 const_cast<MeshLib::Mesh&>(mesh), "gas_pressure_interpolated",
183
184 _process_data.capillary_pressure_interpolated =
185 MeshLib::getOrCreateMeshProperty<double>(
186 const_cast<MeshLib::Mesh&>(mesh), "capillary_pressure_interpolated",
188
189 _process_data.liquid_pressure_interpolated =
190 MeshLib::getOrCreateMeshProperty<double>(
191 const_cast<MeshLib::Mesh&>(mesh), "liquid_pressure_interpolated",
193
194 _process_data.temperature_interpolated =
195 MeshLib::getOrCreateMeshProperty<double>(
196 const_cast<MeshLib::Mesh&>(mesh), "temperature_interpolated",
198
201
202 // Initialize local assemblers after all variables have been set.
206}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
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)
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 executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)

References NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getProperties(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), MeshLib::Node, ProcessLib::setIPDataInitialConditions(), ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables(), and ProcessLib::Deformation::solidMaterialInternalVariablesToIntegrationPointWriter().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::TH2M::TH2MProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 59 of file TH2MProcess.cpp.

60{
61 return false;
62}

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 285 of file TH2MProcess.cpp.

289{
290 DBUG("PostTimestep TH2MProcess.");
291
295 x_prev, t, dt, process_id);
296}
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)

References DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
double const t,
double const dt,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 267 of file TH2MProcess.cpp.

270{
271 DBUG("PreTimestep TH2MProcess.");
272
273 if (hasMechanicalProcess(process_id))
274 {
278 *_local_to_global_index_map, *x[process_id], t, dt);
279 }
280
281 AssemblyMixin<TH2MProcess<DisplacementDim>>::updateActiveElements();
282}
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)

References DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ setInitialConditionsConcreteProcess()

template<int DisplacementDim>
void ProcessLib::TH2M::TH2MProcess< DisplacementDim >::setInitialConditionsConcreteProcess ( std::vector< GlobalVector * > & x,
double const t,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 224 of file TH2MProcess.cpp.

226{
227 if (process_id != 0)
228 {
229 return;
230 }
231
232 DBUG("Set initial conditions of TH2MProcess.");
233
236 local_assemblers_, getDOFTables(x.size()), x, t, process_id);
237}
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)

References DBUG(), and NumLib::SerialExecutor::executeMemberOnDereferenced().

Friends And Related Symbol Documentation

◆ AssemblyMixin< TH2MProcess< DisplacementDim > >

template<int DisplacementDim>
friend class AssemblyMixin< TH2MProcess< DisplacementDim > >
friend

Definition at line 161 of file TH2MProcess.h.

Member Data Documentation

◆ _base_nodes

template<int DisplacementDim>
std::vector<MeshLib::Node*> ProcessLib::TH2M::TH2MProcess< DisplacementDim >::_base_nodes
private

Definition at line 113 of file TH2MProcess.h.

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::TH2M::TH2MProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 121 of file TH2MProcess.h.

◆ _local_to_global_index_map_with_base_nodes

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::TH2M::TH2MProcess< DisplacementDim >::_local_to_global_index_map_with_base_nodes
private

Local to global index mapping for base nodes, which is used for linear interpolation for pressure in the staggered scheme.

Definition at line 126 of file TH2MProcess.h.

◆ _mesh_subset_base_nodes

template<int DisplacementDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::TH2M::TH2MProcess< DisplacementDim >::_mesh_subset_base_nodes
private

Definition at line 114 of file TH2MProcess.h.

◆ _process_data

template<int DisplacementDim>
TH2MProcessData<DisplacementDim> ProcessLib::TH2M::TH2MProcess< DisplacementDim >::_process_data
private

Definition at line 115 of file TH2MProcess.h.

◆ _sparsity_pattern_with_linear_element

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::TH2M::TH2MProcess< DisplacementDim >::_sparsity_pattern_with_linear_element
private

Sparsity pattern for the flow equation, and it is initialized only if the staggered scheme is used.

Definition at line 130 of file TH2MProcess.h.

◆ deformation_process_id

template<int DisplacementDim>
constexpr int ProcessLib::TH2M::TH2MProcess< DisplacementDim >::deformation_process_id = 3
staticconstexprprivate

◆ local_assemblers_

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface<DisplacementDim> > > ProcessLib::TH2M::TH2MProcess< DisplacementDim >::local_assemblers_
private

Definition at line 118 of file TH2MProcess.h.

◆ monolithic_process_id

template<int DisplacementDim>
constexpr int ProcessLib::TH2M::TH2MProcess< DisplacementDim >::monolithic_process_id = 0
staticconstexprprivate

Definition at line 150 of file TH2MProcess.h.

◆ n_capillary_pressure_components

template<int DisplacementDim>
constexpr int ProcessLib::TH2M::TH2MProcess< DisplacementDim >::n_capillary_pressure_components = 1
staticconstexprprivate

Definition at line 155 of file TH2MProcess.h.

◆ n_displacement_components

template<int DisplacementDim>
constexpr int ProcessLib::TH2M::TH2MProcess< DisplacementDim >::n_displacement_components = DisplacementDim
staticconstexprprivate

Definition at line 157 of file TH2MProcess.h.

◆ n_gas_pressure_components

template<int DisplacementDim>
constexpr int ProcessLib::TH2M::TH2MProcess< DisplacementDim >::n_gas_pressure_components = 1
staticconstexprprivate

Definition at line 154 of file TH2MProcess.h.

◆ n_temperature_components

template<int DisplacementDim>
constexpr int ProcessLib::TH2M::TH2MProcess< DisplacementDim >::n_temperature_components = 1
staticconstexprprivate

Definition at line 156 of file TH2MProcess.h.


The documentation for this class was generated from the following files: