OGS
TimeDiscretizedODESystem.cpp
Go to the documentation of this file.
1 
12 
15 #include "NumLib/Exceptions.h"
17 
18 namespace detail
19 {
21 template <typename Solutions, typename Vector>
22 void 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 
43 namespace NumLib
44 {
45 TimeDiscretizedODESystem<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 
65  NonlinearSolverTag::Newton>::~TimeDiscretizedODESystem()
66 {
71 }
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  auto const dxdot_dx = _time_disc.getNewXWeight();
85 
86  std::vector<GlobalVector*> xdot(x_new_timestep.size());
87  _xdot_ids.resize(x_new_timestep.size());
88  for (std::size_t i = 0; i < xdot.size(); i++)
89  {
90  xdot[i] =
92  _time_disc.getXdot(*x_new_timestep[i], *x_prev[i], *xdot[i]);
93  }
94 
95  _M->setZero();
96  _K->setZero();
97  _b->setZero();
98  _Jac->setZero();
99 
100  _ode.preAssemble(t, dt, x_curr);
101  try
102  {
103  _ode.assembleWithJacobian(t, dt, x_new_timestep, xdot, dxdot_dx, 1.0,
104  process_id, *_M, *_K, *_b, *_Jac);
105  }
106  catch (AssemblyException const&)
107  {
108  for (auto& v : xdot)
109  {
111  }
112  throw;
113  }
114 
119 
120  for (auto& v : xdot)
121  {
123  }
124 }
125 
128  NonlinearSolverTag::Newton>::getResidual(GlobalVector const& x_new_timestep,
129  GlobalVector const& x_prev,
130  GlobalVector& res) const
131 {
132  // TODO Maybe the duplicate calculation of xdot here and in assembleJacobian
133  // can be optimized. However, that would make the interface a bit more
134  // fragile.
135  auto& xdot = NumLib::GlobalVectorProvider::provider.getVector(_xdot_id);
136  _time_disc.getXdot(x_new_timestep, x_prev, xdot);
137 
138  _mat_trans->computeResidual(*_M, *_K, *_b, x_new_timestep, xdot, res);
139 
141 }
142 
145  NonlinearSolverTag::Newton>::getJacobian(GlobalMatrix& Jac) const
146 {
147  _mat_trans->computeJacobian(*_Jac, Jac);
148 }
149 
152  NonlinearSolverTag::Newton>::computeKnownSolutions(GlobalVector const& x,
153  int const process_id)
154 {
155  _known_solutions =
156  _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
157 }
158 
162 {
163  ::detail::applyKnownSolutions(_known_solutions, x);
164 }
165 
168  applyKnownSolutionsNewton(GlobalMatrix& Jac, GlobalVector& res,
169  GlobalVector& minus_delta_x) const
170 {
171  if (!_known_solutions)
172  {
173  return;
174  }
175 
177  std::vector<IndexType> ids;
178  for (auto const& bc : *_known_solutions)
179  {
180  ids.insert(end(ids), begin(bc.ids), end(bc.ids));
181  }
182 
183  // For the Newton method the values must be zero
184  std::vector<double> values(ids.size(), 0);
185  MathLib::applyKnownSolution(Jac, res, minus_delta_x, ids, values);
186 }
187 
190  TimeDiscretizedODESystem(const int process_id, ODE& ode,
191  TimeDisc& time_discretization)
192  : _ode(ode),
193  _time_disc(time_discretization),
194  _mat_trans(createMatrixTranslator<ODETag>(time_discretization))
195 {
197  ode.getMatrixSpecifications(process_id), _M_id);
199  ode.getMatrixSpecifications(process_id), _K_id);
201  ode.getMatrixSpecifications(process_id), _b_id);
202 }
203 
206  NonlinearSolverTag::Picard>::~TimeDiscretizedODESystem()
207 {
211 }
212 
215  assemble(std::vector<GlobalVector*> const& x_new_timestep,
216  std::vector<GlobalVector*> const& x_prev,
217  int const process_id)
218 {
219  namespace LinAlg = MathLib::LinAlg;
220 
221  auto const t = _time_disc.getCurrentTime();
222  auto const dt = _time_disc.getCurrentTimeIncrement();
223  auto const& x_curr = *x_new_timestep[process_id];
224  std::vector<GlobalVector*> xdot(x_new_timestep.size());
225  _xdot_ids.resize(x_new_timestep.size());
226 
227  for (std::size_t i = 0; i < xdot.size(); i++)
228  {
229  xdot[i] =
231  _time_disc.getXdot(*x_new_timestep[i], *x_prev[i], *xdot[i]);
232  }
233 
234  _M->setZero();
235  _K->setZero();
236  _b->setZero();
237 
238  _ode.preAssemble(t, dt, x_curr);
239  _ode.assemble(t, dt, x_new_timestep, xdot, process_id, *_M, *_K, *_b);
240 
244 
245  for (auto& v : xdot)
246  {
248  }
249 }
250 
253  NonlinearSolverTag::Picard>::computeKnownSolutions(GlobalVector const& x,
254  int const process_id)
255 {
256  _known_solutions =
257  _ode.getKnownSolutions(_time_disc.getCurrentTime(), x, process_id);
258 }
259 
263 {
264  ::detail::applyKnownSolutions(_known_solutions, x);
265 }
266 
269  applyKnownSolutionsPicard(GlobalMatrix& A,
270  GlobalVector& rhs,
271  GlobalVector& x) const
272 {
273  if (!_known_solutions)
274  {
275  return;
276  }
277 
279  std::vector<IndexType> ids;
280  std::vector<double> values;
281  for (auto const& bc : *_known_solutions)
282  {
283  ids.insert(end(ids), begin(bc.ids), end(bc.ids));
284  values.insert(end(values), begin(bc.values), end(bc.values));
285  }
286  MathLib::applyKnownSolution(A, rhs, x, ids, values);
287 }
288 
289 } // namespace NumLib
Global vector based on Eigen vector.
Definition: EigenVector.h:26
virtual GlobalMatrix & getMatrix(std::size_t &id)=0
Get an uninitialized matrix with the given id.
virtual void releaseMatrix(GlobalMatrix const &A)=0
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)
std::array< double, 3 > Vector
Definition: VariableType.h:28
void finalizeAssembly(PETScMatrix &A)
Definition: LinAlg.cpp:162
void setVector(PETScVector &v, std::initializer_list< double > values)
void applyKnownSolution(EigenMatrix &A, EigenVector &b, EigenVector &, const std::vector< EigenMatrix::IndexType > &vec_knownX_id, const std::vector< double > &vec_knownX_x, double)
Definition: EigenTools.cpp:17
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