OGS  v6.3.3
NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton > Class Referencefinal

Detailed Description

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

See also
ODESystemTag::FirstOrderImplicitQuasilinear

Definition at line 57 of file TimeDiscretizedODESystem.h.

#include <TimeDiscretizedODESystem.h>

Inheritance diagram for NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >:
Collaboration diagram for NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >:

Public Types

using ODE = ODESystem< ODETag, NonlinearSolverTag::Newton >
 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)
 
 ~TimeDiscretizedODESystem () override
 
void assemble (std::vector< GlobalVector * > const &x_new_timestep, std::vector< GlobalVector * > const &x_prev, int const process_id) override
 
void getResidual (GlobalVector const &x_new_timestep, GlobalVector const &x_prev, GlobalVector &res) const override
 
void getJacobian (GlobalMatrix &Jac) const override
 
void computeKnownSolutions (GlobalVector const &x, int const process_id) override
 
void applyKnownSolutions (GlobalVector &x) const override
 
void applyKnownSolutionsNewton (GlobalMatrix &Jac, GlobalVector &res, GlobalVector &minus_delta_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 * _Jac
 the Jacobian of the residual More...
 
GlobalMatrix * _M
 Matrix \( M \). More...
 
GlobalMatrix * _K
 Matrix \( K \). More...
 
GlobalVector * _b
 Matrix \( b \). More...
 
std::size_t _Jac_id = 0u
 ID of the _Jac matrix. 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::size_t _xdot_id = 0u
 ID of the vector storing xdot in intermediate computations. More...
 

Member Typedef Documentation

◆ Index

◆ MatTrans

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

Definition at line 70 of file TimeDiscretizedODESystem.h.

◆ ODE

◆ TimeDisc

Constructor & Destructor Documentation

◆ TimeDiscretizedODESystem()

Constructs a new instance.

Parameters
process_idID of the ODE to be solved.
odethe ODE to be wrapped.
time_discretizationthe time discretization to be used.

Definition at line 46 of file TimeDiscretizedODESystem.cpp.

49  : _ode(ode),
50  _time_disc(time_discretization),
51  _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
52 {
54  _ode.getMatrixSpecifications(process_id), _Jac_id);
56  _ode.getMatrixSpecifications(process_id), _M_id);
58  _ode.getMatrixSpecifications(process_id), _K_id);
60  _ode.getMatrixSpecifications(process_id), _b_id);
61 }
virtual GlobalMatrix & getMatrix()=0
Get an uninitialized matrix.
std::unique_ptr< MatTrans > _mat_trans
the object used to compute the matrix/vector for the nonlinear solver
virtual GlobalVector & getVector()=0
Get an uninitialized vector.
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()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::applyKnownSolutions ( GlobalVector &  x) const
override

Definition at line 159 of file TimeDiscretizedODESystem.cpp.

160 {
162 }
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().

◆ applyKnownSolutionsNewton()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::applyKnownSolutionsNewton ( GlobalMatrix &  Jac,
GlobalVector &  res,
GlobalVector &  minus_delta_x 
) const
override

Definition at line 165 of file TimeDiscretizedODESystem.cpp.

168 {
169  if (!_known_solutions)
170  {
171  return;
172  }
173 
175  std::vector<IndexType> ids;
176  for (auto const& bc : *_known_solutions)
177  {
178  std::copy(bc.ids.cbegin(), bc.ids.cend(), std::back_inserter(ids));
179  }
180 
181  // For the Newton method the values must be zero
182  std::vector<double> values(ids.size(), 0);
183  MathLib::applyKnownSolution(Jac, res, minus_delta_x, ids, values);
184 }
void copy(MatrixOrVector const &x, MatrixOrVector &y)
Copies x to y.
Definition: LinAlg.h:37
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:19

References MathLib::applyKnownSolution(), and MathLib::LinAlg::copy().

◆ assemble()

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

Definition at line 74 of file TimeDiscretizedODESystem.cpp.

78 {
79  namespace LinAlg = MathLib::LinAlg;
80 
81  auto const t = _time_disc.getCurrentTime();
82  auto const dt = _time_disc.getCurrentTimeIncrement();
83  auto const& x_curr = *x_new_timestep[process_id];
84  auto const dxdot_dx = _time_disc.getNewXWeight();
85 
86  std::vector<GlobalVector*> xdot(x_new_timestep.size());
87  for (std::size_t i = 0; i < xdot.size(); i++)
88  {
90  _time_disc.getXdot(*x_new_timestep[i], *x_prev[i], *xdot[i]);
91  }
92 
93  _M->setZero();
94  _K->setZero();
95  _b->setZero();
96  _Jac->setZero();
97 
98  _ode.preAssemble(t, dt, x_curr);
99  try
100  {
101  _ode.assembleWithJacobian(t, dt, x_new_timestep, xdot, dxdot_dx, 1.0,
102  process_id, *_M, *_K, *_b, *_Jac);
103  }
104  catch (AssemblyException const&)
105  {
106  for (auto& v : xdot)
107  {
109  }
110  throw;
111  }
112 
117 
118  for (auto& v : xdot)
119  {
121  }
122 }
void getXdot(GlobalVector const &x_at_new_timestep, GlobalVector const &x_old, GlobalVector &xdot) const
virtual double getCurrentTimeIncrement() const =0
virtual double getNewXWeight() const =0
Returns .
virtual double getCurrentTime() const =0
void finalizeAssembly(Matrix &)

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

◆ computeKnownSolutions()

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

Definition at line 150 of file TimeDiscretizedODESystem.cpp.

152 {
154  _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
155 }

◆ getJacobian()

Definition at line 143 of file TimeDiscretizedODESystem.cpp.

144 {
145  _mat_trans->computeJacobian(*_Jac, Jac);
146 }

◆ getMatrixSpecifications()

Definition at line 116 of file TimeDiscretizedODESystem.h.

118  {
119  return _ode.getMatrixSpecifications(process_id);
120  }

◆ getResidual()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::getResidual ( GlobalVector const &  x_new_timestep,
GlobalVector const &  x_prev,
GlobalVector &  res 
) const
override

Definition at line 126 of file TimeDiscretizedODESystem.cpp.

129 {
130  // TODO Maybe the duplicate calculation of xdot here and in assembleJacobian
131  // can be optimuized. However, that would make the interface a bit more
132  // fragile.
134  _time_disc.getXdot(x_new_timestep, x_prev, xdot);
135 
136  _mat_trans->computeResidual(*_M, *_K, *_b, x_new_timestep, xdot, res);
137 
139 }

References NumLib::VectorProvider::getVector(), NumLib::GlobalVectorProvider::provider, and NumLib::VectorProvider::releaseVector().

◆ getTimeDiscretization()

Exposes the used time discretization scheme.

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

Definition at line 115 of file TimeDiscretizedODESystem.h.

115 { return _time_disc; }

◆ isLinear()

Definition at line 103 of file TimeDiscretizedODESystem.h.

103 { return _ode.isLinear(); }

◆ postIteration()

Definition at line 110 of file TimeDiscretizedODESystem.h.

111  {
112  return _ode.postIteration(x);
113  }

◆ preIteration()

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

Definition at line 105 of file TimeDiscretizedODESystem.h.

106  {
107  _ode.preIteration(iter, x);
108  }

Member Data Documentation

◆ _b

◆ _b_id

ID of the _b vector.

Definition at line 141 of file TimeDiscretizedODESystem.h.

◆ _Jac

the Jacobian of the residual

Definition at line 133 of file TimeDiscretizedODESystem.h.

◆ _Jac_id

ID of the _Jac matrix.

Definition at line 138 of file TimeDiscretizedODESystem.h.

◆ _K

◆ _K_id

ID of the _K matrix.

Definition at line 140 of file TimeDiscretizedODESystem.h.

◆ _known_solutions

Initial value:
=
nullptr

stores precomputed values for known solutions

Definition at line 130 of file TimeDiscretizedODESystem.h.

◆ _M

◆ _M_id

ID of the _M matrix.

Definition at line 139 of file TimeDiscretizedODESystem.h.

◆ _mat_trans

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

Definition at line 127 of file TimeDiscretizedODESystem.h.

◆ _ode

◆ _time_disc

the time discretization to being used

Definition at line 124 of file TimeDiscretizedODESystem.h.

◆ _xdot_id

ID of the vector storing xdot in intermediate computations.

Definition at line 144 of file TimeDiscretizedODESystem.h.

◆ ODETag


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