OGS
|
The flow process is described by
\[ \phi \frac{\partial \rho}{\partial p} \frac{\partial p}{\partial t} + \phi \frac{\partial \rho}{\partial C} \frac{\partial C}{\partial t} - \nabla \cdot \left[\frac{\kappa}{\mu(C)} \rho \nabla \left( p + \rho g z \right)\right] + Q_p = 0, \]
where the storage \(S\) has been substituted by \(\phi \frac{\partial \rho}{\partial p}\), \(\phi\) is the porosity, \(C\) is the concentration, \(p\) is the pressure, \(\kappa\) is permeability, \(\mu\) is viscosity of the fluid, \(\rho\) is the density of the fluid, and \(g\) is the gravitational acceleration.
The mass transport process is described by
\[ \phi R C \frac{\partial \rho}{\partial p} \frac{\partial p}{\partial t} + \phi R \left(\rho + C \frac{\partial \rho}{\partial C}\right) \frac{\partial C}{\partial t} - \nabla \cdot \left[\frac{\kappa}{\mu(C)} \rho C \nabla \left( p + \rho g z \right) + \rho D \nabla C\right] + Q_C + R \vartheta \phi \rho C = 0, \]
where \(R\) is the retardation factor, \(\vec{q} = -\frac{\kappa}{\mu(C)} \nabla \left( p + \rho 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.
The advective term of the concentration equation is given by the confined groundwater flow process, i.e., the concentration distribution depends on Darcy velocity of the groundwater flow process. On the other hand the concentration dependencies of the viscosity and density in the groundwater flow couples the H process to the C process.
Definition at line 95 of file ComponentTransportProcess.h.
#include <ComponentTransportProcess.h>
Public Member Functions | |
ComponentTransportProcess (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, ComponentTransportProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, std::unique_ptr< ProcessLib::SurfaceFluxData > &&surfaceflux, std::unique_ptr< ChemistryLib::ChemicalSolverInterface > &&chemical_solver_interface, bool const is_linear, bool const ls_compute_only_upon_timestep_change) | |
Eigen::Vector3d | getFlux (std::size_t const element_id, MathLib::Point3d const &p, double const t, std::vector< GlobalVector * > const &x) const override |
void | solveReactionEquation (std::vector< GlobalVector * > &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, NumLib::EquationSystem &ode_sys, int const process_id) 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 |
bool | shouldLinearSolverComputeOnlyUponTimestepChange () 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 ¶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 |
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 | 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 | preOutputConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override |
Private Attributes | |
ComponentTransportProcessData | _process_data |
std::vector< std::unique_ptr< ComponentTransportLocalAssemblerInterface > > | _local_assemblers |
std::unique_ptr< ProcessLib::SurfaceFluxData > | _surfaceflux |
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > | _chemical_solver_interface |
std::vector< MeshLib::PropertyVector< double > * > | _residua |
AssembledMatrixCache | _asm_mat_cache |
bool const | _ls_compute_only_upon_timestep_change |
ProcessLib::ComponentTransport::ComponentTransportProcess::ComponentTransportProcess | ( | 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, | ||
ComponentTransportProcessData && | process_data, | ||
SecondaryVariableCollection && | secondary_variables, | ||
bool const | use_monolithic_scheme, | ||
std::unique_ptr< ProcessLib::SurfaceFluxData > && | surfaceflux, | ||
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > && | chemical_solver_interface, | ||
bool const | is_linear, | ||
bool const | ls_compute_only_upon_timestep_change ) |
Definition at line 40 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_process_variables, _residua, ProcessLib::Process::_use_monolithic_scheme, MeshLib::getOrCreateMeshProperty(), MeshLib::Node, OGS_FATAL, and WARN().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 240 of file ComponentTransportProcess.cpp.
References _asm_mat_cache, ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_process_variables, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::AssembledMatrixCache::assemble(), NumLib::computeFluxCorrectedTransport(), DBUG(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getMatrixSpecifications(), ProcessLib::ComponentTransport::ComponentTransportProcessData::hydraulic_process_id, and ProcessLib::ComponentTransport::ComponentTransportProcessData::stabilizer.
Referenced by preOutputConcreteProcess().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 276 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_process_variables, _residua, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::Process::getActiveElementIDs(), and OGS_FATAL.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 438 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, _local_assemblers, ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::Process::getActiveElementIDs().
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 308 of file ComponentTransportProcess.cpp.
References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::getCoupledLocalSolutions(), and NumLib::getRowColumnIndices().
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 127 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, _local_assemblers, _process_data, ProcessLib::Process::_process_variables, ProcessLib::Process::_secondary_variables, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), MeshLib::Cell, ProcessLib::createLocalAssemblers(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtDarcyVelocity(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtMolarFlux(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::ComponentTransport::ComponentTransportProcessData::isothermal, ProcessLib::makeExtrapolator(), ProcessLib::ComponentTransport::ComponentTransportProcessData::mesh_prop_porosity, ProcessLib::ComponentTransport::ComponentTransportProcessData::mesh_prop_velocity, ProcessLib::ComponentTransport::ComponentTransportProcessData::mesh_space_dimension, and ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystemID().
|
inlineoverride |
Definition at line 120 of file ComponentTransportProcess.h.
References _asm_mat_cache, and ProcessLib::AssembledMatrixCache::isLinear().
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 471 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_integration_order, _local_assemblers, ProcessLib::Process::_mesh, _surfaceflux, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postTimestep().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 501 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_process_variables, _residua, ProcessLib::Process::_use_monolithic_scheme, assembleConcreteProcess(), ProcessLib::computeResiduum(), BaseLib::RunTime::elapsed(), ProcessLib::Process::getMatrixSpecifications(), INFO(), and BaseLib::RunTime::start().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 211 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, _local_assemblers, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::initializeChemicalSystem().
|
inlineoverride |
Definition at line 145 of file ComponentTransportProcess.h.
References _ls_compute_only_upon_timestep_change.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 326 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, _local_assemblers, ProcessLib::Process::_mesh, _process_data, ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::assembleReactionEquation(), MathLib::LinAlg::aypx(), ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, MathLib::LinAlg::copy(), BaseLib::RunTime::elapsed(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), MathLib::finalizeMatrixAssembly(), MathLib::finalizeVectorAssembly(), ProcessLib::Process::getActiveElementIDs(), NumLib::MatrixProvider::getMatrix(), ProcessLib::Process::getMatrixSpecifications(), MeshLib::Mesh::getNumberOfNodes(), NumLib::VectorProvider::getVector(), INFO(), ChemistryLib::ChemicalSolverInterface::linear_solver, ProcessLib::ComponentTransport::ComponentTransportProcessData::lookup_table, MathLib::LinAlg::matMult(), NumLib::GlobalMatrixProvider::provider, NumLib::GlobalVectorProvider::provider, NumLib::MatrixProvider::releaseMatrix(), NumLib::VectorProvider::releaseVector(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystem(), MathLib::LinAlg::setLocalAccessibleVector(), MathLib::EigenVector::setZero(), and BaseLib::RunTime::start().
|
private |
Definition at line 188 of file ComponentTransportProcess.h.
Referenced by assembleConcreteProcess(), and isLinear().
|
private |
Definition at line 184 of file ComponentTransportProcess.h.
Referenced by computeSecondaryVariableConcrete(), initializeConcreteProcess(), setInitialConditionsConcreteProcess(), and solveReactionEquation().
|
private |
Definition at line 179 of file ComponentTransportProcess.h.
Referenced by assembleConcreteProcess(), assembleWithJacobianConcreteProcess(), computeSecondaryVariableConcrete(), getFlux(), initializeConcreteProcess(), postTimestepConcreteProcess(), setInitialConditionsConcreteProcess(), and solveReactionEquation().
|
private |
Definition at line 190 of file ComponentTransportProcess.h.
Referenced by shouldLinearSolverComputeOnlyUponTimestepChange().
|
private |
Definition at line 176 of file ComponentTransportProcess.h.
Referenced by assembleConcreteProcess(), initializeConcreteProcess(), and solveReactionEquation().
|
private |
Definition at line 186 of file ComponentTransportProcess.h.
Referenced by ComponentTransportProcess(), assembleWithJacobianConcreteProcess(), and preOutputConcreteProcess().
|
private |
Definition at line 181 of file ComponentTransportProcess.h.
Referenced by postTimestepConcreteProcess().