OGS
|
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>
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 ¶meters, 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 ¶meters, 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::Mesh & | getMesh () 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 |
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.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 130 of file StokesFlowProcess.cpp.
References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberDereferenced().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 151 of file StokesFlowProcess.cpp.
References OGS_FATAL.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 163 of file StokesFlowProcess.cpp.
References ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().
|
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.
References NumLib::BY_LOCATION, and MeshLib::getBaseNodes().
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 214 of file StokesFlowProcess.cpp.
|
override |
Definition at line 117 of file StokesFlowProcess.cpp.
|
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.
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 87 of file StokesFlowProcess.cpp.
References MeshLib::Mesh::getElements(), MeshLib::Mesh::isAxiallySymmetric(), and MeshLib::Node.
|
inlineoverride |
Definition at line 42 of file StokesFlowProcess.h.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 189 of file StokesFlowProcess.cpp.
References NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::LocalAssemblerInterface::postTimestep().
|
private |
Definition at line 85 of file StokesFlowProcess.h.
|
private |
Definition at line 91 of file StokesFlowProcess.h.
|
private |
Definition at line 94 of file StokesFlowProcess.h.
|
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.
|
private |
Definition at line 86 of file StokesFlowProcess.h.
|
private |
Definition at line 88 of file StokesFlowProcess.h.