OpenFPM  5.2.0
Project that contain the implementation of distributed structures
gs-odeint.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(
37  laplacian_type &Lap,
38  verletList_type& verletList
39  ) :
40  Lap(Lap),
41  verletList(verletList)
42  {}
43 
44  void operator()( const state_type_2d_ofp_gpu &X , state_type_2d_ofp_gpu &dxdt , const double t ) const
45  {
46  //Casting the pointers to OpenFPM vector distributions
47  dist_vector_type &Particles= *(dist_vector_type *) PointerDistGlobal;
48 
49  //Aliasing the properties.
50  auto C = getV<0>(Particles);
51  //These expressions only update the bulk values of C.
52  C[x]=X.data.get<0>();
53  C[y]=X.data.get<1>();
54  // Particles.ghost_get<0>(SKIP_LABELLING);
55  // Particles.updateVerlet(verletList,verletList.getRCut());
56 
57  // We do the RHS computations for the Laplacian and reaction term
58  // (Updating bulk only).
59  dxdt.data.get<0>() = d1*Lap(C[x]) - C[x] * C[y] * C[y] + F - F * C[x];
60  dxdt.data.get<1>() = d2*Lap(C[y]) + C[x] * C[y] * C[y] - (F+K) * C[y];
61  //We copy back to the dxdt state_type for Odeint
62  //=dC[x];
63  //=dC[y];
64  }
65 };
66 
67 template<typename dist_vector_type_ker>
68 struct ObserverFunctor {
69 
70  int ctr;
71  double t_old;
72 
73  dist_vector_type_ker &vectorDistKer;
74 
75  //Constructor
76  ObserverFunctor(dist_vector_type_ker &vectorDistKer)
77  : vectorDistKer(vectorDistKer)
78  {
79  //a counter for counting the np. of steps
80  ctr = 0;
81  //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.
82  t_old = -dt;
83  }
84 
85  void operator()(state_type_2d_ofp &X,const double t) {
86  dist_vector_type &Particles= *(dist_vector_type *) PointerDistGlobal;
87  Particles.deviceToHostProp<0,1>();
88  Particles.ghost_get<0>();
89 
90  if (ctr % 300 == 0) {
91  auto C = getV<0>(Particles);
92  C[x]=X.data.get<0>();
93  C[y]=X.data.get<1>();
94  auto &v_cl=create_vcluster();
95  if(v_cl.rank()==0)
96  {
97  std::cout<<"Time: "<<t<<", "<<"dt: "<<t-t_old<<std::endl;
98  }
99  Particles.deleteGhost();
100  Particles.write_frame("PDE_sol",ctr,t);
101  Particles.ghost_get<0>();
102  }
103  t_old=t;
104  ctr++;
105 
106  Particles.hostToDeviceProp<0,1>();
107  }
108 
109 };
110 
111 
112 template <typename stepper_type, typename laplacian_type, typename verletList_type>
113 void run_stepper_const(dist_vector_type &Particles, std::vector<double> &runtime_v, laplacian_type &Lap, verletList_type& verletList) {
114 
115  auto vectorDistKer = vd.toKernel();
117  ObserverFunctor ObserveAndUpdate(vectorDistKer);
118  auto C = getV<0>(Particles);
119  auto InitC = getV<1>(Particles);
120 
122  x0.data.get<x>() = InitC[x];
123  x0.data.get<y>() = InitC[y];
124 
125  Particles.hostToDeviceProp<0,1>();
126  timer timer_integrate;
127  timer_integrate.start();
128  boost::numeric::odeint::integrate_const(stepper_type(), System, x0, 0.0, tf, dt, ObserveAndUpdate);
129 // boost::numeric::odeint::integrate_const(stepper_type(), derivative, x0, 0.0, tf, dt);
130  timer_integrate.stop();
131  double rt=timer_integrate.getwct();
132  auto &v_cl=create_vcluster();
133  v_cl.sum(rt);
134  v_cl.execute();
135 
136  runtime_v.push_back(rt);
137  C[x]=x0.data.get<x>();
138  C[y]=x0.data.get<y>();
139  if(v_cl.rank()==0)std::cout << "Runtime: " << rt << std::endl;
140 }
141 
142 // template <typename stepper_type, typename laplacian_type, typename verletList_type>
143 // void run_stepper_adaptive(dist_vector_type &Particles, std::vector<double> &runtime_v, laplacian_type &Lap, verletList_type& verletList) {
144 
145 // RHSFunctor<laplacian_type, verletList_type> System(Lap, verletList);
146 // ObserverFunctor ObserveAndUpdate;
147 // auto C = getV<0>(Particles);
148 // auto InitC = getV<1>(Particles);
149 
150 // state_type_2d_ofp x0;
151 // x0.data.get<x>() = InitC[x];
152 // x0.data.get<y>() = InitC[y];
153 
154 // timer timer_integrate;
155 // timer_integrate.start();
156 // boost::numeric::odeint::integrate_adaptive(boost::numeric::odeint::make_controlled(1e-3,1e-3,stepper_type()), System , x0 , 0.0 , tf , dt, ObserveAndUpdate);
157 // // boost::numeric::odeint::integrate_adaptive( stepper_type() , derivative , x0 , 0.0 , tf , dt);
158 // timer_integrate.stop();
159 // double rt=timer_integrate.getwct();
160 // auto &v_cl=create_vcluster();
161 // v_cl.sum(rt);
162 // v_cl.execute();
163 // runtime_v.push_back(rt);
164 
165 // C[x]=x0.data.get<x>();
166 // C[y]=x0.data.get<y>();
167 // if(v_cl.rank()==0)std::cout << "Runtime: " << rt << std::endl;
168 // }
169 
170 double average(std::vector<double> &nums) {
171  return std::accumulate(nums.begin(), nums.end(), 0.0) / static_cast<double>(nums.size());
172 }
173 
174 
175 int main(int argc, char *argv[])
176 {
177  // initialize library
178  openfpm_init(&argc, &argv);
179  tf=std::atof(argv[2]);
180  // output
181  std::vector<double> runtime_rk4_const;
182  std::vector<double> runtime_rk5_const;
183  std::vector<double> runtime_rk78_const;
184  std::vector<double> runtime_rk5_adapt;
185  size_t gdsz=std::atof(argv[1]);
186  Box<3,double> box({0.0,0.0,0.0},{2.5,2.5,2.5});
187  size_t sz[3] = {gdsz,gdsz,gdsz};
188  // Define periodicity of the grid
189  size_t bc[3] = {PERIODIC,PERIODIC,PERIODIC};
190  double spacing[3];
191  spacing[0] = 2.5 / (sz[0]);
192  spacing[1] = 2.5 / (sz[1]);
193  spacing[2] = 2.5 / (sz[2]);
194  double rCut = 2.9 * spacing[0];
195  Ghost<3, double> ghost(rCut);
196 
197  dist_vector_type Particles(0, box, bc, ghost);
198  Particles.setPropNames({"Concentration","Initial"});
199 
200  auto it = Particles.getGridIterator(sz);
201  while (it.isNext()) {
202  Particles.add();
203  auto key = it.get();
204  double x = 0.0 + key.get(0) * spacing[0];
205  Particles.getLastPos()[0] = x;
206  double y = 0.0 + key.get(1) * spacing[1];
207  Particles.getLastPos()[1] = y;
208  double z = 0.0 + key.get(2) * spacing[2];
209  Particles.getLastPos()[2] = z;
210  // Here fill the Initial value of the concentration.
211  Particles.template getLastProp<1>()[0] = 1.0;
212  Particles.template getLastProp<1>()[1] = 0.0;
213 
214  if (x > 1.55 && x < 1.85 && y > 1.55 && y < 1.85 && z > 1.55 && z < 1.85) {
215  Particles.template getLastProp<1>()[0] = 0.5 + (((double)std::rand())/RAND_MAX -0.5)/10.0;
216  Particles.template getLastProp<1>()[1] = 0.25 + (((double)std::rand())/RAND_MAX -0.5)/20.0;
217  }
218 
219  ++it;
220  }
221  Particles.map();
222  Particles.ghost_get<0>();
223 
224 
225  // Now we initialize the grid with a filled circle. Outside the circle, the value of Phi_0 will be -1, inside +1.
226  //Now we construct the subsets based on the subset number.
227 
228  //We cast the global pointers to Particles and Particles_bulk as expected by the RHS functor.
229  PointerDistGlobal = (void *) &Particles;
230  auto verletList = Particles.template getVerlet<VL_NON_SYMMETRIC|VL_SKIP_REF_PART>(rCut);
231  //We create the DCPSE Based Laplacian operator.
232 
233  Laplacian<decltype(verletList)> Lap(Particles, verletList, 2, rCut, support_options::RADIUS);
234  auto C = getV<0>(Particles);
235  auto Init = getV<1>(Particles);
236  C=Init;
237  //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.
238  // 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)
239  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;
240  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;
241  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;
242  //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.
243  RHSFunctor<Laplacian<decltype(verletList)>, decltype(verletList)> System(Lap, verletList);
244 
245  //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.
246  ObserverFunctor ObserveAndUpdate;
247 
248 
249  //Furhter, odeint needs data in a state type "state_type_2d_ofp", we create one and fill in the initial condition.
251  //Since we created a 2d state_type we initialize the two fields in the object data using the method get.
252  X.data.get<x>() = C[0];
253  X.data.get<y>() = C[1];
254 
255 
256  std::vector<double> inter_times; // vector to store intermediate time steps taken by odeint.
257  Particles.deleteGhost();
258  Particles.write("Initial");
259  Particles.ghost_get<0>();
260  //for (int i = 0; i < 3; ++i) {
261  run_stepper_const<Odeint_rk4>(Particles, runtime_rk4_const,Lap, verletList);
262  // Particles.deleteGhost();
263  // Particles.write("RK4final");
264  // Particles.ghost_get<0>();
265  /*run_stepper_const<Odeint_rk5>(Particles, runtime_rk5_const,Lap,verletList);
266  Particles.deleteGhost();
267  Particles.write("RK5final");
268  Particles.ghost_get<0>();
269  run_stepper_const<Odeint_rk8>(Particles, runtime_rk78_const,Lap,verletList);
270  Particles.deleteGhost();
271  Particles.write("RK8");
272  Particles.ghost_get<0>();
273  run_stepper_adaptive<Odeint_rk5>(Particles, runtime_rk5_adapt,Lap,verletList);
274  Particles.deleteGhost();
275  Particles.write("AdapRK5");
276  Particles.ghost_get<0>();
277  */
278 
279 
280  //}
281 
282 // size_t steps = boost::numeric::odeint::integrate_const(Odeint_rk4, System, X, 0.0, tf, dt, ObserveAndUpdate);
283 // 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 );
284 
285  auto & vcl = create_vcluster();
286 
287  if (vcl.getProcessUnitID() == 0) {
288 
289  std::ofstream output_runtime;
290 
291  // head line
292  std::string headline = "cores,rk4,dopri5,fehlberg78,dopri5 adaptive";
293  if (access("runtime.csv", F_OK) == -1) {
294  output_runtime.open("runtime.csv");
295  output_runtime << headline << "\n";
296  } else {
297  output_runtime.open("runtime.csv", std::ios::app);
298  }
299 
300  output_runtime << vcl.getProcessingUnits()
301  << "," << average(runtime_rk4_const)
302  << "," << average(runtime_rk5_const)
303  << "," << average(runtime_rk78_const)
304  << "," << average(runtime_rk5_adapt)
305  << "\n";
306  }
307 
308  //Deallocating the operators
309  Lap.deallocate(Particles);
310  openfpm_finalize(); // Finalize openFPM library
311  return 0;
312 } //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 stop()
Stop the timer.
Definition: timer.hpp:119
void start()
Start the timer.
Definition: timer.hpp:90
double getwct()
Return the elapsed real time.
Definition: timer.hpp:130
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 deviceToHostProp()
Move the memory from the device to host memory.
void add()
Add local particle.
void hostToDeviceProp()
Move the memory from the device to host memory.
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.