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 
22 
24 #include "Vector/vector_dist.hpp"
25 
27 
28 int main(int argc, char* argv[])
29 {
30 
47 
49  // initialize the library
50  openfpm_init(&argc,&argv);
51 
52  // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
53  Box<3,float> domain({0.0,0.0,0.0},{22.0,5.0,5.0});
54 
55  // Here we define the boundary conditions of our problem
56  size_t bc[3]={PERIODIC,PERIODIC,PERIODIC};
57 
58  // extended boundary around the domain, and the processor domain
59  Ghost<3,float> g(0.01);
60 
62 
89 
92 
94 
110 
112  auto it = vd.getDomainIterator();
113 
114  while (it.isNext())
115  {
116  auto key = it.get();
117 
118  // we define x, assign a random position between 0.0 and 1.0
119  vd.getPos(key)[0] = 22.0*((float)rand() / RAND_MAX);
120 
121  // we define y, assign a random position between 0.0 and 1.0
122  vd.getPos(key)[1] = 5.0*((float)rand() / RAND_MAX);
123 
124  // we define y, assign a random position between 0.0 and 1.0
125  vd.getPos(key)[1] = 5.0*((float)rand() / RAND_MAX);
126 
127  // next particle
128  ++it;
129  }
130 
132 
147 
149  vd.map();
150 
152 
167 
169  // Get a particle iterator
170  it = vd.getDomainIterator();
171 
172  // For each particle ...
173  while (it.isNext())
174  {
175  // ... p
176  auto p = it.get();
177 
178  // we set the properties of the particle p
179 
180  vd.template getProp<0>(p)[0] = 1.0;
181  vd.template getProp<0>(p)[1] = 1.0;
182  vd.template getProp<0>(p)[2] = 1.0;
183 
184  vd.template getProp<1>(p)[0] = 2.0;
185  vd.template getProp<1>(p)[1] = 2.0;
186  vd.template getProp<1>(p)[2] = 2.0;
187 
188  vd.template getProp<2>(p)[0] = 3.0;
189  vd.template getProp<2>(p)[1] = 3.0;
190  vd.template getProp<2>(p)[2] = 3.0;
191 
192  vd.template getProp<3>(p)[0] = 4.0;
193  vd.template getProp<3>(p)[1] = 4.0;
194  vd.template getProp<3>(p)[2] = 4.0;
195 
196  vd.template getProp<4>(p)[0] = 5.0;
197  vd.template getProp<4>(p)[1] = 5.0;
198  vd.template getProp<4>(p)[2] = 5.0;
199 
200 
201  // next particle
202  ++it;
203  }
204 
206 
230 
232  vd.save("particles_save.hdf5");
233 
235 
262 
265 
266  // load the particles on another vector
267  vd2.load("particles_save.hdf5");
268 
270 
287 
289  vector_dist<3,float, aggregate<float> > vd3(vd.getDecomposition(),0);
290 
291  auto it2 = vd2.getDomainIterator();
292 
293  while (it2.isNext())
294  {
295  auto p = it2.get();
296 
297  float magn_vort = sqrt(vd2.getProp<0>(p)[0]*vd2.getProp<0>(p)[0] +
298  vd2.getProp<0>(p)[1]*vd2.getProp<0>(p)[1] +
299  vd2.getProp<0>(p)[2]*vd2.getProp<0>(p)[2]);
300 
301  if (magn_vort < 0.1)
302  {
303  ++it2;
304  continue;
305  }
306 
307  vd3.add();
308 
309  vd3.getLastPos()[0] = vd2.getPos(p)[0];
310  vd3.getLastPos()[1] = vd2.getPos(p)[1];
311  vd3.getLastPos()[2] = vd2.getPos(p)[2];
312 
313  vd3.template getLastProp<0>() = magn_vort;
314 
315  ++it2;
316  }
317 
318  // We write the vtk file out from vd3
319  vd3.write("output", VTK_WRITER | FORMAT_BINARY );
320  vd3.save("particles_post_process_2");
321 
323 
335 
337  openfpm_finalize();
338 
340 
349 }
Definition: Ghost.hpp:39
bool write(std::string out, int opt=VTK_WRITER)
Output particle position and properties.
void add()
Add local particle.
This class is a trick to indicate the compiler a specific specialization pattern. ...
Definition: memory_c.hpp:201
vector_dist_iterator getDomainIterator() const
Get an iterator that traverse the particles in the domain.
Distributed vector.
void load(const std::string &filename)
Load the distributed vector from an HDF5 file.