OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
main.cpp
1 
63 #include <stddef.h>
65 #include "Vector/vector_dist.hpp"
67 
68 
69 int main(int argc, char* argv[])
70 {
101 
103  // initialize the library
104  openfpm_init(&argc,&argv);
105 
106  // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
107  Box<2,float> domain({0.0,0.0},{1.0,1.0});
108 
109  // Here we define the boundary conditions of our problem
110  size_t bc[2]={PERIODIC,PERIODIC};
111 
112  // extended boundary around the domain, and the processor domain
113  Ghost<2,float> g(0.01);
114 
116 
155 
158 
159  // the scalar is the element at position 0 in the aggregate
160  const int scalar = 0;
161 
162  // the vector is the element at position 1 in the aggregate
163  const int vector = 1;
164 
165  // the tensor is the element at position 2 in the aggregate
166  const int tensor = 2;
167 
169 
192 
194  auto it = vd.getDomainIterator();
195 
196  while (it.isNext())
197  {
198  auto key = it.get();
199 
200  // we define x, assign a random position between 0.0 and 1.0
201  vd.getPos(key)[0] = (float)rand() / RAND_MAX;
202 
203  // we define y, assign a random position between 0.0 and 1.0
204  vd.getPos(key)[1] = (float)rand() / RAND_MAX;
205 
206  // next particle
207  ++it;
208  }
209 
211 
248 
250  vd.map();
251 
253 
267 
269  //Counter we use it later
270  size_t cnt = 0;
271 
272  // Get a particle iterator
273  it = vd.getDomainIterator();
274 
275  // For each particle ...
276  while (it.isNext())
277  {
278  // ... p
279  auto p = it.get();
280 
281  // we set the properties of the particle p
282 
283  // the scalar property
284  vd.template getProp<scalar>(p) = 1.0;
285 
286  vd.template getProp<vector>(p)[0] = 1.0;
287  vd.template getProp<vector>(p)[1] = 1.0;
288  vd.template getProp<vector>(p)[2] = 1.0;
289 
290  vd.template getProp<tensor>(p)[0][0] = 1.0;
291  vd.template getProp<tensor>(p)[0][1] = 1.0;
292  vd.template getProp<tensor>(p)[0][2] = 1.0;
293  vd.template getProp<tensor>(p)[1][0] = 1.0;
294  vd.template getProp<tensor>(p)[1][1] = 1.0;
295  vd.template getProp<tensor>(p)[1][2] = 1.0;
296  vd.template getProp<tensor>(p)[2][0] = 1.0;
297  vd.template getProp<tensor>(p)[2][1] = 1.0;
298  vd.template getProp<tensor>(p)[2][2] = 1.0;
299 
300  // increment the counter
301  cnt++;
302 
303  // next particle
304  ++it;
305  }
306 
308 
323 
325  auto & v_cl = create_vcluster();
326  v_cl.sum(cnt);
327  v_cl.execute();
328 
330 
358 
360  openfpm::vector<std::string> names({"scalar","vector","tensor"});
361  vd.setPropNames(names);
362 
363  // save vtk format (vtk is always the default)
364  vd.write("particles");
365 
366  // save in vtk binary format
367  vd.write("particles_bin",VTK_WRITER | FORMAT_BINARY);
368 
369  // save in vtk format with time
370  vd.write("particles_with_time","time=1.234");
371 
372  // save in vtk format with time
373  vd.write("particles_with_time_bin","time=1.234",VTK_WRITER | FORMAT_BINARY);
374 
376 
388 
390  openfpm_finalize();
391 
393 
402 }
Definition: Ghost.hpp:39
This class represent an N-dimensional box.
Definition: Box.hpp:60
Distributed vector.