OpenFPM_pdata  1.1.0
Project that contain the implementation of distributed structures
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Friends Pages
main.cpp
1 #include "Grid/grid_dist_id.hpp"
2 #include "data_type/aggregate.hpp"
3 #include "timer.hpp"
4 
5 #define FORTRAN_UPDATE
6 
7 #ifndef FORTRAN_UPDATE
8 #include "Vc/Vc"
9 #endif
10 
47 
49 
50 constexpr int x = 0;
51 constexpr int y = 1;
52 constexpr int z = 2;
53 
54 extern "C" void update_new(const int* lo, const int* hi,
55  double* u, const int* ulo, const int* uhi,
56  double* v, const int* vlo, const int* vhi,
57  double* flu, const int* fulo, const int* fuhi,
58  double* flv, const int* fvlo, const int* fvhi,
59  const double * dt, const double * uFactor, const double * vFactor, const double * F,
60  const double * K);
61 
62 
64 
65 void init(grid_dist_id<3,double,aggregate<double> > & OldU,
66  grid_dist_id<3,double,aggregate<double> > & OldV,
67  grid_dist_id<3,double,aggregate<double> > & NewU,
68  grid_dist_id<3,double,aggregate<double> > & NewV,
69  Box<3,double> & domain)
70 {
71  auto it = OldU.getDomainIterator();
72 
73  while (it.isNext())
74  {
75  // Get the local grid key
76  auto key = it.get();
77 
78  // Old values U and V
79  OldU.get(key) = 1.0;
80  OldV.get(key) = 0.0;
81 
82  // Old values U and V
83  NewU.get(key) = 0.0;
84  NewV.get(key) = 0.0;
85 
86  ++it;
87  }
88 
89  long int x_start = OldU.size(0)*1.55f/domain.getHigh(0);
90  long int y_start = OldU.size(1)*1.55f/domain.getHigh(1);
91  long int z_start = OldU.size(1)*1.55f/domain.getHigh(2);
92 
93  long int x_stop = OldU.size(0)*1.85f/domain.getHigh(0);
94  long int y_stop = OldU.size(1)*1.85f/domain.getHigh(1);
95  long int z_stop = OldU.size(1)*1.85f/domain.getHigh(2);
96 
97  grid_key_dx<3> start({x_start,y_start,z_start});
98  grid_key_dx<3> stop ({x_stop,y_stop,z_stop});
99  auto it_init = OldU.getSubDomainIterator(start,stop);
100 
101  while (it_init.isNext())
102  {
103  auto key = it_init.get();
104 
105  OldU.get(key) = 0.5 + (((double)std::rand())/RAND_MAX -0.5)/10.0;
106  OldV.get(key) = 0.25 + (((double)std::rand())/RAND_MAX -0.5)/20.0;
107 
108  ++it_init;
109  }
110 }
111 
112 
114 
115 void step(grid_dist_id<3, double, aggregate<double>> & OldU,
116  grid_dist_id<3, double, aggregate<double>> & OldV,
117  grid_dist_id<3, double, aggregate<double>> & NewU,
118  grid_dist_id<3, double, aggregate<double>> & NewV,
119  grid_key_dx<3> (& star_stencil_3D)[7],
120  double uFactor_s, double vFactor_s, double deltaT, double F, double K)
121 {
122 #ifndef FORTRAN_UPDATE
123 
125 
126  Vc::double uFactor = uFactor_s;
127  Vc::double vFactor = vFactor_s;
128 
129  WHILE_M(OldU,star_stencil_3D)
130  auto & U_old = GET_GRID_M(OldU);
131  auto & V_old = GET_GRID_M(OldV);
132 
133  auto & U_new = GET_GRID_M(NewU);
134  auto & V_new = GET_GRID_M(NewV);
135  ITERATE_3D_M(Vc::double_v::Size)
136 
137  // center point
138  auto Cp = it.getStencil<0>();
139 
140  // plus,minus X,Y,Z
141  auto mx = it.getStencil<1>();
142  auto px = it.getStencil<2>();
143  auto my = it.getStencil<3>();
144  auto py = it.getStencil<4>();
145  auto mz = it.getStencil<5>();
146  auto pz = it.getStencil<6>();
147 
148  //
149  Vc::double_v u_c(&U_old.get<0>(Cp),Vc::Unaligned);
150  Vc::double_v u_mz(&U_old.get<0>(mz),Vc::Unaligned);
151  Vc::double_v u_pz(&U_old.get<0>(pz),Vc::Unaligned);
152  Vc::double_v u_my(&U_old.get<0>(my),Vc::Unaligned);
153  Vc::double_v u_py(&U_old.get<0>(py),Vc::Unaligned);
154  Vc::double_v u_mx(&U_old.get<0>(mx),Vc::Unaligned);
155  Vc::double_v u_px(&U_old.get<0>(px),Vc::Unaligned);
156 
157 
158  Vc::double_v v_c(&V_old.get<0>(Cp),Vc::Unaligned);
159  Vc::double_v v_mz(&V_old.get<0>(mz),Vc::Unaligned);
160  Vc::double_v v_pz(&V_old.get<0>(pz),Vc::Unaligned);
161  Vc::double_v v_my(&V_old.get<0>(my),Vc::Unaligned);
162  Vc::double_v v_py(&V_old.get<0>(py),Vc::Unaligned);
163  Vc::double_v v_mx(&V_old.get<0>(mx),Vc::Unaligned);
164  Vc::double_v v_px(&V_old.get<0>(px),Vc::Unaligned);
165 
166  Vc::double_v out1 = u_c + uFactor * (u_mz + u_pz +
167  u_my + u_py +
168  u_mx + u_px +
169  - 6.0 * u_c) +
170  - deltaT * u_c * v_c * v_c
171  - deltaT * F * (u_c - 1.0);
172 
173  Vc::double_v out2 = v_c + vFactor * (v_mz + v_pz +
174  v_my + v_py +
175  v_mx + v_px +
176  - 6.0 * v_c ) +
177  deltaT * u_c * v_c * v_c +
178  - deltaT * (F+K) * v_c;
179 
180  out1.store(&U_new.get<0>(Cp),Vc::Unaligned);
181  out2.store(&V_new.get<0>(Cp),Vc::Unaligned);
182  END_LOOP_M(Vc::double_v::Size)
183 
185 
186 #else
187 
189 
190  auto & ginfo = OldU.getLocalGridsInfo();
191 
192  for (size_t i = 0 ; i < OldU.getN_loc_grid() ; i++)
193  {
194  auto & U_old = OldU.get_loc_grid(i);
195  auto & V_old = OldV.get_loc_grid(i);
196 
197  auto & U_new = NewU.get_loc_grid(i);
198  auto & V_new = NewV.get_loc_grid(i);
199 
200  int lo[3] = {(int)ginfo.get(i).Dbox.getLow(0),(int)ginfo.get(i).Dbox.getLow(1),(int)ginfo.get(i).Dbox.getLow(2)};
201  int hi[3] = {(int)ginfo.get(i).Dbox.getHigh(0),(int)ginfo.get(i).Dbox.getHigh(1),(int)ginfo.get(i).Dbox.getHigh(2)};
202 
203  int ulo[3] = {0,0,0};
204  int uhi[3] = {(int)ginfo.get(i).GDbox.getHigh(0),(int)ginfo.get(i).GDbox.getHigh(1),(int)ginfo.get(i).GDbox.getHigh(2)};
205  int nulo[3] = {0,0,0};
206  int nuhi[3] = {(int)ginfo.get(i).GDbox.getHigh(0),(int)ginfo.get(i).GDbox.getHigh(1),(int)ginfo.get(i).GDbox.getHigh(2)};
207 
208  int vlo[3] = {0,0,0};
209  int vhi[3] = {(int)ginfo.get(i).GDbox.getHigh(0),(int)ginfo.get(i).GDbox.getHigh(1),(int)ginfo.get(i).GDbox.getHigh(2)};
210  int nvlo[3] = {0,0,0};
211  int nvhi[3] = {(int)ginfo.get(i).GDbox.getHigh(0),(int)ginfo.get(i).GDbox.getHigh(1),(int)ginfo.get(i).GDbox.getHigh(2)};
212 
213  update_new(lo,hi,
214  (double *)U_old.getPointer(),ulo,uhi,
215  (double *)V_old.getPointer(),vlo,vhi,
216  (double *)U_new.getPointer(),nulo,nuhi,
217  (double *)V_new.getPointer(),nulo,nvhi,
218  &deltaT, &uFactor_s, &vFactor_s,&F,&K);
219  }
220 
222 
223 #endif
224 }
225 
227 
228 int main(int argc, char* argv[])
229 {
230  openfpm_init(&argc,&argv);
231 
232  // domain
233  Box<3,double> domain({0.0,0.0},{2.5,2.5,2.5});
234 
235  // grid size
236  size_t sz[3] = {256,256,256};
237 
238  // Define periodicity of the grid
239  periodicity<3> bc = {PERIODIC,PERIODIC,PERIODIC};
240 
241  // Ghost in grid unit
242  Ghost<3,long int> g(1);
243 
244  // deltaT
245  double deltaT = 1;
246 
247  // Diffusion constant for specie U
248  double du = 2*1e-5;
249 
250  // Diffusion constant for specie V
251  double dv = 1*1e-5;
252 
253  // Number of timesteps
254  size_t timeSteps = 5000;
255 
256  // K and F (Physical constant in the equation)
257  double K = 0.053;
258  double F = 0.014;
259 
261 
273 
275  grid_dist_id<3, double, aggregate<double>> OldU(sz,domain,g,bc);
276  grid_dist_id<3, double, aggregate<double>> OldV(OldU.getDecomposition(),sz,g);
277 
278  // New grid with the decomposition of the old grid
279  grid_dist_id<3, double, aggregate<double>> NewU(OldU.getDecomposition(),sz,g);
280  grid_dist_id<3, double, aggregate<double>> NewV(OldV.getDecomposition(),sz,g);
281 
282  // spacing of the grid on x and y
283 
284  double spacing[3] = {OldU.spacing(0),OldU.spacing(1),OldU.spacing(2)};
285 
286  init(OldU,OldV,NewU,NewV,domain);
287 
289 
290  // sync the ghost
291  size_t count = 0;
292  OldU.template ghost_get<0>();
293  OldV.template ghost_get<0>();
294 
295  // because we assume that spacing[x] == spacing[y] we use formula 2
296  // and we calculate the prefactor of Eq 2
297  double uFactor = deltaT * du/(spacing[x]*spacing[x]);
298  double vFactor = deltaT * dv/(spacing[x]*spacing[x]);
299 
300  timer tot_sim;
301  tot_sim.start();
302 
303  static grid_key_dx<3> star_stencil_3D[7] = {{0,0,0},
304  {0,0,-1},
305  {0,0,1},
306  {0,-1,0},
307  {0,1,0},
308  {-1,0,0},
309  {1,0,0}};
310 
311  for (size_t i = 0; i < timeSteps; ++i)
312  {
313  if (i % 300 == 0)
314  std::cout << "STEP: " << i << std::endl;
315 
347 
349  if (i % 2 == 0)
350  {
351  step(OldU,OldV,NewU,NewV,star_stencil_3D,uFactor,vFactor,deltaT,F,K);
352 
353  NewU.ghost_get<0>();
354  NewV.ghost_get<0>();
355  }
356  else
357  {
358  step(NewU,NewV,OldU,OldV,star_stencil_3D,uFactor,vFactor,deltaT,F,K);
359 
360  OldU.ghost_get<0>();
361  OldV.ghost_get<0>();
362  }
363 
365 
375 
377  // Every 2000 time step we output the configuration on hdf5
378  if (i % 2000 == 0)
379  {
380  OldU.save("output_u_" + std::to_string(count));
381  OldV.save("output_v_" + std::to_string(count));
382  count++;
383  }
384 
386  }
387 
388  tot_sim.stop();
389 
390  if (create_vcluster().rank() == 0)
391  {std::cout << "Total simulation: " << tot_sim.getwct() << std::endl;}
392 
393  // We frite the final configuration
394  OldV.write("final");
395 
397 
409 
411  openfpm_finalize();
412 
414 
423 }
424 
425 
grid_key_dx is the key to access any element in the grid
Definition: grid_key.hpp:18
T getHigh(int i) const
get the high interval of the box
Definition: Box.hpp:490
double getwct()
Return the elapsed real time.
Definition: timer.hpp:108
Definition: Ghost.hpp:39
mem_id get(size_t i) const
Get the i index.
Definition: grid_key.hpp:394
This is a distributed grid.
void start()
Start the timer.
Definition: timer.hpp:73
This class represent an N-dimensional box.
Definition: Box.hpp:56
This class is a trick to indicate the compiler a specific specialization pattern. ...
Definition: memory_c.hpp:201
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
Definition: aggregate.hpp:81
Class for cpu time benchmarking.
Definition: timer.hpp:25
St spacing(size_t i) const
Get the spacing of the grid in direction i.
void stop()
Stop the timer.
Definition: timer.hpp:97
[v_transform metafunction]