OGS
NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard > Class Referencefinal

Detailed Description

Time discretized first order implicit quasi-linear ODE; to be solved using the Picard fixpoint iteration method for resolving nonlinearities.

See also
ODESystemTag::FirstOrderImplicitQuasilinear

Definition at line 161 of file TimeDiscretizedODESystem.h.

#include <TimeDiscretizedODESystem.h>

Inheritance diagram for NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >:
[legend]
Collaboration diagram for NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >:
[legend]

Public Types

using ODE = ODESystem< ODETag, NonlinearSolverTag::Picard >
 The type of ODE. More...
 
using MatTrans = MatrixTranslator< ODETag >
 
using TimeDisc = TimeDiscretization
 A shortcut for a general time discretization scheme. More...
 

Public Member Functions

 TimeDiscretizedODESystem (const int process_id, ODE &ode, TimeDisc &time_discretization)
 Constructs a new instance. More...
 
 ~TimeDiscretizedODESystem () override
 
void assemble (std::vector< GlobalVector * > const &x_new_timestep, std::vector< GlobalVector * > const &x_prev, int const process_id) override
 
void getA (GlobalMatrix &A) const override
 
void getRhs (GlobalVector const &x_prev, GlobalVector &rhs) const override
 
void computeKnownSolutions (GlobalVector const &x, int const process_id) override
 
void applyKnownSolutions (GlobalVector &x) const override
 
void applyKnownSolutionsPicard (GlobalMatrix &A, GlobalVector &rhs, GlobalVector &x) const override
 
bool isLinear () const override
 
void preIteration (const unsigned iter, GlobalVector const &x) override
 
IterationResult postIteration (GlobalVector const &x) override
 
TimeDiscgetTimeDiscretization () override
 Exposes the used time discretization scheme. More...
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 

Static Public Attributes

static const ODESystemTag ODETag
 A tag indicating the type of ODE. More...
 

Private Types

using Index = MathLib::MatrixVectorTraits< GlobalMatrix >::Index
 

Private Attributes

ODE_ode
 ode the ODE being wrapped More...
 
TimeDisc_time_disc
 the time discretization to being used More...
 
std::unique_ptr< MatTrans_mat_trans
 the object used to compute the matrix/vector for the nonlinear solver More...
 
std::vector< NumLib::IndexValueVector< Index > > const * _known_solutions
 stores precomputed values for known solutions More...
 
GlobalMatrix_M
 Matrix \( M \). More...
 
GlobalMatrix_K
 Matrix \( K \). More...
 
GlobalVector_b
 Matrix \( b \). More...
 
std::size_t _M_id = 0u
 ID of the _M matrix. More...
 
std::size_t _K_id = 0u
 ID of the _K matrix. More...
 
std::size_t _b_id = 0u
 ID of the _b vector. More...
 
std::vector< std::size_t > _xdot_ids
 

Member Typedef Documentation

◆ Index

◆ MatTrans

The auxiliary class that computes the matrix/vector used by the nonlinear solver.

Definition at line 174 of file TimeDiscretizedODESystem.h.

◆ ODE

◆ TimeDisc

Constructor & Destructor Documentation

◆ TimeDiscretizedODESystem()

Constructs a new instance.

Definition at line 189 of file TimeDiscretizedODESystem.cpp.

192  : _ode(ode),
193  _time_disc(time_discretization),
194  _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
195 {
197  ode.getMatrixSpecifications(process_id), _M_id);
199  ode.getMatrixSpecifications(process_id), _K_id);
201  ode.getMatrixSpecifications(process_id), _b_id);
202 }
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
std::unique_ptr< MatTrans > _mat_trans
the object used to compute the matrix/vector for the nonlinear solver
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider

References NumLib::MatrixProvider::getMatrix(), NumLib::VectorProvider::getVector(), NumLib::GlobalVectorProvider::provider, and NumLib::GlobalMatrixProvider::provider.

◆ ~TimeDiscretizedODESystem()

Member Function Documentation

◆ applyKnownSolutions()

Definition at line 262 of file TimeDiscretizedODESystem.cpp.

263 {
265 }
std::vector< NumLib::IndexValueVector< Index > > const * _known_solutions
stores precomputed values for known solutions
void applyKnownSolutions(std::vector< Solutions > const *const known_solutions, Vector &x)
Applies known solutions to the solution vector x.

References detail::applyKnownSolutions().

◆ applyKnownSolutionsPicard()

Definition at line 268 of file TimeDiscretizedODESystem.cpp.

272 {
273  if (!_known_solutions)
274  {
275  return;
276  }
277 
279  std::vector<IndexType> ids;
280  std::vector<double> values;
281  for (auto const& bc : *_known_solutions)
282  {
283  ids.insert(end(ids), begin(bc.ids), end(bc.ids));
284  values.insert(end(values), begin(bc.values), end(bc.values));
285  }
286  MathLib::applyKnownSolution(A, rhs, x, ids, values);
287 }
void applyKnownSolution(EigenMatrix &A, EigenVector &b, EigenVector &, const std::vector< EigenMatrix::IndexType > &vec_knownX_id, const std::vector< double > &vec_knownX_x, double)
Definition: EigenTools.cpp:17

References MathLib::applyKnownSolution().

◆ assemble()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::assemble ( std::vector< GlobalVector * > const &  x_new_timestep,
std::vector< GlobalVector * > const &  x_prev,
int const  process_id 
)
override

Definition at line 214 of file TimeDiscretizedODESystem.cpp.

218 {
219  namespace LinAlg = MathLib::LinAlg;
220 
221  auto const t = _time_disc.getCurrentTime();
222  auto const dt = _time_disc.getCurrentTimeIncrement();
223  auto const& x_curr = *x_new_timestep[process_id];
224  std::vector<GlobalVector*> xdot(x_new_timestep.size());
225  _xdot_ids.resize(x_new_timestep.size());
226 
227  for (std::size_t i = 0; i < xdot.size(); i++)
228  {
229  xdot[i] =
231  _time_disc.getXdot(*x_new_timestep[i], *x_prev[i], *xdot[i]);
232  }
233 
234  _M->setZero();
235  _K->setZero();
236  _b->setZero();
237 
238  _ode.preAssemble(t, dt, x_curr);
239  _ode.assemble(t, dt, x_new_timestep, xdot, process_id, *_M, *_K, *_b);
240 
244 
245  for (auto& v : xdot)
246  {
248  }
249 }
void setZero()
reset data entries to zero.
Definition: EigenMatrix.h:65
void getXdot(GlobalVector const &x_at_new_timestep, GlobalVector const &x_old, GlobalVector &xdot) const
virtual double getCurrentTimeIncrement() const =0
virtual double getCurrentTime() const =0
void finalizeAssembly(PETScMatrix &A)
Definition: LinAlg.cpp:162

References MathLib::LinAlg::finalizeAssembly(), NumLib::VectorProvider::getVector(), NumLib::GlobalVectorProvider::provider, NumLib::VectorProvider::releaseVector(), and MathLib::EigenVector::setZero().

◆ computeKnownSolutions()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::computeKnownSolutions ( GlobalVector const &  x,
int const  process_id 
)
override

Definition at line 253 of file TimeDiscretizedODESystem.cpp.

255 {
257  _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
258 }

◆ getA()

Definition at line 188 of file TimeDiscretizedODESystem.h.

189  {
190  _mat_trans->computeA(*_M, *_K, A);
191  }

◆ getMatrixSpecifications()

Definition at line 219 of file TimeDiscretizedODESystem.h.

221  {
222  return _ode.getMatrixSpecifications(process_id);
223  }

◆ getRhs()

Definition at line 193 of file TimeDiscretizedODESystem.h.

194  {
195  _mat_trans->computeRhs(*_M, *_K, *_b, x_prev, rhs);
196  }

◆ getTimeDiscretization()

Exposes the used time discretization scheme.

Implements NumLib::TimeDiscretizedODESystemBase< NonlinearSolverTag::Picard >.

Definition at line 218 of file TimeDiscretizedODESystem.h.

218 { return _time_disc; }

◆ isLinear()

Definition at line 206 of file TimeDiscretizedODESystem.h.

206 { return _ode.isLinear(); }

◆ postIteration()

Definition at line 213 of file TimeDiscretizedODESystem.h.

214  {
215  return _ode.postIteration(x);
216  }

◆ preIteration()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::preIteration ( const unsigned  iter,
GlobalVector const &  x 
)
inlineoverride

Definition at line 208 of file TimeDiscretizedODESystem.h.

209  {
210  _ode.preIteration(iter, x);
211  }

Member Data Documentation

◆ _b

◆ _b_id

ID of the _b vector.

Definition at line 242 of file TimeDiscretizedODESystem.h.

◆ _K

◆ _K_id

ID of the _K matrix.

Definition at line 241 of file TimeDiscretizedODESystem.h.

◆ _known_solutions

Initial value:
=
nullptr

stores precomputed values for known solutions

Definition at line 233 of file TimeDiscretizedODESystem.h.

◆ _M

◆ _M_id

ID of the _M matrix.

Definition at line 240 of file TimeDiscretizedODESystem.h.

◆ _mat_trans

the object used to compute the matrix/vector for the nonlinear solver

Definition at line 230 of file TimeDiscretizedODESystem.h.

◆ _ode

◆ _time_disc

the time discretization to being used

Definition at line 227 of file TimeDiscretizedODESystem.h.

◆ _xdot_ids

Definition at line 243 of file TimeDiscretizedODESystem.h.

◆ ODETag


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