Reference documentation for deal.II version 9.1.0-pre
arkode.cc
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 //-----------------------------------------------------------
15 
16 
17 #include <deal.II/base/config.h>
18 
19 #include <deal.II/sundials/arkode.h>
20 
21 #ifdef DEAL_II_WITH_SUNDIALS
22 
23 # include <deal.II/base/utilities.h>
24 
25 # include <deal.II/lac/block_vector.h>
26 # ifdef DEAL_II_WITH_TRILINOS
27 # include <deal.II/lac/trilinos_parallel_block_vector.h>
28 # include <deal.II/lac/trilinos_vector.h>
29 # endif
30 # ifdef DEAL_II_WITH_PETSC
31 # include <deal.II/lac/petsc_block_vector.h>
32 # include <deal.II/lac/petsc_vector.h>
33 # endif
34 # include <deal.II/base/utilities.h>
35 
36 # include <deal.II/sundials/copy.h>
37 
38 # include <arkode/arkode_impl.h>
39 # include <sundials/sundials_config.h>
40 
41 # include <iomanip>
42 # include <iostream>
43 
44 // Make sure we know how to call sundials own ARKode() function
45 const auto &SundialsARKode = ARKode;
46 
47 DEAL_II_NAMESPACE_OPEN
48 
49 namespace SUNDIALS
50 {
51  using namespace internal;
52 
53  namespace
54  {
55  template <typename VectorType>
56  int
57  t_arkode_explicit_function(realtype tt,
58  N_Vector yy,
59  N_Vector yp,
60  void * user_data)
61  {
62  ARKode<VectorType> &solver =
63  *static_cast<ARKode<VectorType> *>(user_data);
65 
66  typename VectorMemory<VectorType>::Pointer src_yy(mem);
67  solver.reinit_vector(*src_yy);
68 
69  typename VectorMemory<VectorType>::Pointer dst_yp(mem);
70  solver.reinit_vector(*dst_yp);
71 
72  copy(*src_yy, yy);
73 
74  int err = solver.explicit_function(tt, *src_yy, *dst_yp);
75 
76  copy(yp, *dst_yp);
77 
78  return err;
79  }
80 
81 
82 
83  template <typename VectorType>
84  int
85  t_arkode_implicit_function(realtype tt,
86  N_Vector yy,
87  N_Vector yp,
88  void * user_data)
89  {
90  ARKode<VectorType> &solver =
91  *static_cast<ARKode<VectorType> *>(user_data);
93 
94  typename VectorMemory<VectorType>::Pointer src_yy(mem);
95  solver.reinit_vector(*src_yy);
96 
97  typename VectorMemory<VectorType>::Pointer dst_yp(mem);
98  solver.reinit_vector(*dst_yp);
99 
100  copy(*src_yy, yy);
101 
102  int err = solver.implicit_function(tt, *src_yy, *dst_yp);
103 
104  copy(yp, *dst_yp);
105 
106  return err;
107  }
108 
109 
110 
111  template <typename VectorType>
112  int
113  t_arkode_setup_jacobian(ARKodeMem arkode_mem,
114  int convfail,
115  N_Vector ypred,
116  N_Vector fpred,
117  booleantype *jcurPtr,
118  N_Vector,
119  N_Vector,
120  N_Vector)
121  {
122  ARKode<VectorType> &solver =
123  *static_cast<ARKode<VectorType> *>(arkode_mem->ark_user_data);
125 
126  typename VectorMemory<VectorType>::Pointer src_ypred(mem);
127  solver.reinit_vector(*src_ypred);
128 
129  typename VectorMemory<VectorType>::Pointer src_fpred(mem);
130  solver.reinit_vector(*src_fpred);
131 
132  copy(*src_ypred, ypred);
133  copy(*src_fpred, fpred);
134 
135  int err = solver.setup_jacobian(convfail,
136  arkode_mem->ark_tn,
137  arkode_mem->ark_gamma,
138  *src_ypred,
139  *src_fpred,
140  (bool &)*jcurPtr);
141 
142  return err;
143  }
144 
145 
146 
147  template <typename VectorType>
148  int
149  t_arkode_solve_jacobian(ARKodeMem arkode_mem,
150  N_Vector b,
151 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
152  N_Vector,
153 # endif
154  N_Vector ycur,
155  N_Vector fcur)
156  {
157  ARKode<VectorType> &solver =
158  *static_cast<ARKode<VectorType> *>(arkode_mem->ark_user_data);
160 
161  typename VectorMemory<VectorType>::Pointer src(mem);
162  solver.reinit_vector(*src);
163 
164  typename VectorMemory<VectorType>::Pointer src_ycur(mem);
165  solver.reinit_vector(*src_ycur);
166 
167  typename VectorMemory<VectorType>::Pointer src_fcur(mem);
168  solver.reinit_vector(*src_fcur);
169 
170  typename VectorMemory<VectorType>::Pointer dst(mem);
171  solver.reinit_vector(*dst);
172 
173  copy(*src, b);
174  copy(*src_ycur, ycur);
175  copy(*src_fcur, fcur);
176 
177  int err = solver.solve_jacobian_system(arkode_mem->ark_tn,
178  arkode_mem->ark_gamma,
179  *src_ycur,
180  *src_fcur,
181  *src,
182  *dst);
183  copy(b, *dst);
184 
185  return err;
186  }
187 
188 
189 
190  template <typename VectorType>
191  int
192  t_arkode_setup_mass(ARKodeMem arkode_mem, N_Vector, N_Vector, N_Vector)
193  {
194  ARKode<VectorType> &solver =
195  *static_cast<ARKode<VectorType> *>(arkode_mem->ark_user_data);
196  int err = solver.setup_mass(arkode_mem->ark_tn);
197  return err;
198  }
199 
200 
201 
202  template <typename VectorType>
203  int
204  t_arkode_solve_mass(ARKodeMem arkode_mem,
205 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
206  N_Vector b,
207  N_Vector
208 # else
209  N_Vector b
210 # endif
211  )
212  {
213  ARKode<VectorType> &solver =
214  *static_cast<ARKode<VectorType> *>(arkode_mem->ark_user_data);
216 
217  typename VectorMemory<VectorType>::Pointer src(mem);
218  solver.reinit_vector(*src);
219 
220  typename VectorMemory<VectorType>::Pointer dst(mem);
221  solver.reinit_vector(*dst);
222 
223  copy(*src, b);
224 
225  int err = solver.solve_mass_system(*src, *dst);
226  copy(b, *dst);
227 
228  return err;
229  }
230  } // namespace
231 
232  template <typename VectorType>
234  const MPI_Comm mpi_comm)
235  : data(data)
236  , arkode_mem(nullptr)
237  , yy(nullptr)
238  , abs_tolls(nullptr)
239  , communicator(is_serial_vector<VectorType>::value ?
240  MPI_COMM_SELF :
241  Utilities::MPI::duplicate_communicator(mpi_comm))
242  {
244  }
245 
246  template <typename VectorType>
248  {
249  if (arkode_mem)
250  ARKodeFree(&arkode_mem);
251 # ifdef DEAL_II_WITH_MPI
253  {
254  const int ierr = MPI_Comm_free(&communicator);
255  (void)ierr;
256  AssertNothrow(ierr == MPI_SUCCESS, ExcMPI(ierr));
257  }
258 # endif
259  }
260 
261 
262 
263  template <typename VectorType>
264  unsigned int
265  ARKode<VectorType>::solve_ode(VectorType &solution)
266  {
267  unsigned int system_size = solution.size();
268 
269  double t = data.initial_time;
270  double h = data.initial_step_size;
271  unsigned int step_number = 0;
272 
273  int status;
274  (void)status;
275 
276  // The solution is stored in
277  // solution. Here we take only a
278  // view of it.
279 # ifdef DEAL_II_WITH_MPI
281  {
282  const IndexSet is = solution.locally_owned_elements();
283  const size_t local_system_size = is.n_elements();
284 
285  yy = N_VNew_Parallel(communicator, local_system_size, system_size);
286 
287  abs_tolls =
288  N_VNew_Parallel(communicator, local_system_size, system_size);
289  }
290  else
291 # endif
292  {
295  "Trying to use a serial code with a parallel vector."));
296  yy = N_VNew_Serial(system_size);
297  abs_tolls = N_VNew_Serial(system_size);
298  }
300 
301  double next_time = data.initial_time;
302 
303  if (output_step)
304  output_step(0, solution, 0);
305 
306  while (t < data.final_time)
307  {
308  next_time += data.output_period;
309 
310  status = SundialsARKode(arkode_mem, next_time, yy, &t, ARK_NORMAL);
311 
312  AssertARKode(status);
313 
314  status = ARKodeGetLastStep(arkode_mem, &h);
315  AssertARKode(status);
316 
317  copy(solution, yy);
318 
319  while (solver_should_restart(t, solution))
320  reset(t, h, solution);
321 
322  step_number++;
323 
324  if (output_step)
325  output_step(t, solution, step_number);
326  }
327 
328  // Free the vectors which are no longer used.
329 # ifdef DEAL_II_WITH_MPI
331  {
332  N_VDestroy_Parallel(yy);
333  N_VDestroy_Parallel(abs_tolls);
334  }
335  else
336 # endif
337  {
338  N_VDestroy_Serial(yy);
339  N_VDestroy_Serial(abs_tolls);
340  }
341 
342  return step_number;
343  }
344 
345  template <typename VectorType>
346  void
347  ARKode<VectorType>::reset(const double & current_time,
348  const double & current_time_step,
349  const VectorType &solution)
350  {
351  unsigned int system_size;
352 
353  if (arkode_mem)
354  ARKodeFree(&arkode_mem);
355 
356  arkode_mem = ARKodeCreate();
357 
358  // Free the vectors which are no longer used.
359  if (yy)
360  {
361 # ifdef DEAL_II_WITH_MPI
363  {
364  N_VDestroy_Parallel(yy);
365  N_VDestroy_Parallel(abs_tolls);
366  }
367  else
368 # endif
369  {
370  N_VDestroy_Serial(yy);
371  N_VDestroy_Serial(abs_tolls);
372  }
373  }
374 
375  int status;
376  (void)status;
377  system_size = solution.size();
378 # ifdef DEAL_II_WITH_MPI
380  {
381  const IndexSet is = solution.locally_owned_elements();
382  const size_t local_system_size = is.n_elements();
383 
384  yy = N_VNew_Parallel(communicator, local_system_size, system_size);
385 
386  abs_tolls =
387  N_VNew_Parallel(communicator, local_system_size, system_size);
388  }
389  else
390 # endif
391  {
392  yy = N_VNew_Serial(system_size);
393  abs_tolls = N_VNew_Serial(system_size);
394  }
395 
396  copy(yy, solution);
397 
399  ExcFunctionNotProvided("explicit_function || implicit_function"));
400 
401  status = ARKodeInit(
402  arkode_mem,
403  explicit_function ? &t_arkode_explicit_function<VectorType> : nullptr,
404  implicit_function ? &t_arkode_implicit_function<VectorType> : nullptr,
405  current_time,
406  yy);
407  AssertARKode(status);
408 
410  {
412  status =
413  ARKodeSVtolerances(arkode_mem, data.relative_tolerance, abs_tolls);
414  AssertARKode(status);
415  }
416  else
417  {
418  status = ARKodeSStolerances(arkode_mem,
421  AssertARKode(status);
422  }
423 
424  status = ARKodeSetInitStep(arkode_mem, current_time_step);
425  AssertARKode(status);
426 
427  status = ARKodeSetUserData(arkode_mem, (void *)this);
428  AssertARKode(status);
429 
430  status = ARKodeSetStopTime(arkode_mem, data.final_time);
431  AssertARKode(status);
432 
433  status =
434  ARKodeSetMaxNonlinIters(arkode_mem, data.maximum_non_linear_iterations);
435  AssertARKode(status);
436 
437  // Initialize solver
438  ARKodeMem ARKode_mem = (ARKodeMem)arkode_mem;
439 
441  {
442  status = ARKodeSetNewton(arkode_mem);
443  AssertARKode(status);
445  {
446  status = ARKodeSetLinear(
447  arkode_mem, data.implicit_function_is_time_independent ? 0 : 1);
448  AssertARKode(status);
449  }
450 
451 
452  ARKode_mem->ark_lsolve = t_arkode_solve_jacobian<VectorType>;
453  if (setup_jacobian)
454  {
455  ARKode_mem->ark_lsetup = t_arkode_setup_jacobian<VectorType>;
456 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
457  ARKode_mem->ark_setupNonNull = true;
458 # endif
459  }
460  }
461  else
462  {
463  status =
464  ARKodeSetFixedPoint(arkode_mem, data.maximum_non_linear_iterations);
465  AssertARKode(status);
466  }
467 
468 
469  if (solve_mass_system)
470  {
471  ARKode_mem->ark_msolve = t_arkode_solve_mass<VectorType>;
472 
473  if (setup_mass)
474  {
475  ARKode_mem->ark_msetup = t_arkode_setup_mass<VectorType>;
476 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
477  ARKode_mem->ark_MassSetupNonNull = true;
478 # endif
479  }
480  }
481 
482  status = ARKodeSetOrder(arkode_mem, data.maximum_order);
483  AssertARKode(status);
484  }
485 
486  template <typename VectorType>
487  void
489  {
490  reinit_vector = [](VectorType &) {
491  AssertThrow(false, ExcFunctionNotProvided("reinit_vector"));
492  };
493 
494  solver_should_restart = [](const double, VectorType &) -> bool {
495  return false;
496  };
497  }
498 
499  template class ARKode<Vector<double>>;
500  template class ARKode<BlockVector<double>>;
501 
502 # ifdef DEAL_II_WITH_MPI
503 
504 # ifdef DEAL_II_WITH_TRILINOS
507 # endif // DEAL_II_WITH_TRILINOS
508 
509 # ifdef DEAL_II_WITH_PETSC
510 # ifndef PETSC_USE_COMPLEX
511  template class ARKode<PETScWrappers::MPI::Vector>;
513 # endif // PETSC_USE_COMPLEX
514 # endif // DEAL_II_WITH_PETSC
515 
516 # endif // DEAL_II_WITH_MPI
517 
518 } // namespace SUNDIALS
519 
520 DEAL_II_NAMESPACE_CLOSE
521 
522 #endif
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1278
size_type n_elements() const
Definition: index_set.h:1743
std::function< int(const VectorType &rhs, VectorType &dst)> solve_mass_system
Definition: arkode.h:787
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
MPI_Comm communicator
Definition: arkode.h:888
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
Definition: arkode.h:806
void reset(const double &t, const double &h, const VectorType &y)
Definition: arkode.cc:347
N_Vector abs_tolls
Definition: arkode.h:881
std::function< VectorType &()> get_local_tolerances
Definition: arkode.h:833
unsigned int maximum_non_linear_iterations
Definition: arkode.h:486
std::function< int(const double t)> setup_mass
Definition: arkode.h:766
std::function< bool(const double t, VectorType &sol)> solver_should_restart
Definition: arkode.h:825
std::function< int(const double t, const double gamma, const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
Definition: arkode.h:730
AdditionalData data
Definition: arkode.h:866
#define Assert(cond, exc)
Definition: exceptions.h:1227
std::function< int(const double t, const VectorType &y, VectorType &res)> implicit_function
Definition: arkode.h:586
Definition: cuda.h:32
std::function< int(const double t, const VectorType &y, VectorType &explicit_f)> explicit_function
Definition: arkode.h:567
std::function< void(VectorType &)> reinit_vector
Definition: arkode.h:547
unsigned int solve_ode(VectorType &solution)
Definition: arkode.cc:265
std::function< int(const int convfail, const double t, const double gamma, const VectorType &ypred, const VectorType &fpred, bool &j_is_current)> setup_jacobian
Definition: arkode.h:678
void set_functions_to_trigger_an_assert()
Definition: arkode.cc:488
ARKode(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
Definition: arkode.cc:233
static::ExceptionBase & ExcInternalError()