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 
24 #include "Vector/vector_dist.hpp"
27 
28 int main(int argc, char* argv[])
29 {
30 
49 
51  // initialize the library
52  openfpm_init(&argc,&argv);
53 
54  // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
55  Box<2,float> domain({0.0,0.0},{1.0,1.0});
56 
57  // Here we define the boundary conditions of our problem
58  size_t bc[2]={PERIODIC,PERIODIC};
59 
60  // extended boundary around the domain, and the processor domain
61  Ghost<2,float> g(0.01);
62 
64 
78 
81 
82  // the scalar is the element at position 0 in the aggregate
83  const int scalar = 0;
84 
85  // the vector is the element at position 1 in the aggregate
86  const int vector = 1;
87 
88  // the tensor is the element at position 2 in the aggregate
89  const int tensor = 2;
90 
92 
106 
108  auto it = vd.getDomainIterator();
109 
110  while (it.isNext())
111  {
112  // Particle p
113  auto p = it.get();
114 
115  // we define x, assign a random position between 0.0 and 1.0
116  vd.getPos(p)[0] = (float)rand() / RAND_MAX;
117 
118  // we define y, assign a random position between 0.0 and 1.0
119  vd.getPos(p)[1] = (float)rand() / RAND_MAX;
120 
121  // the scalar property
122  vd.template getProp<scalar>(p) = vd.getPos(p)[0];
123 
124  // The vector
125  vd.template getProp<vector>(p)[0] = 1.0;
126  vd.template getProp<vector>(p)[1] = 2.0;
127  vd.template getProp<vector>(p)[2] = 3.0;
128 
129  // The tensor
130  vd.template getProp<tensor>(p)[0][0] = 1.0;
131  vd.template getProp<tensor>(p)[0][1] = 2.0;
132  vd.template getProp<tensor>(p)[0][2] = 3.0;
133  vd.template getProp<tensor>(p)[1][0] = 4.0;
134  vd.template getProp<tensor>(p)[1][1] = 5.0;
135  vd.template getProp<tensor>(p)[1][2] = 6.0;
136  vd.template getProp<tensor>(p)[2][0] = 7.0;
137  vd.template getProp<tensor>(p)[2][1] = 8.0;
138  vd.template getProp<tensor>(p)[2][2] = 9.0;
139 
140  // next particle
141  ++it;
142  }
143 
145 
158 
160  vd.map();
161 
163 
193 
195  // write in vtk to see the result before
196  vd.write("before_ghost_get");
197 
198  // Here we synchronize the ghost get only the scalar and the tensor property
199  vd.ghost_get<scalar,tensor>();
200 
201  // write in vtk to see the result after
202  vd.write("after_ghost_get");
203 
205 
251 
253  vd.ghost_get<vector>(KEEP_PROPERTIES);
254 
255  // write in vtk to see the result before
256  vd.write("after_ghost_get_kp");
257 
259 
294 
296  it = vd.getGhostIterator();
297 
298  while (it.isNext())
299  {
300  // Particle p
301  auto p = it.get();
302 
303  // we shift down he particles
304  vd.getPos(p)[0] -= 1.0;
305 
306  // we shift
307  vd.getPos(p)[1] -= 1.0;
308 
309  // The vector
310  vd.template getProp<vector>(p)[0] = -1.0;
311  vd.template getProp<vector>(p)[1] = -2.0;
312  vd.template getProp<vector>(p)[2] = -3.0;
313 
314  // next particle
315  ++it;
316  }
317 
318  vd.ghost_get<vector>(KEEP_PROPERTIES | NO_POSITION);
319 
320  // write in vtk to see the result before
321  vd.write("after_ghost_get_kp_no_pos");
322 
323  vd.ghost_get<vector>(KEEP_PROPERTIES);
324 
325  vd.write("after_ghost_get_kp_no_pos_restore");
326 
328 
363 
365  vd.ghost_put<add_,vector>();
366 
367  // write in vtk to see the result before
368  vd.write("after_ghost_put");
369 
371 
383 
385  openfpm_finalize();
386 
388 
397 }
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:56
Distributed vector.