OGS
ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim > Class Template Referencefinal

Detailed Description

template<int GlobalDim>
class ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >

The Stokes equations are reduced from the Navier-Stokes equations, and include momentum equations and a continuity equation describing the flow in the free-flow region.

Definition at line 25 of file StokesFlowProcess.h.

#include <StokesFlowProcess.h>

Inheritance diagram for ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >:
[legend]

Public Member Functions

 StokesFlowProcess (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, StokesFlowProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
 
bool isLinear () const override
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int) const override
 
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const override
 
void computeSecondaryVariableConcrete (double const, double const, std::vector< GlobalVector * > const &x, GlobalVector const &, int 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
 
- 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 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 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
 

Private Attributes

std::vector< MeshLib::Node * > _base_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
 
StokesFlowProcessData _process_data
 
std::vector< std::unique_ptr< StokesFlowLocalAssemblerInterface > > _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
 

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)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
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

◆ StokesFlowProcess()

template<int GlobalDim>
ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::StokesFlowProcess ( 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,
StokesFlowProcessData && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme )

Definition at line 24 of file StokesFlowProcess.cpp.

35 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
36 integration_order, std::move(process_variables),
37 std::move(secondary_variables), use_monolithic_scheme),
38 _process_data(std::move(process_data))
39{
40}
std::string const name
Definition Process.h:362
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

Member Function Documentation

◆ assembleConcreteProcess()

template<int GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::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 130 of file StokesFlowProcess.cpp.

134{
135 DBUG("Assemble StokesFlowProcess.");
136
137 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
139 {
140 dof_tables.push_back(_local_to_global_index_map.get());
141 }
142
143 // Call global assembler for each local assembly item.
146 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, &M, &K,
147 &b);
148}
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
const bool _use_monolithic_scheme
Definition Process.h:379
std::vector< std::unique_ptr< StokesFlowLocalAssemblerInterface > > _local_assemblers
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 GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::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 151 of file StokesFlowProcess.cpp.

156{
157 OGS_FATAL(
158 "Assembly of Jacobian matrix has not yet been implemented for "
159 "StokesFlowProcess.");
160}
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ computeSecondaryVariableConcrete()

template<int GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::computeSecondaryVariableConcrete ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 163 of file StokesFlowProcess.cpp.

169{
170 if (process_id != 0)
171 {
172 return;
173 }
174
175 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
176 dof_tables.reserve(x.size());
178 {
179 dof_tables.push_back(_local_to_global_index_map.get());
180 }
181
184 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
185 process_id);
186}
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 ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ constructDofTable()

template<int GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::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 43 of file StokesFlowProcess.cpp.

44{
45 // Create single component dof in every of the mesh's nodes.
47 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
48 // Create single component dof in the mesh's base nodes.
51 std::make_unique<MeshLib::MeshSubset>(_mesh, _base_nodes);
52
53 // TODO move the two data members somewhere else.
54 // for extrapolation of secondary variables
55 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
58 std::make_unique<NumLib::LocalToGlobalIndexMap>(
59 std::move(all_mesh_subsets_single_component),
60 // by location order is needed for output
62
64 {
65 // For vector variables, in this case liquid velocity.
66 int const process_id = 0;
67 std::vector<MeshLib::MeshSubset> all_mesh_subsets(
68 getProcessVariables(process_id)[0]
69 .get()
70 .getNumberOfGlobalComponents(),
72
73 // For scalar variables, in this case pressure.
74 all_mesh_subsets.push_back(*_mesh_subset_base_nodes);
75
76 std::vector<int> const vec_n_components{GlobalDim, 1};
77
79 std::make_unique<NumLib::LocalToGlobalIndexMap>(
80 std::move(all_mesh_subsets), vec_n_components,
83 }
84}
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::vector< MeshLib::Node * > _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.
auto & get(Tuples &... ts)
Definition Get.h:59

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

◆ getDOFTable()

template<int GlobalDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::getDOFTable ( const int ) const
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 214 of file StokesFlowProcess.cpp.

216{
218 {
220 }
221}

◆ getMatrixSpecifications()

template<int GlobalDim>
MathLib::MatrixSpecifications ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::getMatrixSpecifications ( const int ) const
override

Definition at line 117 of file StokesFlowProcess.cpp.

119{
121 {
122 auto const& l = *_local_to_global_index_map;
123
124 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
125 &l.getGhostIndices(), &this->_sparsity_pattern};
126 }
127}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:392

◆ initializeBoundaryConditions()

template<int GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::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 104 of file StokesFlowProcess.cpp.

106{
108 {
109 int const process_id = 0;
111 *_local_to_global_index_map, process_id, media);
112 }
113}
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 GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::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 87 of file StokesFlowProcess.cpp.

91{
93 mesh.getElements(), dof_table, _local_assemblers,
94 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
96
99 const_cast<MeshLib::Mesh&>(mesh), "pressure_interpolated",
101}
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
void createLocalAssemblersStokes(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)
MeshLib::PropertyVector< double > * pressure_interpolated

References ProcessLib::createLocalAssemblersStokes(), MeshLib::Mesh::getElements(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::isAxiallySymmetric(), and MeshLib::Node.

◆ isLinear()

template<int GlobalDim>
bool ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::isLinear ( ) const
inlineoverride

Definition at line 42 of file StokesFlowProcess.h.

42{ return false; }

◆ postTimestepConcreteProcess()

template<int GlobalDim>
void ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double dt,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 189 of file StokesFlowProcess.cpp.

195{
196 if (process_id != 0)
197 {
198 return;
199 }
200
201 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
202 dof_tables.reserve(x.size());
204 {
205 dof_tables.push_back(_local_to_global_index_map.get());
206 }
207
210 getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
211}
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 NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::LocalAssemblerInterface::postTimestep().

Member Data Documentation

◆ _base_nodes

template<int GlobalDim>
std::vector<MeshLib::Node*> ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::_base_nodes
private

Definition at line 85 of file StokesFlowProcess.h.

◆ _local_assemblers

template<int GlobalDim>
std::vector<std::unique_ptr<StokesFlowLocalAssemblerInterface> > ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::_local_assemblers
private

Definition at line 91 of file StokesFlowProcess.h.

◆ _local_to_global_index_map_single_component

template<int GlobalDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::_local_to_global_index_map_single_component
private

Definition at line 94 of file StokesFlowProcess.h.

◆ _local_to_global_index_map_with_base_nodes

template<int GlobalDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::_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 99 of file StokesFlowProcess.h.

◆ _mesh_subset_base_nodes

template<int GlobalDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::_mesh_subset_base_nodes
private

Definition at line 86 of file StokesFlowProcess.h.

◆ _process_data

template<int GlobalDim>
StokesFlowProcessData ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::_process_data
private

Definition at line 88 of file StokesFlowProcess.h.


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