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 try
83 {
84 _ode.preAssemble(t, dt, x_curr);
85 _ode.assembleWithJacobian(t, dt, x_new_timestep, x_prev, process_id,
86 *_b, *_Jac);
87 }
88 catch (AssemblyException const&)
89 {
92 throw;
93 }
96}
97
100 getResidual(GlobalVector const& /*x_new_timestep*/,
101 GlobalVector const& /*x_prev*/,
102 GlobalVector& res) const
103{
104 MathLib::LinAlg::copy(*_b, res); // res = b
105 MathLib::LinAlg::scale(res, -1.); // res = -b
106}
107
110 NonlinearSolverTag::Newton>::getJacobian(GlobalMatrix& Jac) const
111{
112 _mat_trans->computeJacobian(*_Jac, Jac);
113}
114
117 NonlinearSolverTag::Newton>::computeKnownSolutions(GlobalVector const& x,
118 int const process_id)
119{
120 _known_solutions =
121 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
122}
123
126 NonlinearSolverTag::Newton>::applyKnownSolutions(GlobalVector& x) const
127{
128 ::detail::applyKnownSolutions(_known_solutions, x);
129}
130
133 applyKnownSolutionsNewton(GlobalMatrix& Jac, GlobalVector& res,
134 GlobalVector const& x,
135 GlobalVector& minus_delta_x) const
136{
137 if (!_known_solutions)
138 {
139 return;
140 }
141
143 std::size_t const size = ranges::accumulate(
144 *_known_solutions | ranges::views::transform([](auto const& bc)
145 { return bc.ids.size(); }),
146 0);
147 std::vector<IndexType> ids;
148 ids.reserve(size);
149 std::vector<double> values;
150 values.reserve(size);
151
152 for (auto const& bc : *_known_solutions)
153 {
154 for (std::size_t i = 0; i < bc.ids.size(); ++i)
155 {
156 auto const id = bc.ids[i];
157 ids.push_back(id);
158 // minus_delta_x will be set to the difference between the current
159 // value and the Dirichlet BC value.
160 values.push_back(x[id] - bc.values[i]);
161 }
162 }
163
164 MathLib::applyKnownSolution(Jac, res, minus_delta_x, ids, values);
165}
166
169 applyKnownSolutionsPETScSNES(GlobalMatrix& Jac, GlobalVector& res,
170 GlobalVector& x) const
171{
172 if (!_known_solutions)
173 {
174 return;
175 }
176
178 std::vector<IndexType> ids;
179 for (auto const& bc : *_known_solutions)
180 {
181 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
182 }
183
184 // For the Newton method the values must be zero
185 std::vector<double> values(ids.size(), 0);
186 MathLib::applyKnownSolution(Jac, res, x, ids, values);
187}
188
191 TimeDiscretizedODESystem(const int process_id, ODE& ode,
192 TimeDisc& time_discretization)
193 : _ode(ode),
194 _time_disc(time_discretization),
195 _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
196{
198 ode.getMatrixSpecifications(process_id), _M_id);
200 ode.getMatrixSpecifications(process_id), _K_id);
202 ode.getMatrixSpecifications(process_id), _b_id);
203}
204
213
216 assemble(std::vector<GlobalVector*> const& x_new_timestep,
217 std::vector<GlobalVector*> const& x_prev,
218 int const process_id)
219{
220 namespace LinAlg = MathLib::LinAlg;
221
222 auto const t = _time_disc.getCurrentTime();
223 auto const dt = _time_disc.getCurrentTimeIncrement();
224 auto const& x_curr = *x_new_timestep[process_id];
225
226 _M->setZero();
227 _K->setZero();
228 _b->setZero();
229
230 _ode.preAssemble(t, dt, x_curr);
231 _ode.assemble(t, dt, x_new_timestep, x_prev, process_id, *_M, *_K, *_b);
232
236}
237
240 NonlinearSolverTag::Picard>::computeKnownSolutions(GlobalVector const& x,
241 int const process_id)
242{
243 _known_solutions =
244 _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
245}
246
249 NonlinearSolverTag::Picard>::applyKnownSolutions(GlobalVector& x) const
250{
251 ::detail::applyKnownSolutions(_known_solutions, x);
252}
253
256 applyKnownSolutionsPicard(GlobalMatrix& A,
257 GlobalVector& rhs,
258 GlobalVector& x) const
259{
260 if (!_known_solutions)
261 {
262 return;
263 }
264
266 std::vector<IndexType> ids;
267 std::vector<double> values;
268 for (auto const& bc : *_known_solutions)
269 {
270 ids.insert(end(ids), begin(bc.ids), end(bc.ids));
271 values.insert(end(values), begin(bc.values), end(bc.values));
272 }
273 MathLib::applyKnownSolution(A, rhs, x, ids, values);
274}
275
276} // 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:198
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