124 #include "interpolation/interpolation.hpp"
125 #include "Grid/grid_dist_id.hpp"
126 #include "Vector/vector_dist.hpp"
127 #include "Matrix/SparseMatrix.hpp"
128 #include "Vector/Vector.hpp"
129 #include "FiniteDifference/FDScheme.hpp"
130 #include "Solvers/petsc_solver.hpp"
131 #include "interpolation/mp4_kernel.hpp"
132 #include "Solvers/petsc_solver_AMG_report.hpp"
147 float sigma = 1.0/3.523;
150 float tgtre = 3000.0;
153 float nu = 1.0/tgtre;
160 const unsigned int nsteps = 10;
162 const unsigned int nsteps = 10001;
166 constexpr
unsigned int vorticity = 0;
167 constexpr
unsigned int velocity = 0;
168 constexpr
unsigned int p_vel = 1;
169 constexpr
int rhs_part = 2;
170 constexpr
unsigned int old_vort = 3;
171 constexpr
unsigned int old_pos = 4;
175 template<
typename gr
id>
void calc_and_print_max_div_and_int(
grid & g_vort)
179 g_vort.template ghost_get<vorticity>();
180 auto it5 = g_vort.getDomainIterator();
182 double max_vort = 0.0;
184 double int_vort[3] = {0.0,0.0,0.0};
188 auto key = it5.get();
192 tmp = 1.0f/g_vort.spacing(x)/2.0f*(g_vort.template get<vorticity>(key.move(x,1))[x] - g_vort.template get<vorticity>(key.move(x,-1))[x] ) +
193 1.0f/g_vort.spacing(y)/2.0f*(g_vort.template get<vorticity>(key.move(y,1))[y] - g_vort.template get<vorticity>(key.move(y,-1))[y] ) +
194 1.0f/g_vort.spacing(z)/2.0f*(g_vort.template get<vorticity>(key.move(z,1))[z] - g_vort.template get<vorticity>(key.move(z,-1))[z] );
196 int_vort[x] += g_vort.template get<vorticity>(key)[x];
197 int_vort[y] += g_vort.template get<vorticity>(key)[y];
198 int_vort[z] += g_vort.template get<vorticity>(key)[z];
208 v_cl.
sum(int_vort[0]);
209 v_cl.
sum(int_vort[1]);
210 v_cl.
sum(int_vort[2]);
214 {std::cout <<
"Max div for vorticity " << max_vort <<
" Integral: " << int_vort[0] <<
" " << int_vort[1] <<
" " << int_vort[2] << std::endl;}
253 constexpr
int nk = 32;
258 for (
size_t i = 0 ; i < nk ; i++)
260 ak[i] = rand()/RAND_MAX;
261 bk[i] = rand()/RAND_MAX;
265 float gamma = nu * tgtre;
266 float rinv2 = 1.0f/(sigma*sigma);
267 float max_vorticity = gamma*rinv2/M_PI;
274 auto key_d = it.get();
275 auto key = it.getGKey(key_d);
277 float tx = (key.get(x)-2)*gr.
spacing(x) + domain.getLow(x);
278 float ty = (key.get(y)-2)*gr.
spacing(y) + domain.getLow(y);
279 float tz = (key.get(z)-2)*gr.
spacing(z) + domain.getLow(z);
280 float theta1 = atan2((ty-domain.getHigh(1)/2.0f),(tz-domain.getHigh(2)/2.0f));
283 float rad1r = sqrt((ty-domain.getHigh(1)/2.0f)*(ty-domain.getHigh(1)/2.0f) + (tz-domain.getHigh(2)/2.0f)*(tz-domain.getHigh(2)/2.0f)) - ringr1;
284 float rad1t = tx - 1.0f;
285 float rad1sq = rad1r*rad1r + rad1t*rad1t;
286 float radstr = -exp(-rad1sq*rinv2)*rinv2*gamma/M_PI;
287 gr.template get<vorticity>(key_d)[x] = 0.0f;
288 gr.template get<vorticity>(key_d)[y] = -radstr * cos(theta1);
289 gr.template get<vorticity>(key_d)[z] = radstr * sin(theta1);
293 float rad1r_ = sqrt((ty-domain.getHigh(1)/2.0f)*(ty-domain.getHigh(1)/2.0f) + (tz-domain.getHigh(2)/2.0f)*(tz-domain.getHigh(2)/2.0f)) + ringr1;
294 float rad1sqTILDA = rad1sq*rinv2;
295 radstr = exp(-rad1sq*rinv2)*rinv2*gamma/M_PI;
296 gr.template get<vorticity>(key_d)[x] = 0.0f;
297 gr.template get<vorticity>(key_d)[y] = -radstr * cos(theta1);
298 gr.template get<vorticity>(key_d)[z] = radstr * sin(theta1);
317 static const unsigned int dims = 3;
319 static const unsigned int nvar = 1;
331 static const int grid_type = NORMAL_GRID;
460 constexpr
unsigned int phi = 0;
476 gr.template ghost_get<vorticity>();
491 psi.template get<phi>(key) = 1.0f/gr.
spacing(x)/2.0f*(gr.template get<vorticity>(key.move(x,1))[x] - gr.template get<vorticity>(key.move(x,-1))[x] ) +
492 1.0f/gr.
spacing(y)/2.0f*(gr.template get<vorticity>(key.move(y,1))[y] - gr.template get<vorticity>(key.move(y,-1))[y] ) +
493 1.0f/gr.
spacing(z)/2.0f*(gr.template get<vorticity>(key.move(z,1))[z] - gr.template get<vorticity>(key.move(z,-1))[z] );
501 calc_and_print_max_div_and_int(gr);
517 fd.template impose_dit<0>(poisson(),psi,psi.getDomainIterator());
551 solver.
solve(fd.getA(),x_,fd.getB());
558 solver.
solve(x_,fd.getB());
563 fd.template copy<phi>(x_,psi);
569 psi.template ghost_get<phi>();
577 auto key = it2.get();
579 gr.template get<vorticity>(key)[x] -= 1.0f/2.0f/psi.
spacing(x)*(psi.template get<phi>(key.move(x,1))-psi.template get<phi>(key.move(x,-1)));
580 gr.template get<vorticity>(key)[y] -= 1.0f/2.0f/psi.
spacing(y)*(psi.template get<phi>(key.move(y,1))-psi.template get<phi>(key.move(y,-1)));
581 gr.template get<vorticity>(key)[z] -= 1.0f/2.0f/psi.
spacing(z)*(psi.template get<phi>(key.move(z,1))-psi.template get<phi>(key.move(z,-1)));
588 calc_and_print_max_div_and_int(gr);
622 auto key = git.get_dist();
633 vd.template getLastProp<vorticity>()[x] = gr.template get<vorticity>(key)[x];
634 vd.template getLastProp<vorticity>()[y] = gr.template get<vorticity>(key)[y];
635 vd.template getLastProp<vorticity>()[z] = gr.template get<vorticity>(key)[z];
705 constexpr
unsigned int phi = 0;
727 calc_and_print_max_div_and_int(g_vort);
730 for (
size_t i = 0 ; i < 3 ; i++)
735 auto it2 = gr_ps.getDomainIterator();
740 auto key = it2.get();
743 gr_ps.get<vorticity>(key) = g_vort.template get<vorticity>(key)[i];
755 fd.template impose_dit<phi>(poisson(),gr_ps,gr_ps.getDomainIterator());
768 solver.
solve(phi_s[i],b);
779 {std::cout <<
"Solved component " << i <<
" Error: " << serr.
err_inf << std::endl;}
782 fd.template copy<phi>(phi_s[i],gr_ps);
788 auto it3 = gr_ps.getDomainIterator();
793 auto key = it3.get();
795 phi_v.get<velocity>(key)[i] = gr_ps.get<phi>(key);
805 phi_v.ghost_get<phi>();
807 float inv_dx = 0.5f / g_vort.
spacing(0);
808 float inv_dy = 0.5f / g_vort.
spacing(1);
809 float inv_dz = 0.5f / g_vort.
spacing(2);
811 auto it3 = phi_v.getDomainIterator();
816 auto key = it3.get();
818 float phi_xy = (phi_v.get<phi>(key.move(y,1))[x] - phi_v.get<phi>(key.move(y,-1))[x])*inv_dy;
819 float phi_xz = (phi_v.get<phi>(key.move(z,1))[x] - phi_v.get<phi>(key.move(z,-1))[x])*inv_dz;
820 float phi_yx = (phi_v.get<phi>(key.move(x,1))[y] - phi_v.get<phi>(key.move(x,-1))[y])*inv_dx;
821 float phi_yz = (phi_v.get<phi>(key.move(z,1))[y] - phi_v.get<phi>(key.move(z,-1))[y])*inv_dz;
822 float phi_zx = (phi_v.get<phi>(key.move(x,1))[z] - phi_v.get<phi>(key.move(x,-1))[z])*inv_dx;
823 float phi_zy = (phi_v.get<phi>(key.move(y,1))[z] - phi_v.get<phi>(key.move(y,-1))[z])*inv_dy;
825 g_vel.template get<velocity>(key)[x] = phi_zy - phi_yz;
826 g_vel.template get<velocity>(key)[y] = phi_xz - phi_zx;
827 g_vel.template get<velocity>(key)[z] = phi_yx - phi_xy;
845 gr.template get<prp>(key)[0] = 0.0;
846 gr.template get<prp>(key)[1] = 0.0;
847 gr.template get<prp>(key)[2] = 0.0;
863 vd.template getProp<prp>(key)[0] = 0.0;
864 vd.template getProp<prp>(key)[1] = 0.0;
865 vd.template getProp<prp>(key)[2] = 0.0;
887 template<
typename gr
id>
void calc_rhs(
grid & g_vort,
grid & g_vel,
grid & g_dwp)
890 constexpr
int rhs = 0;
894 float fac1 = 2.0f*nu/(g_vort.spacing(0)*g_vort.spacing(0));
895 float fac2 = 2.0f*nu/(g_vort.spacing(1)*g_vort.spacing(1));
896 float fac3 = 2.0f*nu/(g_vort.spacing(2)*g_vort.spacing(2));
898 float fac4 = 0.5f/(g_vort.spacing(0));
899 float fac5 = 0.5f/(g_vort.spacing(1));
900 float fac6 = 0.5f/(g_vort.spacing(2));
902 auto it = g_dwp.getDomainIterator();
904 g_vort.template ghost_get<vorticity>();
910 g_dwp.template get<rhs>(key)[x] = fac1*(g_vort.template get<vorticity>(key.move(x,1))[x]+g_vort.template get<vorticity>(key.move(x,-1))[x])+
911 fac2*(g_vort.template get<vorticity>(key.move(y,1))[x]+g_vort.template get<vorticity>(key.move(y,-1))[x])+
912 fac3*(g_vort.template get<vorticity>(key.move(z,1))[x]+g_vort.template get<vorticity>(key.move(z,-1))[x])-
913 2.0f*(fac1+fac2+fac3)*g_vort.template get<vorticity>(key)[x]+
914 fac4*g_vort.template get<vorticity>(key)[x]*
915 (g_vel.template get<velocity>(key.move(x,1))[x] - g_vel.template get<velocity>(key.move(x,-1))[x])+
916 fac5*g_vort.template get<vorticity>(key)[y]*
917 (g_vel.template get<velocity>(key.move(y,1))[x] - g_vel.template get<velocity>(key.move(y,-1))[x])+
918 fac6*g_vort.template get<vorticity>(key)[z]*
919 (g_vel.template get<velocity>(key.move(z,1))[x] - g_vel.template get<velocity>(key.move(z,-1))[x]);
921 g_dwp.template get<rhs>(key)[y] = fac1*(g_vort.template get<vorticity>(key.move(x,1))[y]+g_vort.template get<vorticity>(key.move(x,-1))[y])+
922 fac2*(g_vort.template get<vorticity>(key.move(y,1))[y]+g_vort.template get<vorticity>(key.move(y,-1))[y])+
923 fac3*(g_vort.template get<vorticity>(key.move(z,1))[y]+g_vort.template get<vorticity>(key.move(z,-1))[y])-
924 2.0f*(fac1+fac2+fac3)*g_vort.template get<vorticity>(key)[y]+
925 fac4*g_vort.template get<vorticity>(key)[x]*
926 (g_vel.template get<velocity>(key.move(x,1))[y] - g_vel.template get<velocity>(key.move(x,-1))[y])+
927 fac5*g_vort.template get<vorticity>(key)[y]*
928 (g_vel.template get<velocity>(key.move(y,1))[y] - g_vel.template get<velocity>(key.move(y,-1))[y])+
929 fac6*g_vort.template get<vorticity>(key)[z]*
930 (g_vel.template get<velocity>(key.move(z,1))[y] - g_vel.template get<velocity>(key.move(z,-1))[y]);
933 g_dwp.template get<rhs>(key)[z] = fac1*(g_vort.template get<vorticity>(key.move(x,1))[z]+g_vort.template get<vorticity>(key.move(x,-1))[z])+
934 fac2*(g_vort.template get<vorticity>(key.move(y,1))[z]+g_vort.template get<vorticity>(key.move(y,-1))[z])+
935 fac3*(g_vort.template get<vorticity>(key.move(z,1))[z]+g_vort.template get<vorticity>(key.move(z,-1))[z])-
936 2.0f*(fac1+fac2+fac3)*g_vort.template get<vorticity>(key)[z]+
937 fac4*g_vort.template get<vorticity>(key)[x]*
938 (g_vel.template get<velocity>(key.move(x,1))[z] - g_vel.template get<velocity>(key.move(x,-1))[z])+
939 fac5*g_vort.template get<vorticity>(key)[y]*
940 (g_vel.template get<velocity>(key.move(y,1))[z] - g_vel.template get<velocity>(key.move(y,-1))[z])+
941 fac6*g_vort.template get<vorticity>(key)[z]*
942 (g_vel.template get<velocity>(key.move(z,1))[z] - g_vel.template get<velocity>(key.move(z,-1))[z]);
994 auto key = it2.
get();
1039 auto key = it.
get();
1050 while (it2.isNext())
1052 auto key = it2.
get();
1088 template<
typename gr
id,
typename vector>
void do_step(vector &
particles,
1097 constexpr
int rhs = 0;
1102 inte.template p2m<vorticity,vorticity>(
particles,g_vort);
1104 g_vort.template ghost_put<add_,vorticity>();
1111 comp_vel(domain,g_vort,g_vel,phi_s,solver);
1112 calc_rhs(g_vort,g_vel,g_dvort);
1118 g_dvort.template ghost_get<rhs>();
1119 g_vel.template ghost_get<velocity>();
1122 inte.template m2p<rhs,rhs_part>(g_dvort,
particles);
1123 inte.template m2p<velocity,p_vel>(g_vel,
particles);
1128 template<
typename vector,
typename gr
id>
void check_point_and_save(vector &
particles,
1138 particles.
write(
"part_out_" + std::to_string(i),VTK_WRITER | FORMAT_BINARY);
1139 g_vort.template ghost_get<vorticity>();
1140 g_vel.template ghost_get<velocity>();
1141 g_vel.write_frame(
"grid_velocity",i, VTK_WRITER | FORMAT_BINARY);
1142 g_vort.write_frame(
"grid_vorticity",i, VTK_WRITER | FORMAT_BINARY);
1151 while (it_s.isNext())
1153 auto p = it_s.get();
1155 float vort_magn = sqrt(
particles.template getProp<vorticity>(p)[0] *
particles.template getProp<vorticity>(p)[0] +
1156 particles.template getProp<vorticity>(p)[1] *
particles.template getProp<vorticity>(p)[1] +
1157 particles.template getProp<vorticity>(p)[2] *
particles.template getProp<vorticity>(p)[2]);
1159 if (vort_magn < 0.1)
1171 part_save.template getLastProp<0>() = vort_magn;
1200 int main(
int argc,
char* argv[])
1203 openfpm_init(&argc,&argv);
1216 long int sz[] = {96,96,96};
1217 size_t szu[] = {(size_t)sz[0],(
size_t)sz[1],(size_t)sz[2]};
1222 grid_type g_vel(g_vort.getDecomposition(),szu,g);
1223 grid_type g_dvort(g_vort.getDecomposition(),szu,g);
1236 {std::cout <<
"SPACING " << g_vort.spacing(0) <<
" " << g_vort.spacing(1) <<
" " << g_vort.spacing(2) << std::endl;}
1239 init_ring(g_vort,domain);
1241 x_.
resize(g_vort.size(),g_vort.getLocalDomainSize());
1248 helmotz_hodge_projection(g_vort,domain,solver,x_,
true);
1259 phi_s[0].
resize(g_vort.size(),g_vort.getLocalDomainSize());
1260 phi_s[1].
resize(g_vort.size(),g_vort.getLocalDomainSize());
1261 phi_s[2].
resize(g_vort.size(),g_vort.getLocalDomainSize());
1277 {g_vort.write(
"grid_vorticity_init", VTK_WRITER | FORMAT_BINARY);}
1288 std::string restart(argv[1]);
1291 i = std::stoi(restart);
1299 {std::cout <<
"Restarting from " << i << std::endl;}
1303 for ( ; i < nsteps ; i++)
1306 do_step(
particles,g_vort,g_vel,g_dvort,domain,inte,phi_s,solver);
1312 do_step(
particles,g_vort,g_vel,g_dvort,domain,inte,phi_s,solver);
1319 inte.template p2m<vorticity,vorticity>(
particles,g_vort);
1320 g_vort.template ghost_put<add_,vorticity>();
1323 helmotz_hodge_projection(g_vort,domain,solver,x_,
false);
1329 {std::cout <<
"Step " << i << std::endl;}
1332 if (i % 100 == 0) {check_point_and_save(
particles,g_vort,g_vel,g_dvort,i);}
1342 int main(
int argc,
char* argv[])
Laplacian second order on h (spacing)
Sparse Matrix implementation.
void execute()
Execute all the requests.
void sum(T &num)
Sum the numbers across all processors and get the result.
size_t getProcessUnitID()
Get the process unit id.
size_t getProcessingUnits()
Get the total number of processors.
void max(T &num)
Get the maximum number across all processors (or reduction with infinity norm)
Implementation of VCluster class.
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support.
void resize(size_t row, size_t row_n)
stub resize
This is a distributed grid.
Decomposition & getDecomposition()
Get the object that store the information about the decomposition.
const grid_sm< dim, T > & getGridInfo() const
Get an object containing the grid informations.
grid_dist_iterator< dim, device_grid, decltype(device_grid::type_of_iterator()), FIXED > getDomainGhostIterator() const
It return an iterator that span the grid domain + ghost part.
St spacing(size_t i) const
Get the spacing of the grid in direction i.
grid_dist_iterator< dim, device_grid, decltype(device_grid::type_of_subiterator()), FREE > getDomainIterator() const
It return an iterator that span the full grid domain (each processor span its local domain)
__device__ __host__ const size_t(& getSize() const)[N]
Return the size of the grid as an array.
Main class for interpolation Particle to mest p2m and Mesh to particle m2p.
This class is able to do Matrix inversion in parallel with PETSC solvers.
Vector< T, PETSC_BASE > solve(SparseMatrix< T, int, PETSC_BASE > &A, const Vector< T, PETSC_BASE > &b, bool initial_guess=false)
Here we invert the matrix and solve the system.
void setPreconditionerAMG_cycleType(const std::string &cycle_type, int sweep_up=-1, int sweep_dw=-1, int sweep_crs=-1)
It set the type of cycle and optionally the number of sweep.
solError get_residual_error(SparseMatrix< T, int, PETSC_BASE > &A, const Vector< T, PETSC_BASE > &x, const Vector< T, PETSC_BASE > &b)
It return the resiual error.
void setPreconditionerAMG_nl(int nl)
Set the number of levels for the algebraic-multigrid preconditioner.
void setPreconditionerAMG_coarsenNodalType(int norm)
Set the block coarsening norm type.
void setPreconditioner(PCType type)
Set the preconditioner of the linear solver.
void setPreconditionerAMG_coarsen(const std::string &type)
Set the coarsening method for the algebraic-multi-grid preconditioner.
void setSolver(KSPType type)
Set the Petsc solver.
void setPreconditionerAMG_relax(const std::string &type, int k=REL_ALL)
Set the relaxation method for the algebraic-multi-grid preconditioner.
void setMaxIter(PetscInt n)
Set the maximum number of iteration for Krylov solvers.
void setAbsTol(PetscReal abstol_)
Set the absolute tolerance as stop criteria.
void setPreconditionerAMG_maxit(int nit)
Set the maximum number of V or W cycle for algebraic-multi-grid.
Class for cpu time benchmarking.
void stop()
Stop the timer.
void start()
Start the timer.
vect_dist_key_dx get()
Get the actual key.
size_t size_local() const
return the local size of the vector
auto getProp(vect_dist_key_dx vec_key) -> decltype(vPrp.template get< id >(vec_key.getKey()))
Get the property of an element.
Decomposition & getDecomposition()
Get the decomposition.
grid_dist_id_iterator_dec< Decomposition > getGridIterator(const size_t(&sz)[dim])
auto getPos(vect_dist_key_dx vec_key) -> decltype(vPos.template get< 0 >(vec_key.getKey()))
Get the position of an element.
vector_dist_iterator getDomainIterator() const
Get an iterator that traverse the particles in the domain.
void map(size_t opt=NONE)
It move all the particles that does not belong to the local processor to the respective processor.
void clear()
remove all the elements
bool write(std::string out, int opt=VTK_WRITER)
Output particle position and properties.
void save(const std::string &filename) const
Save the distributed vector on HDF5 file.
void add()
Add local particle.
void load(const std::string &filename)
Load the distributed vector from an HDF5 file.
auto getLastPos() -> decltype(vPos.template get< 0 >(0))
Get the position of the last element.
OutputIteratorT OffsetT ReductionOpT OuputT init
< [in] The initial value of the reduction
Vector< double, PETSC_BASE > Vector_type
Vector to solve the system (PETSC)
float stype
type of the spatial coordinates
SparseMatrix< double, int, PETSC_BASE > SparseMatrix_type
Sparse matrix used to sove the linear system (we use PETSC)
static const bool boundary[]
specify the boundary conditions
static const unsigned int nvar
We have only one scalar unknown.
static const unsigned int dims
Number of dimansions of the problem.
grid_dist_id< 3, float, aggregate< float > > b_grid
grid that store \psi
Meta-function to return a compatible zero-element.
It contain statistic of the error of the calculated solution.
PetscReal err_inf
infinity norm of the error