OGS
|
The flow process is described by
\[ \phi \frac{\partial \rho_w}{\partial p} \frac{\partial p}{\partial t} S - \phi \rho_w \frac{\partial S}{\partial p_c} \frac{\partial p_c}{\partial t} - \nabla \cdot \left[\rho_w \frac{k_{\mathrm{rel}} \kappa}{\mu} \nabla \left( p + \rho_w g z \right)\right] - Q_p = 0, \]
where
Here it is assumed, that
The capillary pressure is given by
\[ p_c = \frac{\rho_w g}{\alpha} \left[S_{\mathrm{eff}}^{-\frac{1}{m}} - 1\right]^{\frac{1}{n}} \]
and the effective saturation by
\[ S_{\mathrm{eff}} = \frac{S-S_r}{S_{\max} - S_r} \]
The mass transport process is described by
\[ \phi R \frac{\partial C}{\partial t} + \nabla \cdot \left(\vec{q} C - D \nabla C \right) + \phi R \vartheta C - Q_C = 0 \]
where
For the hydrodynamic dispersion tensor the relation
\[ D = (\phi D_d + \beta_T \|\vec{q}\|) I + (\beta_L - \beta_T) \frac{\vec{q} \vec{q}^T}{\|\vec{q}\|} \]
is implemented, where \(D_d\) is the molecular diffusion coefficient, \(\beta_L\) the longitudinal dispersivity of chemical species, and \(\beta_T\) the transverse dispersivity of chemical species.
The implementation uses a monolithic approach, i.e., both processes are assembled within one global system of equations.
The advective term of the concentration equation is given by the Richards flow process, i.e., the concentration distribution depends on darcy velocity of the Richards flow process. On the other hand the concentration dependencies of the viscosity and density in the groundwater flow couples the unsaturated H process to the C process.
Definition at line 102 of file RichardsComponentTransportProcess.h.
#include <RichardsComponentTransportProcess.h>
Public Member Functions | |
RichardsComponentTransportProcess (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, RichardsComponentTransportProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme) | |
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 ¶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 |
virtual NumLib::LocalToGlobalIndexMap const & | getDOFTable (const int) const |
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 | initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override |
Process specific initialization called by initialize(). | |
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 | |
RichardsComponentTransportProcessData | _process_data |
std::vector< std::unique_ptr< RichardsComponentTransportLocalAssemblerInterface > > | _local_assemblers |
ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::RichardsComponentTransportProcess | ( | 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, | ||
RichardsComponentTransportProcessData && | process_data, | ||
SecondaryVariableCollection && | secondary_variables, | ||
bool const | use_monolithic_scheme ) |
Definition at line 21 of file RichardsComponentTransportProcess.cpp.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 67 of file RichardsComponentTransportProcess.cpp.
References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 84 of file RichardsComponentTransportProcess.cpp.
References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 39 of file RichardsComponentTransportProcess.cpp.
References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), ProcessLib::createLocalAssemblers(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::Process::getProcessVariables(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().
|
inlineoverride |
Definition at line 122 of file RichardsComponentTransportProcess.h.
|
private |
Definition at line 146 of file RichardsComponentTransportProcess.h.
Referenced by assembleConcreteProcess(), assembleWithJacobianConcreteProcess(), and initializeConcreteProcess().
|
private |
Definition at line 142 of file RichardsComponentTransportProcess.h.
Referenced by initializeConcreteProcess().