OGS
ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >

Linear kinematics poro-mechanical/biphasic (fluid-solid mixture) model.

The mixture momentum balance and the mixture mass balance are solved under fully saturated conditions.

Definition at line 26 of file RichardsMechanicsProcess.h.

#include <RichardsMechanicsProcess.h>

Inheritance diagram for ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >:
[legend]

Public Member Functions

 RichardsMechanicsProcess (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, RichardsMechanicsProcessData< 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 std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Private Types

using LocalAssemblerIF = LocalAssemblerInterface<DisplacementDim>
 

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 postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, const int process_id) override
 
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
 
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const override
 
bool hasMechanicalProcess (int const process_id) const
 

Private Attributes

std::vector< MeshLib::Node * > _base_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
 
RichardsMechanicsProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerIF > > _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
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 
MeshLib::PropertyVector< double > * _hydraulic_flow = nullptr
 

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
 

Member Typedef Documentation

◆ LocalAssemblerIF

template<int DisplacementDim>
using ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::LocalAssemblerIF = LocalAssemblerInterface<DisplacementDim>
private

Definition at line 65 of file RichardsMechanicsProcess.h.

Constructor & Destructor Documentation

◆ RichardsMechanicsProcess()

template<int DisplacementDim>
ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::RichardsMechanicsProcess ( 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,
RichardsMechanicsProcessData< DisplacementDim > && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme )

Definition at line 31 of file RichardsMechanicsProcess.cpp.

42 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
43 integration_order, std::move(process_variables),
44 std::move(secondary_variables), use_monolithic_scheme),
45 _process_data(std::move(process_data))
46{
48 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
49
51 mesh, "MassFlowRate", MeshLib::MeshItemType::Node, 1);
52
55 _integration_point_writer, integration_order,
57}
std::string const name
Definition Process.h:362
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
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)
Definition Process.cpp:44
std::vector< std::unique_ptr< LocalAssemblerIF > > _local_assemblers
RichardsMechanicsProcessData< DisplacementDim > _process_data
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
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::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_hydraulic_flow, ProcessLib::Process::_integration_point_writer, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_local_assemblers, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_nodal_forces, ProcessLib::Reflection::addReflectedIntegrationPointWriters(), MeshLib::getOrCreateMeshProperty(), ProcessLib::RichardsMechanics::LocalAssemblerInterface< DisplacementDim >::getReflectionDataForOutput(), and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 273 of file RichardsMechanicsProcess.cpp.

277{
278 DBUG("Assemble the equations for RichardsMechanics");
279
280 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
282
283 // Call global assembler for each local assembly item.
286 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
287 &b);
288}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const 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)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 291 of file RichardsMechanicsProcess.cpp.

296{
297 // For the monolithic scheme
299 {
300 DBUG(
301 "Assemble the Jacobian of RichardsMechanics for the monolithic"
302 " scheme.");
303 }
304 else
305 {
306 // For the staggered scheme
307 if (process_id == 0)
308 {
309 DBUG(
310 "Assemble the Jacobian equations of liquid fluid process in "
311 "RichardsMechanics for the staggered scheme.");
312 }
313 else
314 {
315 DBUG(
316 "Assemble the Jacobian equations of mechanical process in "
317 "RichardsMechanics for the staggered scheme.");
318 }
319 }
320
321 auto const dof_tables = getDOFTables(x.size());
324 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
325 process_id, &b, &Jac);
326
327 auto copyRhs = [&](int const variable_id, auto& output_vector)
328 {
330 {
331 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
332 output_vector,
333 std::negate<double>());
334 }
335 else
336 {
337 transformVariableFromGlobalVector(b, 0, *dof_tables[process_id],
338 output_vector,
339 std::negate<double>());
340 }
341 };
342 if (_use_monolithic_scheme || process_id == 0)
343 {
344 copyRhs(0, *_hydraulic_flow);
345 }
346 if (_use_monolithic_scheme || process_id == 1)
347 {
348 copyRhs(1, *_nodal_forces);
349 }
350}
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
const bool _use_monolithic_scheme
Definition Process.h:379
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, 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)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().

◆ computeSecondaryVariableConcrete()

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

Reimplemented from ProcessLib::Process.

Definition at line 370 of file RichardsMechanicsProcess.cpp.

375{
376 if (process_id != 0)
377 {
378 return;
379 }
380
381 DBUG("Compute the secondary variables for RichardsMechanicsProcess.");
382
385 getActiveElementIDs(), getDOFTables(x.size()), t, dt, x, x_prev,
386 process_id);
387}
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)
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::RichardsMechanics::RichardsMechanicsProcess< 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 86 of file RichardsMechanicsProcess.cpp.

87{
88 // Create single component dof in every of the mesh's nodes.
90 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
91 // Create single component dof in the mesh's base nodes.
94 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
95
96 // TODO move the two data members somewhere else.
97 // for extrapolation of secondary variables of stress or strain
98 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
101 std::make_unique<NumLib::LocalToGlobalIndexMap>(
102 std::move(all_mesh_subsets_single_component),
103 // by location order is needed for output
105
107 {
108 // For pressure, which is the first
109 std::vector<MeshLib::MeshSubset> all_mesh_subsets{
111
112 // For displacement.
113 const int monolithic_process_id = 0;
114 std::generate_n(std::back_inserter(all_mesh_subsets),
115 getProcessVariables(monolithic_process_id)[1]
116 .get()
117 .getNumberOfGlobalComponents(),
118 [&]() { return *_mesh_subset_all_nodes; });
119
120 std::vector<int> const vec_n_components{1, DisplacementDim};
122 std::make_unique<NumLib::LocalToGlobalIndexMap>(
123 std::move(all_mesh_subsets), vec_n_components,
126 }
127 else
128 {
129 // For displacement equation.
130 const int process_id = 1;
131 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
132 std::generate_n(std::back_inserter(all_mesh_subsets),
133 getProcessVariables(process_id)[0]
134 .get()
135 .getNumberOfGlobalComponents(),
136 [&]() { return *_mesh_subset_all_nodes; });
137
138 std::vector<int> const vec_n_components{DisplacementDim};
140 std::make_unique<NumLib::LocalToGlobalIndexMap>(
141 std::move(all_mesh_subsets), vec_n_components,
143
144 // For pressure equation.
145 // Collect the mesh subsets with base nodes in a vector.
146 std::vector<MeshLib::MeshSubset> all_mesh_subsets_base_nodes{
149 std::make_unique<NumLib::LocalToGlobalIndexMap>(
150 std::move(all_mesh_subsets_base_nodes),
151 // by location order is needed for output
153
156
159 }
160}
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::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_with_base_nodes
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
std::vector< Node * > getBaseNodes(std::vector< Element * > const &elements)
Definition Utils.h:26
@ 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.
auto & get(Tuples &... ts)
Definition Get.h:59

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

◆ getDOFTable()

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

Reimplemented from ProcessLib::Process.

Definition at line 400 of file RichardsMechanicsProcess.cpp.

402{
403 if (hasMechanicalProcess(process_id))
404 {
406 }
407
408 // For the equation of pressure
410}

◆ getDOFTableForExtrapolatorData()

template<int DisplacementDim>
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 391 of file RichardsMechanicsProcess.cpp.

392{
393 const bool manage_storage = false;
394 return std::make_tuple(_local_to_global_index_map_single_component.get(),
395 manage_storage);
396}

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 67 of file RichardsMechanicsProcess.cpp.

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

◆ hasMechanicalProcess()

template<int DisplacementDim>
bool ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 133 of file RichardsMechanicsProcess.h.

134 {
135 return _use_monolithic_scheme || process_id == 1;
136 }

References ProcessLib::Process::_use_monolithic_scheme.

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 230 of file RichardsMechanicsProcess.cpp.

232{
234 {
235 const int monolithic_process_id = 0;
237 *_local_to_global_index_map, monolithic_process_id, media);
238 return;
239 }
240
241 // Staggered scheme:
242 // for the equations of pressure
243 const int hydraulic_process_id = 0;
245 *_local_to_global_index_map_with_base_nodes, hydraulic_process_id,
246 media);
247
248 // for the equations of deformation.
249 const int mechanical_process_id = 1;
251 *_local_to_global_index_map, mechanical_process_id, media);
252}
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

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 163 of file RichardsMechanicsProcess.cpp.

167{
169 RichardsMechanicsLocalAssembler>(
170 mesh.getElements(), dof_table, _local_assemblers,
171 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
173
177
178 auto add_secondary_variable = [&](std::string const& name,
179 int const num_components,
180 auto get_ip_values_function)
181 {
183 name,
184 makeExtrapolator(num_components, getExtrapolator(),
186 std::move(get_ip_values_function)));
187 };
188
189 //
190 // enable output of internal variables defined by material models
191 //
193 LocalAssemblerIF>(_process_data.solid_materials,
194 add_secondary_variable);
195
198 _process_data.solid_materials, _local_assemblers,
199 _integration_point_writer, integration_order);
200
202 const_cast<MeshLib::Mesh&>(mesh), "saturation_avg",
204
206 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
208
210 const_cast<MeshLib::Mesh&>(mesh), "stress_avg",
213 DisplacementDim>::RowsAtCompileTime);
214
215 _process_data.pressure_interpolated =
217 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
219
222
223 // Initialize local assemblers after all variables have been set.
227}
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
LocalAssemblerInterface< DisplacementDim > LocalAssemblerIF
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
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)
void addReflectedSecondaryVariables(ReflData const &reflection_data, SecondaryVariableCollection &secondary_variables, NumLib::Extrapolator &extrapolator, 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)
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 executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)

References ProcessLib::Reflection::addReflectedSecondaryVariables(), MeshLib::Cell, ProcessLib::createLocalAssemblersHM(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), MeshLib::getOrCreateMeshProperty(), 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::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 60 of file RichardsMechanicsProcess.cpp.

61{
62 return false;
63}

◆ postTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 353 of file RichardsMechanicsProcess.cpp.

357{
358 if (hasMechanicalProcess(process_id))
359 {
360 DBUG("PostTimestep RichardsMechanicsProcess.");
361
364 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
365 process_id);
366 }
367}
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().

◆ setInitialConditionsConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 255 of file RichardsMechanicsProcess.cpp.

259{
260 if (process_id != 0)
261 {
262 return;
263 }
264
265 DBUG("SetInitialConditions RichardsMechanicsProcess.");
266
269 getActiveElementIDs(), getDOFTables(x.size()), x, t, process_id);
270}
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::executeSelectedMemberOnDereferenced().

Member Data Documentation

◆ _base_nodes

template<int DisplacementDim>
std::vector<MeshLib::Node*> ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_base_nodes
private

Definition at line 102 of file RichardsMechanicsProcess.h.

◆ _hydraulic_flow

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_hydraulic_flow = nullptr
private

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerIF> > ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_local_assemblers
private

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 109 of file RichardsMechanicsProcess.h.

◆ _local_to_global_index_map_with_base_nodes

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 114 of file RichardsMechanicsProcess.h.

◆ _mesh_subset_base_nodes

template<int DisplacementDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_mesh_subset_base_nodes
private

Definition at line 103 of file RichardsMechanicsProcess.h.

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

template<int DisplacementDim>
RichardsMechanicsProcessData<DisplacementDim> ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::_process_data
private

Definition at line 104 of file RichardsMechanicsProcess.h.

◆ _sparsity_pattern_with_linear_element

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::RichardsMechanics::RichardsMechanicsProcess< 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 118 of file RichardsMechanicsProcess.h.


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