OGS
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 >:
[legend]
Collaboration diagram for NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >:
[legend]

Public Types

using ODE = ODESystem<ODETag, NonlinearSolverTag::Newton>
 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)
 
 ~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 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 const &x, GlobalVector &minus_delta_x) const override
 
void applyKnownSolutionsPETScSNES (GlobalMatrix &Jac, GlobalVector &res, GlobalVector &x) const override
 
void updateConstraints (GlobalVector &lower, GlobalVector &upper, int const process_id) 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.
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) 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_Jac
 the Jacobian of the residual
 
GlobalVector_b
 Matrix \( b \).
 
std::size_t _Jac_id = 0u
 ID of the _Jac 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 70 of file TimeDiscretizedODESystem.h.

◆ ODE

◆ TimeDisc

A shortcut for a general time discretization scheme.

Definition at line 72 of file TimeDiscretizedODESystem.h.

Constructor & Destructor Documentation

◆ TimeDiscretizedODESystem()

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

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

52 : _ode(ode),
53 _time_disc(time_discretization),
54 _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
55{
57 _ode.getMatrixSpecifications(process_id), _Jac_id);
59 _ode.getMatrixSpecifications(process_id), _b_id);
60}
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 119 of file TimeDiscretizedODESystem.cpp.

120{
122}
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 const & x,
GlobalVector & minus_delta_x ) const
override

Definition at line 125 of file TimeDiscretizedODESystem.cpp.

129{
130 if (!_known_solutions)
131 {
132 return;
133 }
134
136 std::size_t const size = ranges::accumulate(
137 *_known_solutions | ranges::views::transform([](auto const& bc)
138 { return bc.ids.size(); }),
139 0);
140 std::vector<IndexType> ids;
141 ids.reserve(size);
142 std::vector<double> values;
143 values.reserve(size);
144
145 for (auto const& bc : *_known_solutions)
146 {
147 for (std::size_t i = 0; i < bc.ids.size(); ++i)
148 {
149 auto const id = bc.ids[i];
150 ids.push_back(id);
151 // minus_delta_x will be set to the difference between the current
152 // value and the Dirichlet BC value.
153 values.push_back(x[id] - bc.values[i]);
154 }
155 }
156
157 MathLib::applyKnownSolution(Jac, res, minus_delta_x, ids, values);
158}
constexpr int size(int const displacement_dim)
Vectorized tensor size for given displacement dimension.
void applyKnownSolution(EigenMatrix &A, EigenVector &b, EigenVector &, const std::vector< EigenMatrix::IndexType > &vec_knownX_id, const std::vector< double > &vec_knownX_x)
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:225

References MathLib::applyKnownSolution().

◆ applyKnownSolutionsPETScSNES()

Definition at line 161 of file TimeDiscretizedODESystem.cpp.

164{
165 if (!_known_solutions)
166 {
167 return;
168 }
169
171 std::vector<IndexType> ids;
172 for (auto const& bc : *_known_solutions)
173 {
174 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
175 }
176
177 // For the Newton method the values must be zero
178 std::vector<double> values(ids.size(), 0);
179 MathLib::applyKnownSolution(Jac, res, x, ids, values);
180}

References MathLib::applyKnownSolution().

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

75{
76 auto const t = _time_disc.getCurrentTime();
77 auto const dt = _time_disc.getCurrentTimeIncrement();
78 auto const& x_curr = *x_new_timestep[process_id];
79
80 _b->setZero();
81 _Jac->setZero();
82
83 _ode.preAssemble(t, dt, x_curr);
84 _ode.assembleWithJacobian(t, dt, x_new_timestep, x_prev, process_id, *_b,
85 *_Jac);
86
89}
void setZero()
reset data entries to zero.
Definition EigenMatrix.h:65
virtual double getCurrentTimeIncrement() const =0
virtual double getCurrentTime() const =0
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:170

References MathLib::LinAlg::finalizeAssembly().

◆ computeKnownSolutions()

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

Definition at line 110 of file TimeDiscretizedODESystem.cpp.

112{
114 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
115}

◆ getIndicesOfResiduumWithoutInitialCompensation()

std::vector< GlobalIndexType > NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::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 92 of file TimeDiscretizedODESystem.h.

93 {
94 return _ode.getIndicesOfResiduumWithoutInitialCompensation();
95 }

◆ getJacobian()

Definition at line 103 of file TimeDiscretizedODESystem.cpp.

104{
105 _mat_trans->computeJacobian(*_Jac, Jac);
106}

◆ getMatrixSpecifications()

Definition at line 134 of file TimeDiscretizedODESystem.h.

136 {
137 return _ode.getMatrixSpecifications(process_id);
138 }

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

96{
97 MathLib::LinAlg::copy(*_b, res); // res = b
98 MathLib::LinAlg::scale(res, -1.); // res = -b
99}
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:44

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

◆ getTimeDiscretization()

Exposes the used time discretization scheme.

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

Definition at line 133 of file TimeDiscretizedODESystem.h.

133{ return _time_disc; }

◆ isLinear()

Definition at line 121 of file TimeDiscretizedODESystem.h.

121{ return _ode.isLinear(); }

◆ postIteration()

Definition at line 128 of file TimeDiscretizedODESystem.h.

129 {
130 return _ode.postIteration(x);
131 }

◆ preIteration()

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

Definition at line 123 of file TimeDiscretizedODESystem.h.

124 {
125 _ode.preIteration(iter, x);
126 }

◆ updateConstraints()

void NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::updateConstraints ( GlobalVector & lower,
GlobalVector & upper,
int const process_id )
inlineoverride

Definition at line 115 of file TimeDiscretizedODESystem.h.

117 {
118 _ode.updateConstraints(lower, upper, process_id);
119 }

Member Data Documentation

◆ _b

◆ _b_id

ID of the _b vector.

Definition at line 155 of file TimeDiscretizedODESystem.h.

◆ _Jac

◆ _Jac_id

ID of the _Jac matrix.

Definition at line 154 of file TimeDiscretizedODESystem.h.

◆ _known_solutions

Initial value:
=
nullptr

stores precomputed values for known solutions

Definition at line 148 of file TimeDiscretizedODESystem.h.

◆ _mat_trans

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

Definition at line 145 of file TimeDiscretizedODESystem.h.

◆ _ode

◆ _time_disc

the time discretization to being used

Definition at line 142 of file TimeDiscretizedODESystem.h.

◆ ODETag


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