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 
56 #include <stddef.h>
58 #include "Vector/vector_dist.hpp"
60 
61 int main(int argc, char* argv[])
62 {
63 
94 
96  // initialize the library
97  openfpm_init(&argc,&argv);
98 
99  // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
100  Box<2,float> domain({0.0,0.0},{1.0,1.0});
101 
102  // Here we define the boundary conditions of our problem
103  size_t bc[2]={PERIODIC,PERIODIC};
104 
105  // extended boundary around the domain, and the processor domain
106  Ghost<2,float> g(0.01);
107 
109 
148 
151 
152  // the scalar is the element at position 0 in the aggregate
153  const int scalar = 0;
154 
155  // the vector is the element at position 1 in the aggregate
156  const int vector = 1;
157 
158  // the tensor is the element at position 2 in the aggregate
159  const int tensor = 2;
160 
162 
185 
187  auto it = vd.getDomainIterator();
188 
189  while (it.isNext())
190  {
191  auto key = it.get();
192 
193  // we define x, assign a random position between 0.0 and 1.0
194  vd.getPos(key)[0] = (float)rand() / RAND_MAX;
195 
196  // we define y, assign a random position between 0.0 and 1.0
197  vd.getPos(key)[1] = (float)rand() / RAND_MAX;
198 
199  // next particle
200  ++it;
201  }
202 
204 
241 
243  vd.map();
244 
246 
260 
262  //Counter we use it later
263  size_t cnt = 0;
264 
265  // Get a particle iterator
266  it = vd.getDomainIterator();
267 
268  // For each particle ...
269  while (it.isNext())
270  {
271  // ... p
272  auto p = it.get();
273 
274  // we set the properties of the particle p
275 
276  // the scalar property
277  vd.template getProp<scalar>(p) = 1.0;
278 
279  vd.template getProp<vector>(p)[0] = 1.0;
280  vd.template getProp<vector>(p)[1] = 1.0;
281  vd.template getProp<vector>(p)[2] = 1.0;
282 
283  vd.template getProp<tensor>(p)[0][0] = 1.0;
284  vd.template getProp<tensor>(p)[0][1] = 1.0;
285  vd.template getProp<tensor>(p)[0][2] = 1.0;
286  vd.template getProp<tensor>(p)[1][0] = 1.0;
287  vd.template getProp<tensor>(p)[1][1] = 1.0;
288  vd.template getProp<tensor>(p)[1][2] = 1.0;
289  vd.template getProp<tensor>(p)[2][0] = 1.0;
290  vd.template getProp<tensor>(p)[2][1] = 1.0;
291  vd.template getProp<tensor>(p)[2][2] = 1.0;
292 
293  // increment the counter
294  cnt++;
295 
296  // next particle
297  ++it;
298  }
299 
301 
316 
318  auto & v_cl = create_vcluster();
319  v_cl.sum(cnt);
320  v_cl.execute();
321 
323 
351 
353  openfpm::vector<std::string> names({"scalar","vector","tensor"});
354  vd.setPropNames(names);
355 
356  // save vtk format (vtk is always the default)
357  vd.write("particles");
358 
359  // save in vtk binary format
360  vd.write("particles_bin",VTK_WRITER | FORMAT_BINARY);
361 
363 
375 
377  openfpm_finalize();
378 
380 
389 }
Definition: Ghost.hpp:39
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
Distributed vector.