OGS
TimeDiscretizedODESystem.cpp
Go to the documentation of this file.
1
12
15#include "NumLib/Exceptions.h"
17
18namespace detail
19{
21template <typename Solutions, typename Vector>
22void applyKnownSolutions(std::vector<Solutions> const* const known_solutions,
23 Vector& x)
24{
25 if (!known_solutions)
26 {
27 return;
28 }
29
30 for (auto const& bc : *known_solutions)
31 {
32 for (std::size_t i = 0; i < bc.ids.size(); ++i)
33 {
34 // TODO that might have bad performance for some Vector types, e.g.,
35 // PETSc.
36 MathLib::setVector(x, bc.ids[i], bc.values[i]);
37 }
38 }
40}
41} // namespace detail
42
43namespace NumLib
44{
45TimeDiscretizedODESystem<ODESystemTag::FirstOrderImplicitQuasilinear,
47 TimeDiscretizedODESystem(const int process_id, ODE& ode,
48 TimeDisc& time_discretization)
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}
62
72
75 assemble(std::vector<GlobalVector*> const& x_new_timestep,
76 std::vector<GlobalVector*> const& x_prev,
77 int const process_id)
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
85 _M->setZero();
86 _K->setZero();
87 _b->setZero();
88 _Jac->setZero();
89
90 _ode.preAssemble(t, dt, x_curr);
91 _ode.assembleWithJacobian(t, dt, x_new_timestep, x_prev, process_id, *_M,
92 *_K, *_b, *_Jac);
93
98}
99
102 NonlinearSolverTag::Newton>::getResidual(GlobalVector const& x_new_timestep,
103 GlobalVector const& x_prev,
104 GlobalVector& res) const
105{
106 double const dt = _time_disc.getCurrentTimeIncrement();
107 _mat_trans->computeResidual(*_M, *_K, *_b, dt, x_new_timestep, x_prev, res);
108}
109
112 NonlinearSolverTag::Newton>::getJacobian(GlobalMatrix& Jac) const
113{
114 _mat_trans->computeJacobian(*_Jac, Jac);
115}
116
119 NonlinearSolverTag::Newton>::computeKnownSolutions(GlobalVector const& x,
120 int const process_id)
121{
122 _known_solutions =
123 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
124}
125
128 NonlinearSolverTag::Newton>::applyKnownSolutions(GlobalVector& x) const
129{
130 ::detail::applyKnownSolutions(_known_solutions, x);
131}
132
135 applyKnownSolutionsNewton(GlobalMatrix& Jac, GlobalVector& res,
136 GlobalVector& minus_delta_x) const
137{
138 if (!_known_solutions)
139 {
140 return;
141 }
142
144 std::vector<IndexType> ids;
145 for (auto const& bc : *_known_solutions)
146 {
147 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
148 }
149
150 // For the Newton method the values must be zero
151 std::vector<double> values(ids.size(), 0);
152 MathLib::applyKnownSolution(Jac, res, minus_delta_x, ids, values);
153}
154
157 TimeDiscretizedODESystem(const int process_id, ODE& ode,
158 TimeDisc& time_discretization)
159 : _ode(ode),
160 _time_disc(time_discretization),
161 _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
162{
164 ode.getMatrixSpecifications(process_id), _M_id);
166 ode.getMatrixSpecifications(process_id), _K_id);
168 ode.getMatrixSpecifications(process_id), _b_id);
169}
170
179
182 assemble(std::vector<GlobalVector*> const& x_new_timestep,
183 std::vector<GlobalVector*> const& x_prev,
184 int const process_id)
185{
186 namespace LinAlg = MathLib::LinAlg;
187
188 auto const t = _time_disc.getCurrentTime();
189 auto const dt = _time_disc.getCurrentTimeIncrement();
190 auto const& x_curr = *x_new_timestep[process_id];
191
192 _M->setZero();
193 _K->setZero();
194 _b->setZero();
195
196 _ode.preAssemble(t, dt, x_curr);
197 _ode.assemble(t, dt, x_new_timestep, x_prev, process_id, *_M, *_K, *_b);
198
202}
203
206 NonlinearSolverTag::Picard>::computeKnownSolutions(GlobalVector const& x,
207 int const process_id)
208{
209 _known_solutions =
210 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
211}
212
215 NonlinearSolverTag::Picard>::applyKnownSolutions(GlobalVector& x) const
216{
217 ::detail::applyKnownSolutions(_known_solutions, x);
218}
219
222 applyKnownSolutionsPicard(GlobalMatrix& A,
223 GlobalVector& rhs,
224 GlobalVector& x) const
225{
226 if (!_known_solutions)
227 {
228 return;
229 }
230
232 std::vector<IndexType> ids;
233 std::vector<double> values;
234 for (auto const& bc : *_known_solutions)
235 {
236 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
237 values.insert(end(values), begin(bc.values), end(bc.values));
238 }
239 MathLib::applyKnownSolution(A, rhs, x, ids, values);
240}
241
242} // namespace NumLib
Global vector based on Eigen vector.
Definition EigenVector.h:25
virtual void releaseMatrix(GlobalMatrix const &A)=0
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
std::unique_ptr< MatrixTranslator< ODETag > > createMatrixTranslator(TimeDiscretization const &timeDisc)
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:170
void applyKnownSolution(EigenMatrix &A, EigenVector &b, EigenVector &, const std::vector< EigenMatrix::IndexType > &vec_knownX_id, const std::vector< double > &vec_knownX_x)
void setVector(PETScVector &v, std::initializer_list< double > values)
void applyKnownSolutions(std::vector< Solutions > const *const known_solutions, Vector &x)
Applies known solutions to the solution vector x.
static NUMLIB_EXPORT MatrixProvider & provider
static NUMLIB_EXPORT VectorProvider & provider