OGS
ProcessLib::ComponentTransport::ComponentTransportProcess Class Referencefinal

Detailed Description

ComponentTransport process

Governing equations

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.

Process Coupling

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.

Note
At the moment there is not any coupling by source or sink terms, i.e., the coupling is implemented only by density and viscosity changes due to concentration changes as well as by the temporal derivatives of each variable.

Definition at line 88 of file ComponentTransportProcess.h.

#include <ComponentTransportProcess.h>

Inheritance diagram for ProcessLib::ComponentTransport::ComponentTransportProcess:
[legend]
Collaboration diagram for ProcessLib::ComponentTransport::ComponentTransportProcess:
[legend]

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 &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)
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 &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
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 >

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

◆ 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.

50 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
51 integration_order, std::move(process_variables),
52 std::move(secondary_variables), use_monolithic_scheme),
54 use_monolithic_scheme},
55 _process_data(std::move(process_data)),
56 _surfaceflux(std::move(surfaceflux)),
57 _chemical_solver_interface(std::move(chemical_solver_interface)),
59 ls_compute_only_upon_timestep_change}
60{
61 this->_jacobian_assembler->checkPerturbationSize(_process_variables.size());
62
63 if (ls_compute_only_upon_timestep_change)
64 {
65 // TODO move this feature to some common location for all processes.
66 if (!is_linear)
67 {
69 "Using the linear solver compute() method only upon timestep "
70 "change only makes sense for linear model equations.");
71 }
72
73 WARN(
74 "You specified that the ComponentTransport linear solver will do "
75 "the compute() step only upon timestep change. This is an expert "
76 "option. It is your responsibility to ensure that "
77 "the conditions for the correct use of this feature are met! "
78 "Otherwise OGS might compute garbage without being recognized. "
79 "There is no safety net!");
80
81 WARN(
82 "You specified that the ComponentTransport linear solver will do "
83 "the compute() step only upon timestep change. This option will "
84 "only be used by the Picard non-linear solver. The Newton-Raphson "
85 "non-linear solver will silently ignore this setting, i.e., it "
86 "won't do any harm, there, but you won't observe the speedup you "
87 "probably expect.");
88 }
89
91 {
92 // For numerical Jacobian:
93 std::vector<int> non_deformation_component_ids(
94 _process_variables.size());
95 std::iota(non_deformation_component_ids.begin(),
96 non_deformation_component_ids.end(), 0);
97 this->_jacobian_assembler->setNonDeformationComponentIDsNoSizeCheck(
98 non_deformation_component_ids);
99 }
100}
#define OGS_FATAL(...)
Definition Error.h:19
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > _chemical_solver_interface
std::unique_ptr< ProcessLib::SurfaceFluxData > _surfaceflux
std::string const name
Definition Process.h:361
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:399
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
const bool _use_monolithic_scheme
Definition Process.h:378

References ComponentTransportProcess(), ProcessLib::Process::Process(), ProcessLib::Process::_jacobian_assembler, and ProcessLib::Process::name.

Referenced by ComponentTransportProcess(), and AssemblyMixin< ComponentTransportProcess >.

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::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 222 of file ComponentTransportProcess.cpp.

226{
227 DBUG("Assemble ComponentTransportProcess.");
228
230 {
231 this->_jacobian_assembler->setNonDeformationComponentIDsNoSizeCheck(
232 {process_id});
233 }
234
236 process_id, M, K, b);
237
238 // TODO (naumov) What about temperature? A test with FCT would reveal any
239 // problems.
240 if (process_id == _process_data.hydraulic_process_id)
241 {
242 return;
243 }
244 auto const matrix_specification = getMatrixSpecifications(process_id);
245
247 x_prev, process_id,
248 matrix_specification, M, K, b);
249}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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)
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:203
void computeFluxCorrectedTransport(NumericalStabilization const &stabilizer, const double t, const double dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, const MathLib::MatrixSpecifications &matrix_specification, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)

References ProcessLib::Process::_jacobian_assembler, _process_data, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::AssemblyMixin< Process >::assemble(), NumLib::computeFluxCorrectedTransport(), DBUG(), and ProcessLib::Process::getMatrixSpecifications().

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::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 251 of file ComponentTransportProcess.cpp.

255{
256 DBUG("AssembleWithJacobian ComponentTransportProcess.");
257
259 {
260 OGS_FATAL(
261 "The AssembleWithJacobian() for ComponentTransportProcess is "
262 "implemented only for staggered scheme.");
263 }
264
266 t, dt, x, x_prev, process_id, b, Jac);
267}
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)

References ProcessLib::Process::_use_monolithic_scheme, ProcessLib::AssemblyMixin< Process >::assembleWithJacobian(), DBUG(), and OGS_FATAL.

◆ computeSecondaryVariableConcrete()

void ProcessLib::ComponentTransport::ComponentTransportProcess::computeSecondaryVariableConcrete ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 401 of file ComponentTransportProcess.cpp.

407{
408 if (process_id != 0)
409 {
410 return;
411 }
412
413 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
414 dof_tables.reserve(x.size());
415 std::generate_n(std::back_inserter(dof_tables), x.size(),
416 [&]() { return _local_to_global_index_map.get(); });
417
420 local_assemblers_, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
421 process_id);
422
424 {
425 return;
426 }
427
429 &ComponentTransportLocalAssemblerInterface::
430 computeReactionRelatedSecondaryVariable,
431 local_assemblers_, _chemical_solver_interface->activeElementIDs());
432}
std::vector< std::unique_ptr< ComponentTransportLocalAssemblerInterface > > local_assemblers_
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References _chemical_solver_interface, ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), and local_assemblers_.

◆ getFlux()

Eigen::Vector3d ProcessLib::ComponentTransport::ComponentTransportProcess::getFlux ( std::size_t const element_id,
MathLib::Point3d const & p,
double const t,
std::vector< GlobalVector * > const & x ) const
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 269 of file ComponentTransportProcess.cpp.

274{
275 std::vector<GlobalIndexType> indices_cache;
276 auto const r_c_indices = NumLib::getRowColumnIndices(
277 element_id, *_local_to_global_index_map, indices_cache);
278
279 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes{
280 x.size(), r_c_indices.rows};
281 auto const local_xs =
282 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
283
284 return local_assemblers_[element_id]->getFlux(p, t, local_xs);
285}
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
NumLib::LocalToGlobalIndexMap::RowColumnIndices getRowColumnIndices(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< GlobalIndexType > &indices)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)

References ProcessLib::Process::_local_to_global_index_map, ProcessLib::getCoupledLocalSolutions(), NumLib::getRowColumnIndices(), and local_assemblers_.

◆ initializeAssemblyOnSubmeshes()

std::vector< std::vector< std::string > > ProcessLib::ComponentTransport::ComponentTransportProcess::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & meshes)
overridevirtual

Initializes the assembly on submeshes

Parameters
meshesthe submeshes on whom the assembly shall proceed.
Attention
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!
Returns
The names of the residuum vectors that will be assembled for each process: outer vector of size 1 for monolithic schemes and greater for staggered schemes.

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 435 of file ComponentTransportProcess.cpp.

437{
438 DBUG("CT process initializeSubmeshOutput().");
439
440 namespace R = ranges;
441 namespace RV = ranges::views;
442
443 auto get_residuum_name =
444 [](ProcessLib::ProcessVariable const& pv) -> std::string
445 {
446 std::string const& pv_name = pv.getName();
447 if (pv_name == "pressure")
448 {
449 return "LiquidMassFlowRate";
450 }
451 if (pv_name == "temperature")
452 {
453 return "HeatFlowRate";
454 }
455 return pv_name + "FlowRate";
456 };
457
458 auto get_residuum_names = [get_residuum_name](auto const& pvs)
459 { return pvs | RV::transform(get_residuum_name); };
460
461 auto residuum_names = _process_variables |
462 RV::transform(get_residuum_names) |
463 R::to<std::vector<std::vector<std::string>>>();
464
466 meshes, residuum_names);
467
468 return residuum_names;
469}
void initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)

References ProcessLib::Process::_process_variables, DBUG(), and ProcessLib::AssemblyMixin< Process >::initializeAssemblyOnSubmeshes().

◆ initializeConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::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 102 of file ComponentTransportProcess.cpp.

106{
107 int const mesh_space_dimension = _process_data.mesh_space_dimension;
108
110 const_cast<MeshLib::Mesh&>(mesh), "velocity",
111 MeshLib::MeshItemType::Cell, mesh_space_dimension);
112
114 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
116
117 std::vector<std::reference_wrapper<ProcessLib::ProcessVariable>>
118 transport_process_variables;
120 {
121 const int process_id = 0;
122 for (auto pv_iter = std::next(_process_variables[process_id].begin());
123 pv_iter != _process_variables[process_id].end();
124 ++pv_iter)
125 {
126 transport_process_variables.push_back(*pv_iter);
127 }
128 }
129 else
130 {
131 // All process variables but the pressure and optionally the temperature
132 // are transport variables.
133 for (auto const& pv :
135 ranges::views::drop(_process_data.isothermal ? 1 : 2))
136 {
137 transport_process_variables.push_back(pv[0]);
138 }
139 }
140
142 mesh_space_dimension, mesh.getElements(), dof_table, local_assemblers_,
143 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
144 _process_data, transport_process_variables);
145
147 {
150 local_assemblers_, _chemical_solver_interface->activeElementIDs());
151
152 _chemical_solver_interface->initialize();
153 }
154
155 _secondary_variables.addSecondaryVariable(
156 "liquid_density",
160
161 _secondary_variables.addSecondaryVariable(
162 "darcy_velocity",
164 mesh_space_dimension, getExtrapolator(), local_assemblers_,
166
167 for (std::size_t component_id = 0;
168 component_id < transport_process_variables.size();
169 ++component_id)
170 {
171 auto const& pv = transport_process_variables[component_id].get();
172
173 auto integration_point_values_accessor =
174 [component_id](
175 ComponentTransportLocalAssemblerInterface const& loc_asm,
176 const double t,
177 std::vector<GlobalVector*> const& x,
178 std::vector<NumLib::LocalToGlobalIndexMap const*> const&
179 dof_tables,
180 std::vector<double>& cache) -> auto const&
181 {
182 return loc_asm.getIntPtMolarFlux(t, x, dof_tables, cache,
183 component_id);
184 };
185
186 _secondary_variables.addSecondaryVariable(
187 pv.getName() + "Flux",
188 makeExtrapolator(mesh_space_dimension, getExtrapolator(),
190 integration_point_values_accessor));
191 }
192}
virtual std::vector< double > const & getIntPtLiquidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
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 _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().

◆ isLinear()

bool ProcessLib::ComponentTransport::ComponentTransportProcess::isLinear ( ) const
inlineoverride

◆ postTimestepConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double dt,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 483 of file ComponentTransportProcess.cpp.

489{
490 if (process_id != 0)
491 {
492 return;
493 }
494
495 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
496 dof_tables.reserve(x.size());
497 std::generate_n(std::back_inserter(dof_tables), x.size(),
498 [&]() { return _local_to_global_index_map.get(); });
499
502 local_assemblers_, getActiveElementIDs(), dof_tables, x, x_prev, t, dt,
503 process_id);
504
505 if (!_surfaceflux) // computing the surfaceflux is optional
506 {
507 return;
508 }
509 _surfaceflux->integrate(x, t, *this, process_id, _integration_order, _mesh,
511}
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
MeshLib::Mesh & _mesh
Definition Process.h:364
unsigned const _integration_order
Definition Process.h:383

References ProcessLib::Process::_integration_order, ProcessLib::Process::_mesh, _surfaceflux, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), local_assemblers_, and ProcessLib::LocalAssemblerInterface::postTimestep().

◆ preOutputConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::preOutputConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 513 of file ComponentTransportProcess.cpp.

519{
521 process_id);
522}
void preOutput(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)

References ProcessLib::AssemblyMixin< Process >::preOutput().

◆ preTimestepConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & ,
const double ,
const double ,
const int process_id )
overridevirtual

◆ setInitialConditionsConcreteProcess()

void ProcessLib::ComponentTransport::ComponentTransportProcess::setInitialConditionsConcreteProcess ( std::vector< GlobalVector * > & x,
double const t,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 194 of file ComponentTransportProcess.cpp.

196{
198 {
199 return;
200 }
201
202 if (process_id != static_cast<int>(x.size() - 1))
203 {
204 return;
205 }
206
207 std::for_each(
208 x.begin(), x.end(), [](auto const process_solution)
209 { MathLib::LinAlg::setLocalAccessibleVector(*process_solution); });
210
211 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
212 dof_tables.reserve(x.size());
213 std::generate_n(std::back_inserter(dof_tables), x.size(),
214 [&]() { return _local_to_global_index_map.get(); });
215
219 dof_tables, x, t);
220}
void initializeChemicalSystem(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t)

References _chemical_solver_interface, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::initializeChemicalSystem(), and local_assemblers_.

◆ shouldLinearSolverComputeOnlyUponTimestepChange()

bool ProcessLib::ComponentTransport::ComponentTransportProcess::shouldLinearSolverComputeOnlyUponTimestepChange ( ) const
inlineoverride

Definition at line 154 of file ComponentTransportProcess.h.

155 {
156 // TODO unify for all processes?
158 }

References _ls_compute_only_upon_timestep_change.

◆ solveReactionEquation()

void ProcessLib::ComponentTransport::ComponentTransportProcess::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 )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 287 of file ComponentTransportProcess.cpp.

291{
292 // todo (renchao): move chemical calculation to elsewhere.
293 if (_process_data.lookup_table && process_id == 0)
294 {
295 INFO("Update process solutions via the look-up table approach");
296 _process_data.lookup_table->lookup(x, x_prev, _mesh.getNumberOfNodes());
297
298 return;
299 }
300
302 {
303 return;
304 }
305
306 // Sequential non-iterative approach applied here to split the reactive
307 // transport process into the transport stage followed by the reaction
308 // stage.
309 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
310 dof_tables.reserve(x.size());
311 std::generate_n(std::back_inserter(dof_tables), x.size(),
312 [&]() { return _local_to_global_index_map.get(); });
313
314 if (process_id == 0)
315 {
319 dof_tables, x, t, dt);
320
321 BaseLib::RunTime time_phreeqc;
322 time_phreeqc.start();
323
324 _chemical_solver_interface->setAqueousSolutionsPrevFromDumpFile();
325
326 _chemical_solver_interface->executeSpeciationCalculation(dt);
327
328 INFO("[time] Phreeqc took {:g} s.", time_phreeqc.elapsed());
329
331 &ComponentTransportLocalAssemblerInterface::
332 postSpeciationCalculation,
334 t, dt);
335
336 return;
337 }
338
339 auto const matrix_specification = getMatrixSpecifications(process_id);
340
341 std::size_t matrix_id = 0u;
343 matrix_specification, matrix_id);
345 matrix_specification, matrix_id);
346 auto& b =
348
349 M.setZero();
350 K.setZero();
351 b.setZero();
352
355 local_assemblers_, getActiveElementIDs(), dof_tables, x, t, dt, M, K, b,
356 process_id);
357
358 // todo (renchao): incorporate Neumann boundary condition
362
364 matrix_specification, matrix_id);
365 auto& rhs =
367
368 A.setZero();
369 rhs.setZero();
370
373
374 // compute A
376 MathLib::LinAlg::aypx(A, 1.0 / dt, K);
377
378 // compute rhs
379 MathLib::LinAlg::matMult(M, *x[process_id], rhs);
380 MathLib::LinAlg::aypx(rhs, 1.0 / dt, b);
381
382 using Tag = NumLib::NonlinearSolverTag;
383 using EqSys = NumLib::NonlinearSystem<Tag::Picard>;
384 auto& equation_system = static_cast<EqSys&>(ode_sys);
385 equation_system.applyKnownSolutionsPicard(
386 A, rhs, *x[process_id],
388
389 auto& linear_solver =
390 _process_data.chemical_solver_interface->linear_solver;
392 linear_solver.solve(A, rhs, *x[process_id]);
393
399}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:31
void start()
Start the timer.
Definition RunTime.h:21
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
void setChemicalSystem(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
void assembleReactionEquation(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, int const process_id)
NonlinearSolverTag
Tag used to specify which nonlinear solver will be used.
Definition Types.h:13
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:20
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:142
void aypx(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:43
@ COMPLETE_MATRIX_UPDATE
Both A and b fully updated.
Definition LinAlgEnums.h:34
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
bool finalizeMatrixAssembly(MAT_T &)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider

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

◆ AssemblyMixin< ComponentTransportProcess >

Member Data Documentation

◆ _chemical_solver_interface

std::unique_ptr<ChemistryLib::ChemicalSolverInterface> ProcessLib::ComponentTransport::ComponentTransportProcess::_chemical_solver_interface
private

◆ _ls_compute_only_upon_timestep_change

bool const ProcessLib::ComponentTransport::ComponentTransportProcess::_ls_compute_only_upon_timestep_change
private

◆ _process_data

ComponentTransportProcessData ProcessLib::ComponentTransport::ComponentTransportProcess::_process_data
private

◆ _surfaceflux

std::unique_ptr<ProcessLib::SurfaceFluxData> ProcessLib::ComponentTransport::ComponentTransportProcess::_surfaceflux
private

Definition at line 191 of file ComponentTransportProcess.h.

Referenced by postTimestepConcreteProcess().

◆ local_assemblers_

std::vector<std::unique_ptr<ComponentTransportLocalAssemblerInterface> > ProcessLib::ComponentTransport::ComponentTransportProcess::local_assemblers_
private

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