OGS
ComponentTransportProcess.cpp
Go to the documentation of this file.
1 
12 
13 #include <cassert>
14 
15 #include "BaseLib/RunTime.h"
20 #include "MathLib/LinAlg/LinAlg.h"
25 
26 namespace ProcessLib
27 {
28 namespace ComponentTransport
29 {
31  std::string name,
32  MeshLib::Mesh& mesh,
33  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
34  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
35  unsigned const integration_order,
36  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
37  process_variables,
38  ComponentTransportProcessData&& process_data,
39  SecondaryVariableCollection&& secondary_variables,
40  bool const use_monolithic_scheme,
41  std::unique_ptr<ProcessLib::SurfaceFluxData>&& surfaceflux,
42  std::unique_ptr<ChemistryLib::ChemicalSolverInterface>&&
43  chemical_solver_interface)
44  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
45  integration_order, std::move(process_variables),
46  std::move(secondary_variables), use_monolithic_scheme),
47  _process_data(std::move(process_data)),
48  _surfaceflux(std::move(surfaceflux)),
49  _chemical_solver_interface(std::move(chemical_solver_interface))
50 {
51 }
52 
54  NumLib::LocalToGlobalIndexMap const& dof_table,
55  MeshLib::Mesh const& mesh,
56  unsigned const integration_order)
57 {
58  _process_data.mesh_prop_velocity = MeshLib::getOrCreateMeshProperty<double>(
59  const_cast<MeshLib::Mesh&>(mesh), "velocity",
61 
62  const int process_id = 0;
63  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
64 
65  std::vector<std::reference_wrapper<ProcessLib::ProcessVariable>>
66  transport_process_variables;
68  {
69  for (auto pv_iter = std::next(_process_variables[process_id].begin());
70  pv_iter != _process_variables[process_id].end();
71  ++pv_iter)
72  {
73  transport_process_variables.push_back(*pv_iter);
74  }
75  }
76  else
77  {
78  for (auto pv_iter = std::next(_process_variables.begin());
79  pv_iter != _process_variables.end();
80  ++pv_iter)
81  {
82  transport_process_variables.push_back((*pv_iter)[0]);
83  }
84  }
85 
86  ProcessLib::createLocalAssemblers<LocalAssemblerData>(
87  mesh.getDimension(), mesh.getElements(), dof_table,
89  mesh.isAxiallySymmetric(), integration_order, _process_data,
90  transport_process_variables);
91 
93  {
97 
98  _chemical_solver_interface->initialize();
99  }
100 
102  "darcy_velocity",
106 }
107 
109  std::vector<GlobalVector*>& x, double const t, int const process_id)
110 {
112  {
113  return;
114  }
115 
116  if (process_id != static_cast<int>(x.size() - 1))
117  {
118  return;
119  }
120 
121  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
122 
123  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
124  dof_tables.reserve(x.size());
125  std::generate_n(std::back_inserter(dof_tables), x.size(),
126  [&]() { return _local_to_global_index_map.get(); });
127 
130  _local_assemblers, pv.getActiveElementIDs(), dof_tables, x, t);
131 
132  BaseLib::RunTime time_phreeqc;
133  time_phreeqc.start();
134 
135  _chemical_solver_interface->executeSpeciationCalculation(0. /*dt*/);
136 
138  t, _chemical_solver_interface->getIntPtProcessSolutions(), x);
139 
140  INFO("[time] Phreeqc took {:g} s.", time_phreeqc.elapsed());
141 }
142 
144  const double t, double const dt, std::vector<GlobalVector*> const& x,
145  std::vector<GlobalVector*> const& xdot, int const process_id,
147 {
148  DBUG("Assemble ComponentTransportProcess.");
149 
150  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
151 
152  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
153  dof_tables;
155  {
156  dof_tables.push_back(std::ref(*_local_to_global_index_map));
157  }
158  else
159  {
160  std::generate_n(
161  std::back_inserter(dof_tables), _process_variables.size(),
162  [&]() { return std::ref(*_local_to_global_index_map); });
163  }
164  // Call global assembler for each local assembly item.
167  pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot, process_id, M, K,
168  b);
169 }
170 
172  const double t, double const dt, std::vector<GlobalVector*> const& x,
173  std::vector<GlobalVector*> const& xdot, const double dxdot_dx,
174  const double dx_dx, int const process_id, GlobalMatrix& M, GlobalMatrix& K,
175  GlobalVector& b, GlobalMatrix& Jac)
176 {
177  DBUG("AssembleWithJacobian ComponentTransportProcess.");
178 
179  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
180  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
181  dof_table = {std::ref(*_local_to_global_index_map)};
182  // Call global assembler for each local assembly item.
185  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
186  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
187 }
188 
190  std::size_t const element_id,
191  MathLib::Point3d const& p,
192  double const t,
193  std::vector<GlobalVector*> const& x) const
194 {
195  std::vector<GlobalIndexType> indices_cache;
196  auto const r_c_indices = NumLib::getRowColumnIndices(
197  element_id, *_local_to_global_index_map, indices_cache);
198 
199  std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes{
200  x.size(), r_c_indices.rows};
201  auto const local_xs =
202  getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
203 
204  return _local_assemblers[element_id]->getFlux(p, t, local_xs);
205 }
206 
209 {
210  DBUG("Set the coupled term for the staggered scheme to local assembers.");
211 
212  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
215  setStaggeredCoupledSolutions,
217 }
218 
220  std::vector<GlobalVector*>& x, std::vector<GlobalVector*> const& x_prev,
221  double const t, double const dt, NumLib::EquationSystem& ode_sys,
222  int const process_id)
223 {
224  // todo (renchao): move chemical calculation to elsewhere.
225  if (_process_data.lookup_table && process_id == 0)
226  {
227  INFO("Update process solutions via the look-up table approach");
228  _process_data.lookup_table->lookup(x, x_prev, _mesh.getNumberOfNodes());
229 
230  return;
231  }
232 
234  {
235  return;
236  }
237 
238  // Sequential non-iterative approach applied here to split the reactive
239  // transport process into the transport stage followed by the reaction
240  // stage.
241  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
242 
243  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
244  dof_tables.reserve(x.size());
245  std::generate_n(std::back_inserter(dof_tables), x.size(),
246  [&]() { return _local_to_global_index_map.get(); });
247 
248  if (process_id == 0)
249  {
252  _local_assemblers, pv.getActiveElementIDs(), dof_tables, x, t, dt);
253 
254  BaseLib::RunTime time_phreeqc;
255  time_phreeqc.start();
256 
257  _chemical_solver_interface->setAqueousSolutionsPrevFromDumpFile();
258 
259  _chemical_solver_interface->executeSpeciationCalculation(dt);
260 
261  INFO("[time] Phreeqc took {:g} s.", time_phreeqc.elapsed());
262 
265  postSpeciationCalculation,
267 
268  return;
269  }
270 
271  auto const matrix_specification =
272  ode_sys.getMatrixSpecifications(process_id);
273 
274  std::size_t matrix_id = 0u;
276  matrix_specification, matrix_id);
278  matrix_specification, matrix_id);
279  auto& b =
281 
282  M.setZero();
283  K.setZero();
284  b.setZero();
285 
288  _local_assemblers, pv.getActiveElementIDs(), dof_tables, x, t, dt, M, K,
289  b, process_id);
290 
291  // todo (renchao): incorporate Neumann boundary condition
295 
297  matrix_specification, matrix_id);
298  auto& rhs =
300 
301  A.setZero();
302  rhs.setZero();
303 
304  // compute A
305  MathLib::LinAlg::copy(M, A);
306  MathLib::LinAlg::aypx(A, 1.0 / dt, K);
307 
308  // compute rhs
309  MathLib::LinAlg::matMult(M, *x[process_id], rhs);
310  MathLib::LinAlg::aypx(rhs, 1.0 / dt, b);
311 
312  using Tag = NumLib::NonlinearSolverTag;
314  auto& equation_system = static_cast<EqSys&>(ode_sys);
315  equation_system.applyKnownSolutionsPicard(A, rhs, *x[process_id]);
316 
317  auto& linear_solver =
319  linear_solver.solve(A, rhs, *x[process_id]);
320 
325 }
326 
328  const double t,
329  std::vector<GlobalVector*> const& integration_point_values_vectors,
330  std::vector<GlobalVector*>& nodal_values_vectors)
331 {
332  auto& extrapolator = getExtrapolator();
333  auto const extrapolatables =
336  getInterpolatedLocalSolution);
337 
338  for (unsigned transport_process_id = 0;
339  transport_process_id < integration_point_values_vectors.size();
340  ++transport_process_id)
341  {
342  auto const& pv = _process_variables[transport_process_id + 1][0].get();
343  auto const& int_pt_C =
344  integration_point_values_vectors[transport_process_id];
345 
346  extrapolator.extrapolate(pv.getNumberOfGlobalComponents(),
347  extrapolatables, t, {int_pt_C},
348  {_local_to_global_index_map.get()});
349 
350  auto const& nodal_values = extrapolator.getNodalValues();
351  MathLib::LinAlg::copy(nodal_values,
352  *nodal_values_vectors[transport_process_id + 1]);
354  *nodal_values_vectors[transport_process_id + 1]);
355  }
356 }
357 
359  double const t,
360  double const dt,
361  std::vector<GlobalVector*> const& x,
362  GlobalVector const& x_dot,
363  int const process_id)
364 {
365  if (process_id != 0)
366  {
367  return;
368  }
369 
370  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
371  dof_tables.reserve(x.size());
372  std::generate_n(std::back_inserter(dof_tables), x.size(),
373  [&]() { return _local_to_global_index_map.get(); });
374 
375  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
378  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
379  x_dot, process_id);
380 }
381 
383  std::vector<GlobalVector*> const& x,
384  const double t,
385  const double dt,
386  int const process_id)
387 {
388  if (process_id != 0)
389  {
390  return;
391  }
392 
393  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
394  dof_tables.reserve(x.size());
395  std::generate_n(std::back_inserter(dof_tables), x.size(),
396  [&]() { return _local_to_global_index_map.get(); });
397 
398  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
401  _local_assemblers, pv.getActiveElementIDs(), dof_tables, x, t, dt);
402 
403  if (!_surfaceflux) // computing the surfaceflux is optional
404  {
405  return;
406  }
407  _surfaceflux->integrate(x, t, *this, process_id, _integration_order, _mesh,
408  pv.getActiveElementIDs());
409 }
410 
411 } // namespace ComponentTransport
412 } // namespace ProcessLib
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Definition of the RunTime class.
Count the running time.
Definition: RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition: RunTime.h:42
void start()
Start the timer.
Definition: RunTime.h:32
bool solve(EigenMatrix &A, EigenVector &b, EigenVector &x)
Global vector based on Eigen vector.
Definition: EigenVector.h:26
bool isAxiallySymmetric() const
Definition: Mesh.h:126
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition: Mesh.h:89
virtual MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const =0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
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)
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)
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
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > _chemical_solver_interface
void extrapolateIntegrationPointValuesToNodes(const double t, std::vector< GlobalVector * > const &integration_point_values_vectors, std::vector< GlobalVector * > &nodal_values_vectors) override
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &x, GlobalVector const &, int const) override
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) 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 postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
void setCoupledTermForTheStaggeredSchemeToLocalAssemblers(int const process_id) override
std::vector< std::unique_ptr< ComponentTransportLocalAssemblerInterface > > _local_assemblers
std::unique_ptr< ProcessLib::SurfaceFluxData > _surfaceflux
Eigen::Vector3d getFlux(std::size_t const element_id, MathLib::Point3d const &p, double const t, std::vector< GlobalVector * > const &x) const override
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)
virtual void postTimestep(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)
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_dot, int const process_id)
unsigned getShapeFunctionOrder() const
std::vector< std::size_t > const & getActiveElementIDs() const
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
MeshLib::Mesh & _mesh
Definition: Process.h:326
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition: Process.h:360
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
CoupledSolutionsForStaggeredScheme * _coupled_solutions
Definition: Process.h:339
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
unsigned const _integration_order
Definition: Process.h:344
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
const bool _use_monolithic_scheme
Definition: Process.h:335
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
NonlinearSolverTag
Tag used to specify which nonlinear solver will be used.
Definition: Types.h:20
void finalizeAssembly(PETScMatrix &A)
Definition: LinAlg.cpp:162
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:141
void aypx(PETScVector &y, double const a, PETScVector const &x)
Definition: LinAlg.cpp:50
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
bool finalizeMatrixAssembly(MAT_T &)
NumLib::LocalToGlobalIndexMap::RowColumnIndices getRowColumnIndices(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< GlobalIndexType > &indices)
ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection > makeExtrapolatable(LocalAssemblerCollection const &local_assemblers, IntegrationPointValuesMethod integration_point_values_method)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType >> const &indices)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
ChemistryLib::ChemicalSolverInterface *const chemical_solver_interface