![]() |
OGS
|
|
The flow process is described by
\[\left(\rho \beta_s + \phi \frac{\partial \rho}{\partial p} \right)\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, \]
\(\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, \(g\) is the gravitational acceleration, and \(\beta_s = \frac{\partial \phi}{\partial p}\).
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 88 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 |
| std::vector< std::vector< std::string > > | initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override |
| void | preTimestepConcreteProcess (std::vector< GlobalVector * > const &, const double, const double, const int) 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 | ~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 Member Functions inherited from ProcessLib::AssemblyMixin< ComponentTransportProcess > | |
| void | initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names) |
| void | updateActiveElements () |
| void | assemble (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, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false) |
| void | assembleWithJacobian (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false) |
| void | preOutput (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) |
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 |
| bool const | _ls_compute_only_upon_timestep_change |
Friends | |
| class | AssemblyMixin< ComponentTransportProcess > |
| 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 34 of file ComponentTransportProcess.cpp.
References ComponentTransportProcess(), ProcessLib::Process::Process(), ProcessLib::Process::_jacobian_assembler, and ProcessLib::Process::name.
Referenced by ComponentTransportProcess(), and AssemblyMixin< ComponentTransportProcess >.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 222 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_jacobian_assembler, _process_data, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::AssemblyMixin< Process >::assemble(), NumLib::computeFluxCorrectedTransport(), DBUG(), and ProcessLib::Process::getMatrixSpecifications().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 251 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_use_monolithic_scheme, ProcessLib::AssemblyMixin< Process >::assembleWithJacobian(), DBUG(), and OGS_FATAL.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 401 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), and local_assemblers_.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 269 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_local_to_global_index_map, ProcessLib::getCoupledLocalSolutions(), NumLib::getRowColumnIndices(), and local_assemblers_.
|
overridevirtual |
Initializes the assembly on submeshes
| meshes | the submeshes on whom the assembly shall proceed. |
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!Reimplemented from ProcessLib::SubmeshAssemblySupport.
Definition at line 435 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_process_variables, DBUG(), and ProcessLib::AssemblyMixin< Process >::initializeAssemblyOnSubmeshes().
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 102 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, _process_data, ProcessLib::Process::_process_variables, ProcessLib::Process::_secondary_variables, ProcessLib::Process::_use_monolithic_scheme, MeshLib::Cell, ProcessLib::createLocalAssemblers(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtDarcyVelocity(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtLiquidDensity(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtMolarFlux(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::isAxiallySymmetric(), local_assemblers_, ProcessLib::makeExtrapolator(), and ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystemID().
|
inlineoverride |
Definition at line 117 of file ComponentTransportProcess.h.
References ProcessLib::AssemblyMixinBase::isLinear().
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 483 of file ComponentTransportProcess.cpp.
References ProcessLib::Process::_integration_order, ProcessLib::Process::_mesh, _surfaceflux, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), local_assemblers_, and ProcessLib::LocalAssemblerInterface::postTimestep().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 513 of file ComponentTransportProcess.cpp.
References ProcessLib::AssemblyMixin< Process >::preOutput().
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 471 of file ComponentTransportProcess.cpp.
References ProcessLib::AssemblyMixin< Process >::updateActiveElements().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 194 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::initializeChemicalSystem(), and local_assemblers_.
|
inlineoverride |
Definition at line 154 of file ComponentTransportProcess.h.
References _ls_compute_only_upon_timestep_change.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 287 of file ComponentTransportProcess.cpp.
References _chemical_solver_interface, ProcessLib::Process::_mesh, _process_data, ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::assembleReactionEquation(), MathLib::LinAlg::aypx(), MathLib::COMPLETE_MATRIX_UPDATE, MathLib::LinAlg::copy(), BaseLib::RunTime::elapsed(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), MathLib::finalizeMatrixAssembly(), MathLib::finalizeVectorAssembly(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getMatrixSpecifications(), INFO(), local_assemblers_, MathLib::LinAlg::matMult(), NumLib::GlobalMatrixProvider::provider, NumLib::GlobalVectorProvider::provider, ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystem(), MathLib::LinAlg::setLocalAccessibleVector(), and BaseLib::RunTime::start().
|
friend |
Definition at line 88 of file ComponentTransportProcess.h.
References ComponentTransportProcess(), ProcessLib::Process::Process(), and ProcessLib::Process::name.
|
private |
Definition at line 194 of file ComponentTransportProcess.h.
Referenced by computeSecondaryVariableConcrete(), initializeConcreteProcess(), setInitialConditionsConcreteProcess(), and solveReactionEquation().
|
private |
Definition at line 196 of file ComponentTransportProcess.h.
Referenced by shouldLinearSolverComputeOnlyUponTimestepChange().
|
private |
Definition at line 186 of file ComponentTransportProcess.h.
Referenced by assembleConcreteProcess(), initializeConcreteProcess(), and solveReactionEquation().
|
private |
Definition at line 191 of file ComponentTransportProcess.h.
Referenced by postTimestepConcreteProcess().
|
private |
Definition at line 189 of file ComponentTransportProcess.h.
Referenced by computeSecondaryVariableConcrete(), getFlux(), initializeConcreteProcess(), postTimestepConcreteProcess(), setInitialConditionsConcreteProcess(), and solveReactionEquation().