OGS
TimeDiscretizedODESystem.cpp
Go to the documentation of this file.
1
12
13#include <range/v3/numeric/accumulate.hpp>
14#include <range/v3/view/transform.hpp>
15
18#include "NumLib/Exceptions.h"
20
21namespace detail
22{
24template <typename Solutions, typename Vector>
25void applyKnownSolutions(std::vector<Solutions> const* const known_solutions,
26 Vector& x)
27{
28 if (!known_solutions)
29 {
30 return;
31 }
32
33 for (auto const& bc : *known_solutions)
34 {
35 for (std::size_t i = 0; i < bc.ids.size(); ++i)
36 {
37 // TODO that might have bad performance for some Vector types, e.g.,
38 // PETSc.
39 MathLib::setVector(x, bc.ids[i], bc.values[i]);
40 }
41 }
43}
44} // namespace detail
45
46namespace NumLib
47{
48TimeDiscretizedODESystem<ODESystemTag::FirstOrderImplicitQuasilinear,
50 TimeDiscretizedODESystem(const int process_id, ODE& ode,
51 TimeDisc& time_discretization)
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}
61
69
72 assemble(std::vector<GlobalVector*> const& x_new_timestep,
73 std::vector<GlobalVector*> const& x_prev,
74 int const process_id)
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}
90
93 getResidual(GlobalVector const& /*x_new_timestep*/,
94 GlobalVector const& /*x_prev*/,
95 GlobalVector& res) const
96{
97 MathLib::LinAlg::copy(*_b, res); // res = b
98 MathLib::LinAlg::scale(res, -1.); // res = -b
99}
100
103 NonlinearSolverTag::Newton>::getJacobian(GlobalMatrix& Jac) const
104{
105 _mat_trans->computeJacobian(*_Jac, Jac);
106}
107
110 NonlinearSolverTag::Newton>::computeKnownSolutions(GlobalVector const& x,
111 int const process_id)
112{
113 _known_solutions =
114 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
115}
116
119 NonlinearSolverTag::Newton>::applyKnownSolutions(GlobalVector& x) const
120{
121 ::detail::applyKnownSolutions(_known_solutions, x);
122}
123
126 applyKnownSolutionsNewton(GlobalMatrix& Jac, GlobalVector& res,
127 GlobalVector const& x,
128 GlobalVector& minus_delta_x) const
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}
159
162 applyKnownSolutionsPETScSNES(GlobalMatrix& Jac, GlobalVector& res,
163 GlobalVector& x) const
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}
181
184 TimeDiscretizedODESystem(const int process_id, ODE& ode,
185 TimeDisc& time_discretization)
186 : _ode(ode),
187 _time_disc(time_discretization),
188 _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
189{
191 ode.getMatrixSpecifications(process_id), _M_id);
193 ode.getMatrixSpecifications(process_id), _K_id);
195 ode.getMatrixSpecifications(process_id), _b_id);
196}
197
206
209 assemble(std::vector<GlobalVector*> const& x_new_timestep,
210 std::vector<GlobalVector*> const& x_prev,
211 int const process_id)
212{
213 namespace LinAlg = MathLib::LinAlg;
214
215 auto const t = _time_disc.getCurrentTime();
216 auto const dt = _time_disc.getCurrentTimeIncrement();
217 auto const& x_curr = *x_new_timestep[process_id];
218
219 _M->setZero();
220 _K->setZero();
221 _b->setZero();
222
223 _ode.preAssemble(t, dt, x_curr);
224 _ode.assemble(t, dt, x_new_timestep, x_prev, process_id, *_M, *_K, *_b);
225
229}
230
233 NonlinearSolverTag::Picard>::computeKnownSolutions(GlobalVector const& x,
234 int const process_id)
235{
236 _known_solutions =
237 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
238}
239
242 NonlinearSolverTag::Picard>::applyKnownSolutions(GlobalVector& x) const
243{
244 ::detail::applyKnownSolutions(_known_solutions, x);
245}
246
249 applyKnownSolutionsPicard(GlobalMatrix& A,
250 GlobalVector& rhs,
251 GlobalVector& x) const
252{
253 if (!_known_solutions)
254 {
255 return;
256 }
257
259 std::vector<IndexType> ids;
260 std::vector<double> values;
261 for (auto const& bc : *_known_solutions)
262 {
263 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
264 values.insert(end(values), begin(bc.values), end(bc.values));
265 }
266 MathLib::applyKnownSolution(A, rhs, x, ids, values);
267}
268
269} // 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 copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:44
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