OpenFPM  5.2.0
Project that contain the implementation of distributed structures
gs-odeint-cuda.cu
1 
2 // Include Vector Expression,Vector Expressions for Subset,DCPSE,Odeint header files
3 #include "Operators/Vector/vector_dist_operators.hpp"
4 #include "DCPSE/DCPSE_op/DCPSE_op.hpp"
5 #include "OdeIntegrators/OdeIntegrators.hpp"
6 #include <string>
8 
9 constexpr int x = 0;
10 constexpr int y = 1;
11 
12 double dt=1.0,tf=5000.0;
13 
14 void *PointerDistGlobal;
15 
18 
19 
20 template<typename laplacian_type, typename verletList_type>
21 struct RHSFunctor
22 {
23 
24  //Intializing the operators
25  laplacian_type &Lap;
26  verletList_type &verletList;
27 
28  // Physical contants
29  double K = 0.053;
30  double F = 0.014;
31 
32  double d1 = 2*1e-4;
33  double d2 = 1*1e-4;
34 
35  //Constructor
36  RHSFunctor(laplacian_type &Lap, verletList_type& verletList) : Lap(Lap), verletList(verletList)
37  {}
38 
39  void operator()( const state_type_2d_ofp &X , state_type_2d_ofp &dxdt , const double t ) const
40  {
41  //Casting the pointers to OpenFPM vector distributions
42  dist_vector_type &Particles= *(dist_vector_type *) PointerDistGlobal;
43 
44  //Aliasing the properties.
45  auto C = getV<0>(Particles);
46  //These expressions only update the bulk values of C.
47  C[x]=X.data.get<0>();
48  C[y]=X.data.get<1>();
49  Particles.ghost_get<0>(SKIP_LABELLING);
50  // Particles.updateVerlet(verletList,verletList.getRCut());
51 
52  // We do the RHS computations for the Laplacian and reaction term
53  // (Updating bulk only).
54  dxdt.data.get<0>() = d1*Lap(C[x]) - C[x] * C[y] * C[y] + F - F * C[x];
55  dxdt.data.get<1>() = d2*Lap(C[y]) + C[x] * C[y] * C[y] - (F+K) * C[y];
56  //We copy back to the dxdt state_type for Odeint
57  //=dC[x];
58  //=dC[y];
59  }
60 };
61 
62 struct ObserverFunctor {
63 
64  int ctr;
65  double t_old;
66 
67  //Constructor
68  ObserverFunctor() {
69  //a counter for counting the np. of steps
70  ctr = 0;
71  //Starting with t=0, we compute the step size take by t-t_old. So for the first observed step size is what we provide. Which will be 0-(-dt)=dt.
72  t_old = -dt;
73  }
74 
75  void operator()(state_type_2d_ofp &X,const double t) {
76  if (ctr % 300 == 0) {
77  dist_vector_type &Particles= *(dist_vector_type *) PointerDistGlobal;
78  auto C = getV<0>(Particles);
79  C[x]=X.data.get<0>();
80  C[y]=X.data.get<1>();
81  auto &v_cl=create_vcluster();
82  if(v_cl.rank()==0)
83  {
84  std::cout<<"Time: "<<t<<", "<<"dt: "<<t-t_old<<std::endl;
85  }
86  Particles.deleteGhost();
87  Particles.write_frame("PDE_sol",ctr,t);
88  Particles.ghost_get<0>();
89  }
90  t_old=t;
91  ctr++;
92  }
93 
94 };
95 
96 
97 template <typename stepper_type, typename laplacian_type, typename verletList_type>
98 void run_stepper_const(dist_vector_type &Particles, std::vector<double> &runtime_v, laplacian_type &Lap, verletList_type& verletList) {
99 
101  ObserverFunctor ObserveAndUpdate;
102  auto C = getV<0>(Particles);
103  auto InitC = getV<1>(Particles);
104 
106  x0.data.get<x>() = InitC[x];
107  x0.data.get<y>() = InitC[y];
108 
109  timer timer_integrate;
110  timer_integrate.start();
111  boost::numeric::odeint::integrate_const(stepper_type(), System, x0, 0.0, tf, dt, ObserveAndUpdate);
112 // boost::numeric::odeint::integrate_const(stepper_type(), derivative, x0, 0.0, tf, dt);
113  timer_integrate.stop();
114  double rt=timer_integrate.getwct();
115  auto &v_cl=create_vcluster();
116  v_cl.sum(rt);
117  v_cl.execute();
118 
119  runtime_v.push_back(rt);
120  C[x]=x0.data.get<x>();
121  C[y]=x0.data.get<y>();
122  if(v_cl.rank()==0)std::cout << "Runtime: " << rt << std::endl;
123 }
124 
125 template <typename stepper_type, typename laplacian_type, typename verletList_type>
126 void run_stepper_adaptive(dist_vector_type &Particles, std::vector<double> &runtime_v, laplacian_type &Lap, verletList_type& verletList) {
127 
129  ObserverFunctor ObserveAndUpdate;
130  auto C = getV<0>(Particles);
131  auto InitC = getV<1>(Particles);
132 
134  x0.data.get<x>() = InitC[x];
135  x0.data.get<y>() = InitC[y];
136 
137  timer timer_integrate;
138  timer_integrate.start();
139  boost::numeric::odeint::integrate_adaptive(boost::numeric::odeint::make_controlled(1e-3,1e-3,stepper_type()), System , x0 , 0.0 , tf , dt, ObserveAndUpdate);
140 // boost::numeric::odeint::integrate_adaptive( stepper_type() , derivative , x0 , 0.0 , tf , dt);
141  timer_integrate.stop();
142  double rt=timer_integrate.getwct();
143  auto &v_cl=create_vcluster();
144  v_cl.sum(rt);
145  v_cl.execute();
146  runtime_v.push_back(rt);
147 
148  C[x]=x0.data.get<x>();
149  C[y]=x0.data.get<y>();
150  if(v_cl.rank()==0)std::cout << "Runtime: " << rt << std::endl;
151 }
152 
153 double average(std::vector<double> &nums) {
154  return std::accumulate(nums.begin(), nums.end(), 0.0) / static_cast<double>(nums.size());
155 }
156 
157 
158 int main(int argc, char *argv[])
159 {
160  // initialize library
161  openfpm_init(&argc, &argv);
162  tf=std::atof(argv[2]);
163  // output
164  std::vector<double> runtime_rk4_const;
165  std::vector<double> runtime_rk5_const;
166  std::vector<double> runtime_rk78_const;
167  std::vector<double> runtime_rk5_adapt;
168  size_t gdsz=std::atof(argv[1]);
169  Box<3,double> box({0.0,0.0,0.0},{2.5,2.5,2.5});
170  size_t sz[3] = {gdsz,gdsz,gdsz};
171  // Define periodicity of the grid
172  size_t bc[3] = {PERIODIC,PERIODIC,PERIODIC};
173  double spacing[3];
174  spacing[0] = 2.5 / (sz[0]);
175  spacing[1] = 2.5 / (sz[1]);
176  spacing[2] = 2.5 / (sz[2]);
177  double rCut = 2.9 * spacing[0];
178  Ghost<3, double> ghost(rCut);
179 
180  dist_vector_type Particles(0, box, bc, ghost);
181  Particles.setPropNames({"Concentration","Initial"});
182 
183  auto it = Particles.getGridIterator(sz);
184  while (it.isNext()) {
185  Particles.add();
186  auto key = it.get();
187  double x = 0.0 + key.get(0) * spacing[0];
188  Particles.getLastPos()[0] = x;
189  double y = 0.0 + key.get(1) * spacing[1];
190  Particles.getLastPos()[1] = y;
191  double z = 0.0 + key.get(2) * spacing[2];
192  Particles.getLastPos()[2] = z;
193  // Here fill the Initial value of the concentration.
194  Particles.template getLastProp<1>()[0] = 1.0;
195  Particles.template getLastProp<1>()[1] = 0.0;
196 
197  if (x > 1.55 && x < 1.85 && y > 1.55 && y < 1.85 && z > 1.55 && z < 1.85) {
198  Particles.template getLastProp<1>()[0] = 0.5 + (((double)std::rand())/RAND_MAX -0.5)/10.0;
199  Particles.template getLastProp<1>()[1] = 0.25 + (((double)std::rand())/RAND_MAX -0.5)/20.0;
200  }
201 
202  ++it;
203  }
204  Particles.map();
205  Particles.ghost_get<0>();
206 
207 
208  // Now we initialize the grid with a filled circle. Outside the circle, the value of Phi_0 will be -1, inside +1.
209  //Now we construct the subsets based on the subset number.
210 
211  //We cast the global pointers to Particles and Particles_bulk as expected by the RHS functor.
212  PointerDistGlobal = (void *) &Particles;
213  auto verletList = Particles.template getVerlet<VL_NON_SYMMETRIC|VL_SKIP_REF_PART>(rCut);
214  //We create the DCPSE Based Laplacian operator.
215 
216  Laplacian_gpu<decltype(verletList)> Lap(Particles, verletList, 2, rCut, support_options::RADIUS);
217  // Laplacian<decltype(verletList)> Lap(Particles, verletList, 2, rCut, support_options::RADIUS);
218  auto C = getV<0>(Particles);
219  auto Init = getV<1>(Particles);
220  C=Init;
221  //Now we create a odeint stepper object (RK4). Since we are in 2d, we are going to use "state_type_2d_ofp". Which is a structure or state_type compatible with odeint. We further pass all the parameters including "boost::numeric::odeint::vector_space_algebra_ofp",which tell odeint to use openfpm algebra.
222  // The template parameters are: state_type_2d_ofp (state type of X), double (type of the value inside the state), state_type_2d_ofp (state type of DxDt), double (type of the time), boost::numeric::odeint::vector_space_algebra_ofp (our algebra)
223  typedef boost::numeric::odeint::runge_kutta4<state_type_2d_ofp, double, state_type_2d_ofp, double, boost::numeric::odeint::vector_space_algebra_ofp> Odeint_rk4;
224  typedef boost::numeric::odeint::runge_kutta_cash_karp54< state_type_2d_ofp,double,state_type_2d_ofp,double,boost::numeric::odeint::vector_space_algebra_ofp> Odeint_rk5;
225  typedef boost::numeric::odeint::runge_kutta_fehlberg78< state_type_2d_ofp,double,state_type_2d_ofp,double,boost::numeric::odeint::vector_space_algebra_ofp> Odeint_rk8;
226  //The method Odeint_rk4 from Odeint, requires system (a function which computes RHS of the PDE), an instance of the Compute RHS functor. We create the System with the correct types and parameteres for the operators as declared before.
227  RHSFunctor<Laplacian_gpu<decltype(verletList)>, decltype(verletList)> System(Lap, verletList);
228 
229  //Since we are using Odeint to control the time steps, we also create a observer instance. Which also updates the position via an euler step for moving thr particles.
230  ObserverFunctor ObserveAndUpdate;
231 
232 
233  //Furhter, odeint needs data in a state type "state_type_2d_ofp", we create one and fill in the initial condition.
235  //Since we created a 2d state_type we initialize the two fields in the object data using the method get.
236  X.data.get<x>() = C[0];
237  X.data.get<y>() = C[1];
238 
239 
240  std::vector<double> inter_times; // vector to store intermediate time steps taken by odeint.
241  Particles.deleteGhost();
242  Particles.write("Initial");
243  Particles.ghost_get<0>();
244  //for (int i = 0; i < 3; ++i) {
245  run_stepper_const<Odeint_rk4>(Particles, runtime_rk4_const,Lap, verletList);
246  Particles.deleteGhost();
247  Particles.write("RK4final");
248  Particles.ghost_get<0>();
249  /*run_stepper_const<Odeint_rk5>(Particles, runtime_rk5_const,Lap,verletList);
250  Particles.deleteGhost();
251  Particles.write("RK5final");
252  Particles.ghost_get<0>();
253  run_stepper_const<Odeint_rk8>(Particles, runtime_rk78_const,Lap,verletList);
254  Particles.deleteGhost();
255  Particles.write("RK8");
256  Particles.ghost_get<0>();
257  run_stepper_adaptive<Odeint_rk5>(Particles, runtime_rk5_adapt,Lap,verletList);
258  Particles.deleteGhost();
259  Particles.write("AdapRK5");
260  Particles.ghost_get<0>();
261  */
262 
263 
264  //}
265 
266 // size_t steps = boost::numeric::odeint::integrate_const(Odeint_rk4, System, X, 0.0, tf, dt, ObserveAndUpdate);
267 // size_t steps = boost::numeric::odeint::integrate_adaptive( boost::numeric::odeint::make_controlled( 1.0e-7 , 1.0e-7 , Odeint_rk5()) , System , X , 0.0 , tf , dt, ObserveAndUpdate );
268 
269  auto & vcl = create_vcluster();
270 
271  if (vcl.getProcessUnitID() == 0) {
272 
273  std::ofstream output_runtime;
274 
275  // head line
276  std::string headline = "cores,rk4,dopri5,fehlberg78,dopri5 adaptive";
277  if (access("runtime.csv", F_OK) == -1) {
278  output_runtime.open("runtime.csv");
279  output_runtime << headline << "\n";
280  } else {
281  output_runtime.open("runtime.csv", std::ios::app);
282  }
283 
284  output_runtime << vcl.getProcessingUnits()
285  << "," << average(runtime_rk4_const)
286  << "," << average(runtime_rk5_const)
287  << "," << average(runtime_rk78_const)
288  << "," << average(runtime_rk5_adapt)
289  << "\n";
290  }
291 
292  //Deallocating the operators
293  Lap.deallocate(Particles);
294  openfpm_finalize(); // Finalize openFPM library
295  return 0;
296 } //main end
This class represent an N-dimensional box.
Definition: Box.hpp:60
Definition: Ghost.hpp:40
Laplacian second order on h (spacing)
Definition: Laplacian.hpp:23
This class implement the point shape in an N-dimensional space.
Definition: Point.hpp:28
Class for cpu time benchmarking.
Definition: timer.hpp:28
void start()
Start the timer.
Definition: timer.hpp:90
Distributed vector.
void setPropNames(const openfpm::vector< std::string > &names)
Set the properties names.
grid_dist_id_iterator_dec< Decomposition > getGridIterator(const size_t(&sz)[dim])
void ghost_get(size_t opt=WITH_POSITION)
It synchronize the properties and position of the ghost particles.
void map(size_t opt=NONE)
It move all the particles that does not belong to the local processor to the respective processor.
void add()
Add local particle.
void deleteGhost()
Delete the particles on the ghost.
auto getLastPos() -> decltype(vPos.template get< 0 >(0))
Get the position of the last element.
[v_transform metafunction]
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
Definition: aggregate.hpp:221
A 2d Odeint and Openfpm compatible structure.