OGS
ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess Class Referencefinal

Detailed Description

RichardsComponentTransport process

Governing Equations

Richards Flow

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

  • \(\phi\) is the porosity,
  • \(S\) is the saturation,
  • \(p\) is the pressure,
  • \(k_{\mathrm{rel}}\) is the relative permeability (depending on \(S\)),
  • \(\kappa\) is the specific permeability,
  • \(\mu\) is viscosity of the fluid,
  • \(\rho_w\) is the mass density of the fluid, and
  • \(g\) is the gravitational acceleration.

Here it is assumed, that

  • the porosity is constant and
  • the pressure of the gas phase is zero.

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} \]

Mass Transport

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

  • \(R\) is the retardation factor,
  • \(C\) is the concentration,
  • \(\vec{q} = \frac{k_{\mathrm{rel}} \kappa}{\mu(C)} \nabla \left( p + \rho_w g z \right)\) is the Darcy velocity,
  • \(D\) is the hydrodynamic dispersion tensor,
  • \(\vartheta\) is the decay rate.

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.

Process Coupling

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.

Note
At the moment there is not any coupling by source or sink terms, i.e., the coupling is implemented only by density changes due to concentration changes in the buoyance term in the groundwater flow. This coupling schema is referred to as the Boussinesq approximation.

Definition at line 102 of file RichardsComponentTransportProcess.h.

#include <RichardsComponentTransportProcess.h>

Inheritance diagram for ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess:
[legend]
Collaboration diagram for ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess:
[legend]

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 &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)
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
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
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 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

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)
virtual void constructDofTable ()
void constructMonolithicProcessDofTable ()
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
void setReleaseNodalForces (GlobalVector const *r_neq, int const process_id) 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

◆ RichardsComponentTransportProcess()

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.

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

References ProcessLib::Process::Process(), _process_data, and ProcessLib::Process::name.

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::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 67 of file RichardsComponentTransportProcess.cpp.

71{
72 DBUG("Assemble RichardsComponentTransportProcess.");
73
74 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
76
77 // Call global assembler for each local assembly item.
80 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
81 &b);
82}
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:383
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:374
std::vector< std::unique_ptr< RichardsComponentTransportLocalAssemblerInterface > > _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::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::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 84 of file RichardsComponentTransportProcess.cpp.

88{
89 DBUG("AssembleWithJacobian RichardsComponentTransportProcess.");
90
91 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
93
94 // Call global assembler for each local assembly item.
97 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
98 process_id, &b, &Jac);
99}
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)

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().

◆ initializeConcreteProcess()

void ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::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 39 of file RichardsComponentTransportProcess.cpp.

43{
44 const int monolithic_process_id = 0;
45 ProcessLib::ProcessVariable const& pv =
46 getProcessVariables(monolithic_process_id)[0];
47
49 mesh.getDimension(), mesh.getElements(), dof_table, _local_assemblers,
50 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
51 _process_data, pv);
52
53 _secondary_variables.addSecondaryVariable(
54 "darcy_velocity",
55 makeExtrapolator(mesh.getDimension(), getExtrapolator(),
57 &RichardsComponentTransportLocalAssemblerInterface::
58 getIntPtDarcyVelocity));
59
60 _secondary_variables.addSecondaryVariable(
61 "saturation",
63 &RichardsComponentTransportLocalAssemblerInterface::
64 getIntPtSaturation));
65}
SecondaryVariableCollection _secondary_variables
Definition Process.h:376
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables() const
Definition Process.h:156
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void createLocalAssemblers(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)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)

References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::createLocalAssemblers(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::Process::getProcessVariables(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

bool ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::isLinear ( ) const
inlineoverride

Definition at line 122 of file RichardsComponentTransportProcess.h.

122{ return false; }

Member Data Documentation

◆ _local_assemblers

std::vector< std::unique_ptr<RichardsComponentTransportLocalAssemblerInterface> > ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::_local_assemblers
private

◆ _process_data

RichardsComponentTransportProcessData ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::_process_data
private

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