OGS
ProcessLib::ComponentTransport::ComponentTransportProcess Class Referencefinal

Detailed Description

ComponentTransport process

Governing equations

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.

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 95 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
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 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

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 40 of file ComponentTransportProcess.cpp.

56 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
57 integration_order, std::move(process_variables),
58 std::move(secondary_variables), use_monolithic_scheme),
59 _process_data(std::move(process_data)),
60 _surfaceflux(std::move(surfaceflux)),
61 _chemical_solver_interface(std::move(chemical_solver_interface)),
62 _asm_mat_cache{is_linear, use_monolithic_scheme},
64 ls_compute_only_upon_timestep_change}
65{
66 if (ls_compute_only_upon_timestep_change)
67 {
68 // TODO move this feature to some common location for all processes.
69 if (!is_linear)
70 {
72 "Using the linear solver compute() method only upon timestep "
73 "change only makes sense for linear model equations.");
74 }
75
76 WARN(
77 "You specified that the ComponentTransport linear solver will do "
78 "the compute() step only upon timestep change. This is an expert "
79 "option. It is your responsibility to ensure that "
80 "the conditions for the correct use of this feature are met! "
81 "Otherwise OGS might compute garbage without being recognized. "
82 "There is no safety net!");
83
84 WARN(
85 "You specified that the ComponentTransport linear solver will do "
86 "the compute() step only upon timestep change. This option will "
87 "only be used by the Picard non-linear solver. The Newton-Raphson "
88 "non-linear solver will silently ignore this setting, i.e., it "
89 "won't do any harm, there, but you won't observe the speedup you "
90 "probably expect.");
91 }
92
93 auto residuum_name = [&](auto const& pv) -> std::string
94 {
95 std::string const& pv_name = pv.getName();
96 if (pv_name == "pressure")
97 {
98 return "LiquidMassFlowRate";
99 }
100 if (pv_name == "temperature")
101 {
102 return "HeatFlowRate";
103 }
104 return pv_name + "FlowRate";
105 };
106
108 {
109 int const process_id = 0;
110 for (auto const& pv : _process_variables[process_id])
111 {
113 mesh, residuum_name(pv.get()), MeshLib::MeshItemType::Node, 1));
114 }
115 }
116 else
117 {
118 for (auto const& pv : _process_variables)
119 {
121 mesh, residuum_name(pv[0].get()), MeshLib::MeshItemType::Node,
122 1));
123 }
124 }
125}
#define OGS_FATAL(...)
Definition Error.h:26
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:42
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > _chemical_solver_interface
std::vector< MeshLib::PropertyVector< double > * > _residua
std::unique_ptr< ProcessLib::SurfaceFluxData > _surfaceflux
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
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:406
const bool _use_monolithic_scheme
Definition Process.h:385
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
auto & get(Tuples &... ts)
Definition Get.h:59

References ProcessLib::Process::Process(), _asm_mat_cache, _chemical_solver_interface, _ls_compute_only_upon_timestep_change, _process_data, ProcessLib::Process::_process_variables, _residua, _surfaceflux, ProcessLib::Process::_use_monolithic_scheme, MeshLib::getOrCreateMeshProperty(), ProcessLib::Process::name, MeshLib::Node, OGS_FATAL, and WARN().

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 240 of file ComponentTransportProcess.cpp.

244{
245 DBUG("Assemble ComponentTransportProcess.");
246
247 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
249 {
250 dof_tables.push_back(_local_to_global_index_map.get());
251 }
252 else
253 {
254 std::generate_n(std::back_inserter(dof_tables),
255 _process_variables.size(),
256 [&]() { return _local_to_global_index_map.get(); });
257 }
258
259 _asm_mat_cache.assemble(t, dt, x, x_prev, process_id, &M, &K, &b,
262
263 // TODO (naumov) What about temperature? A test with FCT would reveal any
264 // problems.
265 if (process_id == _process_data.hydraulic_process_id)
266 {
267 return;
268 }
269 auto const matrix_specification = getMatrixSpecifications(process_id);
270
272 x_prev, process_id,
273 matrix_specification, M, K, b);
274}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::unique_ptr< ComponentTransportLocalAssemblerInterface > > _local_assemblers
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
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:210
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 _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, NumLib::computeFluxCorrectedTransport(), DBUG(), ProcessLib::Process::getActiveElementIDs(), and ProcessLib::Process::getMatrixSpecifications().

Referenced by preOutputConcreteProcess().

◆ 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 276 of file ComponentTransportProcess.cpp.

280{
281 DBUG("AssembleWithJacobian ComponentTransportProcess.");
282
284 {
285 OGS_FATAL(
286 "The AssembleWithJacobian() for ComponentTransportProcess is "
287 "implemented only for staggered scheme.");
288 }
289
290 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
291
292 std::generate_n(std::back_inserter(dof_tables), _process_variables.size(),
293 [&]() { return _local_to_global_index_map.get(); });
294
295 // Call global assembler for each local assembly item.
298 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
299 process_id, &b, &Jac);
300
301 // b is the negated residumm used in the Newton's method.
302 // Here negating b is to recover the primitive residuum.
304 *_residua[process_id],
305 std::negate<double>());
306}
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)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)
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::Process::_process_variables, _residua, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::Process::getActiveElementIDs(), 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 440 of file ComponentTransportProcess.cpp.

446{
447 if (process_id != 0)
448 {
449 return;
450 }
451
452 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
453 dof_tables.reserve(x.size());
454 std::generate_n(std::back_inserter(dof_tables), x.size(),
455 [&]() { return _local_to_global_index_map.get(); });
456
459 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
460 process_id);
461
463 {
464 return;
465 }
466
468 &ComponentTransportLocalAssemblerInterface::
469 computeReactionRelatedSecondaryVariable,
470 _local_assemblers, _chemical_solver_interface->activeElementIDs());
471}
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)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

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

◆ 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 308 of file ComponentTransportProcess.cpp.

313{
314 std::vector<GlobalIndexType> indices_cache;
315 auto const r_c_indices = NumLib::getRowColumnIndices(
316 element_id, *_local_to_global_index_map, indices_cache);
317
318 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes{
319 x.size(), r_c_indices.rows};
320 auto const local_xs =
321 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
322
323 return _local_assemblers[element_id]->getFlux(p, t, local_xs);
324}
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 _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::getCoupledLocalSolutions(), and NumLib::getRowColumnIndices().

◆ 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 127 of file ComponentTransportProcess.cpp.

131{
132 int const mesh_space_dimension = _process_data.mesh_space_dimension;
133
135 const_cast<MeshLib::Mesh&>(mesh), "velocity",
136 MeshLib::MeshItemType::Cell, mesh_space_dimension);
137
139 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
141
142 std::vector<std::reference_wrapper<ProcessLib::ProcessVariable>>
143 transport_process_variables;
145 {
146 const int process_id = 0;
147 for (auto pv_iter = std::next(_process_variables[process_id].begin());
148 pv_iter != _process_variables[process_id].end();
149 ++pv_iter)
150 {
151 transport_process_variables.push_back(*pv_iter);
152 }
153 }
154 else
155 {
156 // All process variables but the pressure and optionally the temperature
157 // are transport variables.
158 for (auto const& pv :
160 ranges::views::drop(_process_data.isothermal ? 1 : 2))
161 {
162 transport_process_variables.push_back(pv[0]);
163 }
164 }
165
167 mesh_space_dimension, mesh.getElements(), dof_table, _local_assemblers,
168 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
169 _process_data, transport_process_variables);
170
172 {
175 _local_assemblers, _chemical_solver_interface->activeElementIDs());
176
177 _chemical_solver_interface->initialize();
178 }
179
180 _secondary_variables.addSecondaryVariable(
181 "darcy_velocity",
183 mesh_space_dimension, getExtrapolator(), _local_assemblers,
185
186 for (std::size_t component_id = 0;
187 component_id < transport_process_variables.size();
188 ++component_id)
189 {
190 auto const& pv = transport_process_variables[component_id].get();
191
192 auto integration_point_values_accessor = [component_id](
193 ComponentTransportLocalAssemblerInterface const& loc_asm,
194 const double t,
195 std::vector<GlobalVector*> const& x,
196 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
197 std::vector<double>& cache) -> auto const&
198 {
199 return loc_asm.getIntPtMolarFlux(t, x, dof_tables, cache,
200 component_id);
201 };
202
203 _secondary_variables.addSecondaryVariable(
204 pv.getName() + "Flux",
205 makeExtrapolator(mesh_space_dimension, getExtrapolator(),
207 integration_point_values_accessor));
208 }
209}
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:376
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 _chemical_solver_interface, _local_assemblers, _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::getIntPtMolarFlux(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystemID().

◆ isLinear()

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

Definition at line 120 of file ComponentTransportProcess.h.

120{ return _asm_mat_cache.isLinear(); }

References _asm_mat_cache.

◆ 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 473 of file ComponentTransportProcess.cpp.

479{
480 if (process_id != 0)
481 {
482 return;
483 }
484
485 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
486 dof_tables.reserve(x.size());
487 std::generate_n(std::back_inserter(dof_tables), x.size(),
488 [&]() { return _local_to_global_index_map.get(); });
489
492 _local_assemblers, getActiveElementIDs(), dof_tables, x, x_prev, t, dt,
493 process_id);
494
495 if (!_surfaceflux) // computing the surfaceflux is optional
496 {
497 return;
498 }
499 _surfaceflux->integrate(x, t, *this, process_id, _integration_order, _mesh,
501}
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:371
unsigned const _integration_order
Definition Process.h:390

References ProcessLib::Process::_integration_order, _local_assemblers, ProcessLib::Process::_mesh, _surfaceflux, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), 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 503 of file ComponentTransportProcess.cpp.

509{
510 auto const matrix_specification = getMatrixSpecifications(process_id);
511
512 auto M = MathLib::MatrixVectorTraits<GlobalMatrix>::newInstance(
513 matrix_specification);
514 auto K = MathLib::MatrixVectorTraits<GlobalMatrix>::newInstance(
515 matrix_specification);
516 auto b = MathLib::MatrixVectorTraits<GlobalVector>::newInstance(
517 matrix_specification);
518
519 M->setZero();
520 K->setZero();
521 b->setZero();
522
523 assembleConcreteProcess(t, dt, x, x_prev, process_id, *M, *K, *b);
524
525 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
527 {
528 dof_tables.push_back(_local_to_global_index_map.get());
529 }
530 else
531 {
532 std::generate_n(std::back_inserter(dof_tables),
533 _process_variables.size(),
534 [&]() { return _local_to_global_index_map.get(); });
535 }
536
537 BaseLib::RunTime time_residuum;
538 time_residuum.start();
539
541 {
542 auto const residuum =
543 computeResiduum(dt, *x[0], *x_prev[0], *M, *K, *b);
544 for (std::size_t variable_id = 0; variable_id < _residua.size();
545 ++variable_id)
546 {
548 residuum, variable_id, *dof_tables[0], *_residua[variable_id],
549 std::negate<double>());
550 }
551 }
552 else
553 {
554 auto const residuum = computeResiduum(dt, *x[process_id],
555 *x_prev[process_id], *M, *K, *b);
556 transformVariableFromGlobalVector(residuum, 0, *dof_tables[process_id],
557 *_residua[process_id],
558 std::negate<double>());
559 }
560
561 INFO("[time] Computing residuum flow rates took {:g} s",
562 time_residuum.elapsed());
563}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:36
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
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
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)

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

◆ 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 211 of file ComponentTransportProcess.cpp.

213{
215 {
216 return;
217 }
218
219 if (process_id != static_cast<int>(x.size() - 1))
220 {
221 return;
222 }
223
224 std::for_each(
225 x.begin(), x.end(),
226 [](auto const process_solution)
227 { MathLib::LinAlg::setLocalAccessibleVector(*process_solution); });
228
229 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
230 dof_tables.reserve(x.size());
231 std::generate_n(std::back_inserter(dof_tables), x.size(),
232 [&]() { return _local_to_global_index_map.get(); });
233
237 dof_tables, x, t);
238}
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, _local_assemblers, NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::initializeChemicalSystem().

◆ shouldLinearSolverComputeOnlyUponTimestepChange()

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

◆ 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 326 of file ComponentTransportProcess.cpp.

330{
331 // todo (renchao): move chemical calculation to elsewhere.
332 if (_process_data.lookup_table && process_id == 0)
333 {
334 INFO("Update process solutions via the look-up table approach");
335 _process_data.lookup_table->lookup(x, x_prev, _mesh.getNumberOfNodes());
336
337 return;
338 }
339
341 {
342 return;
343 }
344
345 // Sequential non-iterative approach applied here to split the reactive
346 // transport process into the transport stage followed by the reaction
347 // stage.
348 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
349 dof_tables.reserve(x.size());
350 std::generate_n(std::back_inserter(dof_tables), x.size(),
351 [&]() { return _local_to_global_index_map.get(); });
352
353 if (process_id == 0)
354 {
358 dof_tables, x, t, dt);
359
360 BaseLib::RunTime time_phreeqc;
361 time_phreeqc.start();
362
363 _chemical_solver_interface->setAqueousSolutionsPrevFromDumpFile();
364
365 _chemical_solver_interface->executeSpeciationCalculation(dt);
366
367 INFO("[time] Phreeqc took {:g} s.", time_phreeqc.elapsed());
368
370 &ComponentTransportLocalAssemblerInterface::
371 postSpeciationCalculation,
373 t, dt);
374
375 return;
376 }
377
378 auto const matrix_specification = getMatrixSpecifications(process_id);
379
380 std::size_t matrix_id = 0u;
382 matrix_specification, matrix_id);
384 matrix_specification, matrix_id);
385 auto& b =
387
388 M.setZero();
389 K.setZero();
390 b.setZero();
391
394 _local_assemblers, getActiveElementIDs(), dof_tables, x, t, dt, M, K, b,
395 process_id);
396
397 // todo (renchao): incorporate Neumann boundary condition
401
403 matrix_specification, matrix_id);
404 auto& rhs =
406
407 A.setZero();
408 rhs.setZero();
409
412
413 // compute A
415 MathLib::LinAlg::aypx(A, 1.0 / dt, K);
416
417 // compute rhs
418 MathLib::LinAlg::matMult(M, *x[process_id], rhs);
419 MathLib::LinAlg::aypx(rhs, 1.0 / dt, b);
420
421 using Tag = NumLib::NonlinearSolverTag;
422 using EqSys = NumLib::NonlinearSystem<Tag::Picard>;
423 auto& equation_system = static_cast<EqSys&>(ode_sys);
424 equation_system.applyKnownSolutionsPicard(
425 A, rhs, *x[process_id],
427
428 auto& linear_solver =
429 _process_data.chemical_solver_interface->linear_solver;
431 linear_solver.solve(A, rhs, *x[process_id]);
432
438}
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:20
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:149
void aypx(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50
@ COMPLETE_MATRIX_UPDATE
Both A and b fully updated.
Definition LinAlgEnums.h:43
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, _local_assemblers, 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(), MathLib::LinAlg::matMult(), NumLib::GlobalMatrixProvider::provider, NumLib::GlobalVectorProvider::provider, ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystem(), MathLib::LinAlg::setLocalAccessibleVector(), and BaseLib::RunTime::start().

Member Data Documentation

◆ _asm_mat_cache

AssembledMatrixCache ProcessLib::ComponentTransport::ComponentTransportProcess::_asm_mat_cache
private

◆ _chemical_solver_interface

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

◆ _local_assemblers

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

◆ _residua

std::vector<MeshLib::PropertyVector<double>*> ProcessLib::ComponentTransport::ComponentTransportProcess::_residua
private

◆ _surfaceflux

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

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