accelerInt  v0.1
solver_rk78.cpp
Go to the documentation of this file.
1 
11 //wrapper code
12 #include "rk78_typedefs.hpp"
13 
14 //boost includes
15 #include <boost/numeric/odeint/integrate/integrate_adaptive.hpp>
16 extern "C" {
17 #include "solver.h"
18 }
19 
20 #ifdef GENERATE_DOCS
21 namespace rk78 {
22 #endif
23 
24 extern std::vector<state_type*> state_vectors;
25 extern std::vector<rhs_eval*> evaluators;
26 extern std::vector<controller> controllers;
27 
28 extern "C" void intDriver(const int, const double, const double, const double*, double*);
29 
30 #ifdef STIFFNESS_MEASURE
31  extern std::vector<double> max_stepsize;
32  controlled_step_result test_step(int index, const state_type& y, const double t, state_type& y_out, const double dt)
33  {
34  try
35  {
36  double t_copy = t;
37  double dt_copy = dt;
38  return controllers[index]->try_step(*evaluators[index], y, t_copy, y_out, dt_copy);
39  }
40  catch(...)
41  {
42  return fail;
43  }
44  }
45 #endif
46 
57 void intDriver (const int NUM, const double t, const double t_end,
58  const double *pr_global, double *y_global)
59 {
60  #ifdef STIFFNESS_MEASURE
61  max_stepsize.clear();
62  max_stepsize.resize(NUM, 0.0);
63  #endif
64 
65  int tid = 0;
66 #ifdef STIFFNESS_MEASURE
67  #pragma omp parallel for shared(state_vectors, evaluators, controllers, max_stepsize) private(tid)
68 #else
69  #pragma omp parallel for shared(state_vectors, evaluators, controllers) private(tid)
70 #endif
71  for (tid = 0; tid < NUM; ++tid) {
72  int index = omp_get_thread_num();
73 
74  // local array with initial values
75  state_type& vec = *state_vectors[index];
76  evaluators[index]->set_state_var(pr_global[tid]);
77 
78  // load local array with initial values from global array
79  for (int i = 0; i < NSP; i++)
80  {
81  vec[i] = y_global[tid + i * NUM];
82  }
83 
84 #ifndef STIFFNESS_MEASURE
85  integrate_adaptive(controllers[index],
86  *evaluators[index], vec, t, t_end, t_end - t);
87 #else
88  double tol = 1e-15;
89  state_type y_copy(vec);
90  //do a binary search to find the maximum stepsize
91  double left_size = 1.0;
92  while (test_step(index, vec, t, y_copy, left_size) == success)
93  {
94  left_size *= 10.0;
95  }
96  double right_size = 1e-20;
97  while (test_step(index, vec, t, y_copy, right_size) == fail)
98  {
99  right_size /= 10.0;
100  }
101  double delta = 1.0;
102  double mid = 0;
103  while (delta > tol) {
104  mid = (left_size + right_size) / 2.0;
105  controlled_step_result result = test_step(index, vec, t, y_copy, mid);
106  if (result == fail) {
107  //mid becomes the new left
108  delta = fabs(left_size - mid) / left_size;
109  left_size = mid;
110  }
111  else{
112  delta = fabs(right_size - mid) / right_size;
113  right_size = mid;
114  }
115  }
116  max_stepsize[tid] = mid;
117 #endif
118 
119  // update global array with integrated values
120  for (int i = 0; i < NSP; i++)
121  {
122  y_global[tid + i * NUM] = vec[i];
123  }
124 
125  }
126 }
127 
128 #ifdef GENERATE_DOCS
129 }
130 #endif
std::vector< controller > controllers
ODE controllers.
Definition: rk78_init.cpp:24
std::vector< state_type * > state_vectors
State vector containers for boost.
Definition: rk78_init.cpp:18
std::vector< rhs_eval * > evaluators
RHS wrappers for boost.
Definition: rk78_init.cpp:20
#define NSP
The IVP system size.
Definition: header.cuh:20
std::vector< double > state_type
state vector
Contains skeleton of all methods that need to be defined on a per solver basis.
Defines an interface for boost&#39;s runge_kutta_fehlberg78 solver.
void intDriver(const int, const double, const double, const double *, double *)
Integration driver for the CPU integrators.
Definition: solver_rk78.cpp:57