OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
main.cpp
1 
23 #include "Vector/vector_dist.hpp"
26 
27 int main(int argc, char* argv[])
28 {
29 
48 
50  // initialize the library
51  openfpm_init(&argc,&argv);
52 
53  // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
54  Box<2,float> domain({0.0,0.0},{1.0,1.0});
55 
56  // Here we define the boundary conditions of our problem
57  size_t bc[2]={PERIODIC,PERIODIC};
58 
59  // extended boundary around the domain, and the processor domain
60  Ghost<2,float> g(0.01);
61 
63 
77 
80 
81  // the scalar is the element at position 0 in the aggregate
82  const int scalar = 0;
83 
84  // the vector is the element at position 1 in the aggregate
85  const int vector = 1;
86 
87  // the tensor is the element at position 2 in the aggregate
88  const int tensor = 2;
89 
91 
105 
107  auto it = vd.getDomainIterator();
108 
109  while (it.isNext())
110  {
111  // Particle p
112  auto p = it.get();
113 
114  // we define x, assign a random position between 0.0 and 1.0
115  vd.getPos(p)[0] = (float)rand() / RAND_MAX;
116 
117  // we define y, assign a random position between 0.0 and 1.0
118  vd.getPos(p)[1] = (float)rand() / RAND_MAX;
119 
120  // the scalar property
121  vd.template getProp<scalar>(p) = vd.getPos(p)[0];
122 
123  // The vector
124  vd.template getProp<vector>(p)[0] = 1.0;
125  vd.template getProp<vector>(p)[1] = 2.0;
126  vd.template getProp<vector>(p)[2] = 3.0;
127 
128  // The tensor
129  vd.template getProp<tensor>(p)[0][0] = 1.0;
130  vd.template getProp<tensor>(p)[0][1] = 2.0;
131  vd.template getProp<tensor>(p)[0][2] = 3.0;
132  vd.template getProp<tensor>(p)[1][0] = 4.0;
133  vd.template getProp<tensor>(p)[1][1] = 5.0;
134  vd.template getProp<tensor>(p)[1][2] = 6.0;
135  vd.template getProp<tensor>(p)[2][0] = 7.0;
136  vd.template getProp<tensor>(p)[2][1] = 8.0;
137  vd.template getProp<tensor>(p)[2][2] = 9.0;
138 
139  // next particle
140  ++it;
141  }
142 
144 
157 
159  vd.map();
160 
162 
220 
222  // write in vtk to see the result before
223  vd.write("before_ghost_get");
224 
225  // Here we synchronize the ghost get only the scalar and the tensor property
226  vd.ghost_get<scalar,tensor>();
227 
228  // write in vtk to see the result after
229  vd.write("after_ghost_get");
230 
232 
278 
280  vd.ghost_get<vector>(KEEP_PROPERTIES);
281 
282  // write in vtk to see the result before
283  vd.write("after_ghost_get_kp");
284 
286 
321 
323  it = vd.getGhostIterator();
324 
325  while (it.isNext())
326  {
327  // Particle p
328  auto p = it.get();
329 
330  // we shift down he particles
331  vd.getPos(p)[0] -= 1.0;
332 
333  // we shift
334  vd.getPos(p)[1] -= 1.0;
335 
336  // The vector
337  vd.template getProp<vector>(p)[0] = -1.0;
338  vd.template getProp<vector>(p)[1] = -2.0;
339  vd.template getProp<vector>(p)[2] = -3.0;
340 
341  // next particle
342  ++it;
343  }
344 
345  vd.ghost_get<vector>(KEEP_PROPERTIES | NO_POSITION);
346 
347  // write in vtk to see the result before
348  vd.write("after_ghost_get_kp_no_pos");
349 
350  vd.ghost_get<vector>(KEEP_PROPERTIES);
351 
352  vd.write("after_ghost_get_kp_no_pos_restore");
353 
355 
390 
392  vd.ghost_put<add_,vector>();
393 
394  // write in vtk to see the result before
395  vd.write("after_ghost_put");
396 
398 
410 
412  openfpm_finalize();
413 
415 
424 }
This structure define the operation add to use with copy general.
Definition: Ghost.hpp:39
This class represent an N-dimensional box.
Definition: Box.hpp:60
Distributed vector.