OpenFPM  5.2.0
Project that contain the implementation of distributed structures
Vortex in Cell 3D (Optimization)

Vortex in Cell 3D ring with ringlets optimization

In this example we solve the Navier-Stokes equation in the vortex formulation in 3D for an incompressible fluid. This example has the following changes compared to Vortex in Cell 3D

  • Constructing grid is expensive in particular with a lot of cores. For this reason we create the grid in the main function rather than in the comp_vel function and helmotz_hodge_projection
  • Constructing also FDScheme is expensive so we construct it once in the main. We set the left hand side to the poisson operator, and inside the functions comp_vel and helmotz_hodge_projection just write the right hand side with impose_dit_b
grid_type g_vort(szu,domain,g,bc);
grid_type g_vel(g_vort.getDecomposition(),szu,g);
grid_type g_dvort(g_vort.getDecomposition(),szu,g);
particles_type particles(g_vort.getDecomposition(),0);
// Construct an FDScheme is heavy so we construct it here
// We also construct distributed temporal grids here because
// they are expensive
// Here we create temporal distributed grid, create a distributed grid is expensive we do it once outside
// And the vectorial phi_v
grid_type_s psi(g_vort.getDecomposition(),g_vort.getGridInfo().getSize(),g);
grid_type phi_v(g_vort.getDecomposition(),g_vort.getGridInfo().getSize(),g);
// In order to create a matrix that represent the poisson equation we have to indicate
// we have to indicate the maximum extension of the stencil and we we need an extra layer
// of points in case we have to impose particular boundary conditions.
Ghost<3,long int> stencil_max(2);
// We define a field phi_f
// We assemble the poisson equation doing the
// poisson of the Laplacian of phi using Laplacian
// central scheme (where the both derivative use central scheme, not
// forward composed backward like usually)
FDScheme<poisson_nn_helm> fd(stencil_max, domain, psi);
fd.template impose_dit<0>(poisson(),psi,psi.getDomainIterator());
Finite Differences.
Definition: FDScheme.hpp:127
Definition: eq.hpp:83
Definition: Ghost.hpp:40
Laplacian second order on h (spacing)
Definition: Laplacian.hpp:23
This is a distributed grid.
Distributed vector.
fd.new_b();
fd.template impose_dit_b<0>(psi,psi.getDomainIterator());

Another optimization that we do is to use a Hilbert space-filling curve as sub-sub-domain distribution strategy