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 176 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.
using MatTrans = MatrixTranslator<ODETag>
using TimeDisc = TimeDiscretization
 A shortcut for a general time discretization scheme.

Public Member Functions

 TimeDiscretizedODESystem (const int process_id, ODE &ode, TimeDisc &time_discretization)
 Constructs a new instance.
 ~TimeDiscretizedODESystem () override
void assemble (std::vector< GlobalVector * > const &x_new_timestep, std::vector< GlobalVector * > const &x_prev, int const process_id) override
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
void setReleaseNodalForces (GlobalVector const *r_neq, int const process_id) override
void getA (GlobalMatrix &A) const override
void getRhs (GlobalVector const &x_prev, GlobalVector &rhs) const override
void getAandRhsNormalized (GlobalMatrix &A, 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, MathLib::DirichletBCApplicationMode const mode) const override
bool isLinear () const override
bool requiresNormalization () 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.
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
MathLib::LinearSolverBehaviour linearSolverNeedsToCompute () const override

Static Public Attributes

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

Private Types

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

Private Attributes

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

Member Typedef Documentation

◆ Index

◆ MatTrans

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

Definition at line 189 of file TimeDiscretizedODESystem.h.

◆ ODE

◆ TimeDisc

A shortcut for a general time discretization scheme.

Definition at line 191 of file TimeDiscretizedODESystem.h.

Constructor & Destructor Documentation

◆ TimeDiscretizedODESystem()

NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::TimeDiscretizedODESystem ( const int process_id,
ODE & ode,
TimeDisc & time_discretization )
explicit

Constructs a new instance.

Definition at line 194 of file TimeDiscretizedODESystem.cpp.

197 : _ode(ode),
200{
202 ode.getMatrixSpecifications(process_id), _M_id);
204 ode.getMatrixSpecifications(process_id), _K_id);
206 ode.getMatrixSpecifications(process_id), _b_id);
207}
std::unique_ptr< MatTrans > _mat_trans
the object used to compute the matrix/vector for the nonlinear solver

References _b, _b_id, _K, _K_id, _M, _M_id, _mat_trans, _ode, _time_disc, NumLib::createMatrixTranslator(), ODETag, NumLib::Picard, NumLib::GlobalMatrixProvider::provider, and NumLib::GlobalVectorProvider::provider.

◆ ~TimeDiscretizedODESystem()

Member Function Documentation

◆ applyKnownSolutions()

Definition at line 253 of file TimeDiscretizedODESystem.cpp.

254{
256}
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 _known_solutions, detail::applyKnownSolutions(), applyKnownSolutions(), and NumLib::Picard.

Referenced by applyKnownSolutions().

◆ applyKnownSolutionsPicard()

Definition at line 259 of file TimeDiscretizedODESystem.cpp.

263{
264 if (!_known_solutions)
265 {
266 return;
267 }
268
272 for (auto const& bc : *_known_solutions)
273 {
274 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
275 values.insert(end(values), begin(bc.values), end(bc.values));
276 }
278}
void applyKnownSolution(EigenMatrix &A, EigenVector &b, EigenVector &, const std::vector< EigenMatrix::IndexType > &vec_knownX_id, const std::vector< double > &vec_knownX_x, DirichletBCApplicationMode const mode)

References _known_solutions, MathLib::applyKnownSolution(), and NumLib::Picard.

◆ 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 219 of file TimeDiscretizedODESystem.cpp.

223{
224 namespace LinAlg = MathLib::LinAlg;
225
226 auto const t = _time_disc.getCurrentTime();
227 auto const dt = _time_disc.getCurrentTimeIncrement();
228 auto const& x_curr = *x_new_timestep[process_id];
229
230 _M->setZero();
231 _K->setZero();
232 _b->setZero();
233
234 _ode.preAssemble(t, dt, x_curr);
235 _ode.assemble(t, dt, x_new_timestep, x_prev, process_id, *_M, *_K, *_b);
236
240}
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:198

References _b, _K, _M, _ode, _time_disc, MathLib::LinAlg::finalizeAssembly(), and NumLib::Picard.

◆ computeKnownSolutions()

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

Definition at line 244 of file TimeDiscretizedODESystem.cpp.

246{
248 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
249}

References _known_solutions, _ode, _time_disc, computeKnownSolutions(), and NumLib::Picard.

Referenced by computeKnownSolutions().

◆ getA()

Definition at line 217 of file TimeDiscretizedODESystem.h.

218 {
219 _mat_trans->computeA(*_M, *_K, A);
220 }

References _K, _M, and _mat_trans.

◆ getAandRhsNormalized()

Definition at line 227 of file TimeDiscretizedODESystem.h.

228 {
229 _mat_trans->normalizeAandRhs(A, rhs);
230 }

References _mat_trans.

◆ getIndicesOfResiduumWithoutInitialCompensation()

std::vector< GlobalIndexType > NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::getIndicesOfResiduumWithoutInitialCompensation ( ) const
inlineoverride
Returns
The global indices for the entries of the global residuum vector that do not need initial non-equilibrium compensation.

Definition at line 206 of file TimeDiscretizedODESystem.h.

207 {
208 return _ode.getIndicesOfResiduumWithoutInitialCompensation();
209 }

References _ode.

◆ getMatrixSpecifications()

Definition at line 259 of file TimeDiscretizedODESystem.h.

261 {
262 return _ode.getMatrixSpecifications(process_id);
263 }

References _ode.

◆ getRhs()

Definition at line 222 of file TimeDiscretizedODESystem.h.

223 {
224 _mat_trans->computeRhs(*_M, *_K, *_b, x_prev, rhs);
225 }

References _b, _K, _M, and _mat_trans.

◆ getTimeDiscretization()

Exposes the used time discretization scheme.

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

Definition at line 258 of file TimeDiscretizedODESystem.h.

258{ return _time_disc; }

References _time_disc.

◆ isLinear()

Definition at line 241 of file TimeDiscretizedODESystem.h.

241{ return _ode.isLinear(); }

References _ode.

◆ linearSolverNeedsToCompute()

Definition at line 265 of file TimeDiscretizedODESystem.h.

266 {
267 if (_ode.shouldLinearSolverComputeOnlyUponTimestepChange() &&
268 _time_disc.getCurrentTimeIncrement() !=
269 _time_disc.getPreviousTimeIncrement())
270 {
272 }
273 if (_ode.shouldLinearSolverComputeOnlyUponTimestepChange() &&
274 _time_disc.getCurrentTimeIncrement() ==
275 _time_disc.getPreviousTimeIncrement())
276 {
278 }
280 }

References _ode, _time_disc, MathLib::RECOMPUTE, MathLib::RECOMPUTE_AND_STORE, and MathLib::REUSE.

◆ postIteration()

Definition at line 253 of file TimeDiscretizedODESystem.h.

254 {
255 return _ode.postIteration(x);
256 }

References _ode.

◆ preIteration()

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

Definition at line 248 of file TimeDiscretizedODESystem.h.

249 {
250 _ode.preIteration(iter, x);
251 }

References _ode.

◆ requiresNormalization()

Definition at line 243 of file TimeDiscretizedODESystem.h.

244 {
245 return _ode.requiresNormalization();
246 }

References _ode.

◆ setReleaseNodalForces()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::setReleaseNodalForces ( GlobalVector const * r_neq,
int const process_id )
inlineoverride

Definition at line 211 of file TimeDiscretizedODESystem.h.

213 {
214 _ode.setReleaseNodalForces(r_neq, process_id);
215 }

References _ode.

Member Data Documentation

◆ _b

◆ _b_id

ID of the _b vector.

Definition at line 299 of file TimeDiscretizedODESystem.h.

Referenced by TimeDiscretizedODESystem().

◆ _K

◆ _K_id

ID of the _K matrix.

Definition at line 298 of file TimeDiscretizedODESystem.h.

Referenced by TimeDiscretizedODESystem().

◆ _known_solutions

Initial value:
=
nullptr

stores precomputed values for known solutions

Definition at line 290 of file TimeDiscretizedODESystem.h.

Referenced by applyKnownSolutions(), applyKnownSolutionsPicard(), and computeKnownSolutions().

◆ _M

◆ _M_id

ID of the _M matrix.

Definition at line 297 of file TimeDiscretizedODESystem.h.

Referenced by TimeDiscretizedODESystem().

◆ _mat_trans

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

Definition at line 287 of file TimeDiscretizedODESystem.h.

Referenced by TimeDiscretizedODESystem(), getA(), getAandRhsNormalized(), and getRhs().

◆ _ode

◆ _time_disc

◆ ODETag


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