53#define PRINT_STACKTRACE
62#include "Vector/vector_dist.hpp"
64#include "Draw/DrawParticles.hpp"
68typedef float real_number;
77const real_number dp = 0.00425;
80real_number h_swl = 0.0;
83const real_number coeff_sound = 20.0;
86const real_number gamma_ = 7.0;
89const real_number H = 0.00736121593217;
92const real_number Eta2 = 0.01 * H*H;
94const real_number FourH2 = 4.0 * H*H;
97const real_number visco = 0.1;
100real_number cbar = 0.0;
103const real_number MassFluid = 0.0000767656;
106const real_number MassBound = 0.0000767656;
112const real_number t_end = 0.001;
114const real_number t_end = 1.5;
118const real_number gravity = 9.81;
121const real_number rho_zero = 1000.0;
127const real_number CFLnumber = 0.2;
130const real_number DtMin = 0.00001;
133const real_number RhoMin = 700.0;
136const real_number RhoMax = 1300.0;
139real_number max_fluid_height = 0.0;
144const size_t type = 0;
150const int rho_prev = 2;
153const int Pressure = 3;
162const int velocity = 6;
165const int velocity_prev = 7;
172typedef vector_dist_gpu<3,real_number,aggregate<unsigned int,real_number, real_number, real_number, real_number, real_number[3], real_number[3], real_number[3], real_number, real_number>>
particles;
181 template<
typename Decomposition,
typename vector>
inline void addComputation(
Decomposition & dec,
186 if (vd.template getProp<type>(p) == FLUID)
187 dec.addComputationCost(v,4);
189 dec.addComputationCost(v,3);
192 template<
typename Decomposition>
inline void applyModel(
Decomposition & dec,
size_t v)
194 dec.setSubSubDomainComputationCost(v, dec.getSubSubDomainComputationCost(v) * dec.getSubSubDomainComputationCost(v));
197 real_number distributionTol()
203template<
typename vd_type>
204__global__
void EqState_gpu(vd_type vd, real_number B)
206 auto a = GET_PARTICLE(vd);
208 real_number rho_a = vd.template getProp<rho>(a);
209 real_number rho_frac = rho_a / rho_zero;
211 vd.template getProp<Pressure>(a) = B*( rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac - 1.0);
216 auto it = vd.getDomainIteratorGPU();
218 CUDA_LAUNCH(EqState_gpu,it,vd.toKernel(),B);
222const real_number a2 = 1.0/M_PI/H/H/H;
224inline __device__ __host__ real_number Wab(real_number r)
229 return (1.0 - 3.0/2.0*r*r + 3.0/4.0*r*r*r)*a2;
231 return (1.0/4.0*(2.0 - r*r)*(2.0 - r*r)*(2.0 - r*r))*a2;
237const real_number c1 = -3.0/M_PI/H/H/H/H;
238const real_number d1 = 9.0/4.0/M_PI/H/H/H/H;
239const real_number c2 = -3.0/4.0/M_PI/H/H/H/H;
240const real_number a2_4 = 0.25*a2;
242real_number W_dap = 0.0;
246 const real_number qq=r/H;
248 real_number qq2 = qq * qq;
249 real_number fac1 = (c1*qq + d1*qq2)/r;
250 real_number b1 = (qq < 1.0f)?1.0f:0.0f;
252 real_number wqq = (2.0f - qq);
253 real_number fac2 = c2 * wqq * wqq / r;
254 real_number b2 = (qq >= 1.0f && qq < 2.0f)?1.0f:0.0f;
256 real_number factor = (b1*fac1 + b2*fac2);
258 DW.
get(0) = factor * dx.
get(0);
259 DW.
get(1) = factor * dx.
get(1);
260 DW.
get(2) = factor * dx.
get(2);
264inline __device__ __host__ real_number Tensile(real_number r, real_number rhoa, real_number rhob, real_number prs1, real_number prs2, real_number W_dap)
266 const real_number qq=r/H;
271 real_number wqq1=2.0f-qq;
272 real_number wqq2=wqq1*wqq1;
274 wab=a2_4*(wqq2*wqq1);
278 real_number wqq2=qq*qq;
279 real_number wqq3=wqq2*qq;
281 wab=a2*(1.0f-1.5f*wqq2+0.75f*wqq3);
285 real_number fab=wab*W_dap;
287 const real_number tensilp1=(prs1/(rhoa*rhoa))*(prs1>0.0f? 0.01f: -0.2f);
288 const real_number tensilp2=(prs2/(rhob*rhob))*(prs2>0.0f? 0.01f: -0.2f);
290 return (fab*(tensilp1+tensilp2));
294inline __device__ __host__ real_number Pi(
const Point<3,real_number> & dr, real_number rr2,
Point<3,real_number> & dv, real_number rhoa, real_number rhob, real_number massb, real_number cbar, real_number & visc)
296 const real_number dot = dr.
get(0)*dv.
get(0) + dr.
get(1)*dv.
get(1) + dr.
get(2)*dv.
get(2);
297 const real_number dot_rr2 = dot/(rr2+Eta2);
298 visc=(dot_rr2 < visc)?visc:dot_rr2;
302 const float amubar=H*dot_rr2;
303 const float robar=(rhoa+rhob)*0.5f;
304 const float pi_visc=(-visco*cbar*amubar/robar);
312template<
typename particles_type,
typename flu
id_
ids_type,
typename NN_type>
313__global__
void calc_forces_fluid_gpu(
particles_type vd, fluid_ids_type fids, NN_type NN, real_number W_dap, real_number cbar)
318 GET_PARTICLE_BY_ID(a,fids);
320 real_number max_visc = 0.0f;
326 unsigned int typea = vd.template getProp<type>(a);
329 real_number rhoa = vd.template getProp<rho>(a);
332 real_number Pa = vd.template getProp<Pressure>(a);
338 force_.
get(0) = 0.0f;
339 force_.
get(1) = 0.0f;
340 force_.
get(2) = -gravity;
341 real_number drho_ = 0.0f;
344 auto Np = NN.getNNIteratorBox(NN.getCell(xa));
347 while (Np.isNext() ==
true)
350 auto b = Np.get_sort();
358 unsigned int typeb = vd.template getProp<type>(b);
360 real_number massb = (typeb == FLUID)?MassFluid:MassBound;
362 real_number Pb = vd.template getProp<Pressure>(b);
363 real_number rhob = vd.template getProp<rho>(b);
369 real_number r2 = norm2(dr);
372 if (r2 < FourH2 && r2 >= 1e-16)
374 real_number r = sqrtf(r2);
379 real_number factor = - massb*((Pa + Pb) / (rhoa * rhob) + Tensile(r,rhoa,rhob,Pa,Pb,W_dap) + Pi(dr,r2,v_rel,rhoa,rhob,massb,cbar,max_visc));
382 factor = (typea == BOUNDARY && typeb == BOUNDARY)?0.0f:factor;
384 force_.
get(0) += factor * DW.
get(0);
385 force_.
get(1) += factor * DW.
get(1);
386 force_.
get(2) += factor * DW.
get(2);
388 real_number scal = massb*(v_rel.
get(0)*DW.
get(0)+v_rel.
get(1)*DW.
get(1)+v_rel.
get(2)*DW.
get(2));
389 scal = (typea == BOUNDARY && typeb == BOUNDARY)?0.0f:scal;
397 vd.template getProp<red>(a) = max_visc;
399 vd.template getProp<force>(a)[0] = force_.
get(0);
400 vd.template getProp<force>(a)[1] = force_.
get(1);
401 vd.template getProp<force>(a)[2] = force_.
get(2);
402 vd.template getProp<drho>(a) = drho_;
405template<
typename particles_type,
typename flu
id_
ids_type,
typename NN_type>
406__global__
void calc_forces_border_gpu(
particles_type vd, fluid_ids_type fbord, NN_type NN, real_number W_dap, real_number cbar)
411 GET_PARTICLE_BY_ID(a,fbord);
413 real_number max_visc = 0.0f;
419 unsigned int typea = vd.template getProp<type>(a);
424 real_number drho_ = 0.0f;
427 auto Np = NN.getNNIteratorBox(NN.getCell(xa));
430 while (Np.isNext() ==
true)
433 auto b = Np.get_sort();
441 unsigned int typeb = vd.template getProp<type>(b);
443 real_number massb = (typeb == FLUID)?MassFluid:MassBound;
450 real_number r2 = norm2(dr);
453 if (r2 < FourH2 && r2 >= 1e-16)
455 real_number r = sqrtf(r2);
460 real_number scal = massb*(v_rel.
get(0)*DW.
get(0)+v_rel.
get(1)*DW.
get(1)+v_rel.
get(2)*DW.
get(2));
461 scal = (typea == BOUNDARY && typeb == BOUNDARY)?0.0f:scal;
469 vd.template getProp<red>(a) = max_visc;
471 vd.template getProp<drho>(a) = drho_;
476 __device__
static bool check(
int c)
484 __device__
static bool check(
int c)
486 return c == BOUNDARY;
503 auto part = fluid_ids.getGPUIterator(96);
504 CUDA_LAUNCH(calc_forces_fluid_gpu,part,vd.toKernel_sorted(),fluid_ids.toKernel(),NN.toKernel(),W_dap,cbar);
506 part = border_ids.getGPUIterator(96);
507 CUDA_LAUNCH(calc_forces_border_gpu,part,vd.toKernel_sorted(),border_ids.toKernel(),NN.toKernel(),W_dap,cbar);
511 vd.merge_sort<force,drho,
red>(NN);
513 max_visc = reduce_local<red,_max_>(vd);
516template<
typename vector_type>
517__global__
void max_acceleration_and_velocity_gpu(
vector_type vd)
519 auto a = GET_PARTICLE(vd);
522 vd.template getProp<red>(a) = norm(acc);
525 vd.template getProp<red2>(a) = norm(vel);
528void max_acceleration_and_velocity(
particles & vd, real_number & max_acc, real_number & max_vel)
531 auto part = vd.getDomainIteratorGPU();
533 CUDA_LAUNCH(max_acceleration_and_velocity_gpu,part,vd.toKernel());
535 max_acc = reduce_local<red,_max_>(vd);
536 max_vel = reduce_local<red2,_max_>(vd);
545real_number calc_deltaT(
particles & vd, real_number ViscDtMax)
547 real_number Maxacc = 0.0;
548 real_number Maxvel = 0.0;
549 max_acceleration_and_velocity(vd,Maxacc,Maxvel);
552 const real_number dt_f = (Maxacc)?sqrt(H/Maxacc):std::numeric_limits<float>::max();
555 const real_number dt_cv = H/(std::max(cbar,Maxvel*10.f) + H*ViscDtMax);
558 real_number dt=real_number(CFLnumber)*std::min(dt_f,dt_cv);
559 if(dt<real_number(DtMin))
560 {dt=real_number(DtMin);}
565template<
typename vector_dist_type>
566__global__
void verlet_int_gpu(vector_dist_type vd, real_number dt, real_number dt2, real_number dt205)
569 auto a = GET_PARTICLE(vd);
572 if (vd.template getProp<type>(a) == BOUNDARY)
575 real_number rhop = vd.template getProp<rho>(a);
578 vd.template getProp<velocity>(a)[0] = 0.0;
579 vd.template getProp<velocity>(a)[1] = 0.0;
580 vd.template getProp<velocity>(a)[2] = 0.0;
581 real_number rhonew = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
582 vd.template getProp<rho>(a) = (rhonew < rho_zero)?rho_zero:rhonew;
584 vd.template getProp<rho_prev>(a) = rhop;
586 vd.template getProp<red>(a) = 0;
592 real_number dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
593 real_number dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
594 real_number dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
596 vd.getPos(a)[0] += dx;
597 vd.getPos(a)[1] += dy;
598 vd.getPos(a)[2] += dz;
600 real_number velX = vd.template getProp<velocity>(a)[0];
601 real_number velY = vd.template getProp<velocity>(a)[1];
602 real_number velZ = vd.template getProp<velocity>(a)[2];
604 real_number rhop = vd.template getProp<rho>(a);
606 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity_prev>(a)[0] + vd.template getProp<force>(a)[0]*dt2;
607 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity_prev>(a)[1] + vd.template getProp<force>(a)[1]*dt2;
608 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity_prev>(a)[2] + vd.template getProp<force>(a)[2]*dt2;
609 vd.template getProp<rho>(a) = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
612 if (vd.getPos(a)[0] < 0.0 || vd.getPos(a)[1] < 0.0 || vd.getPos(a)[2] < 0.0 ||
613 vd.getPos(a)[0] > 1.61 || vd.getPos(a)[1] > 0.68 || vd.getPos(a)[2] > 0.50 ||
614 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
616 vd.template getProp<red>(a) = 1;
620 vd.template getProp<red>(a) = 0;
624 vd.template getProp<velocity_prev>(a)[0] = velX;
625 vd.template getProp<velocity_prev>(a)[1] = velY;
626 vd.template getProp<velocity_prev>(a)[2] = velZ;
627 vd.template getProp<rho_prev>(a) = rhop;
632void verlet_int(
particles & vd, real_number dt)
635 auto part = vd.getDomainIteratorGPU();
637 real_number dt205 = dt*dt*0.5;
638 real_number dt2 = dt*2.0;
640 CUDA_LAUNCH(verlet_int_gpu,part,vd.toKernel(),dt,dt2,dt205);
643 remove_marked<red>(vd);
649template<
typename vector_type>
650__global__
void euler_int_gpu(
vector_type vd,real_number dt, real_number dt205)
653 auto a = GET_PARTICLE(vd);
656 if (vd.template getProp<type>(a) == BOUNDARY)
659 real_number rhop = vd.template getProp<rho>(a);
662 vd.template getProp<velocity>(a)[0] = 0.0;
663 vd.template getProp<velocity>(a)[1] = 0.0;
664 vd.template getProp<velocity>(a)[2] = 0.0;
665 real_number rhonew = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
666 vd.template getProp<rho>(a) = (rhonew < rho_zero)?rho_zero:rhonew;
668 vd.template getProp<rho_prev>(a) = rhop;
670 vd.template getProp<red>(a) = 0;
676 real_number dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
677 real_number dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
678 real_number dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
684 real_number velX = vd.template getProp<velocity>(a)[0];
685 real_number velY = vd.template getProp<velocity>(a)[1];
686 real_number velZ = vd.template getProp<velocity>(a)[2];
687 real_number rhop = vd.template getProp<rho>(a);
689 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity>(a)[0] + vd.template getProp<force>(a)[0]*dt;
690 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity>(a)[1] + vd.template getProp<force>(a)[1]*dt;
691 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity>(a)[2] + vd.template getProp<force>(a)[2]*dt;
692 vd.template getProp<rho>(a) = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
697 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
698 {vd.template getProp<red>(a) = 1;}
700 {vd.template getProp<red>(a) = 0;}
702 vd.template getProp<velocity_prev>(a)[0] = velX;
703 vd.template getProp<velocity_prev>(a)[1] = velY;
704 vd.template getProp<velocity_prev>(a)[2] = velZ;
705 vd.template getProp<rho_prev>(a) = rhop;
708void euler_int(
particles & vd, real_number dt)
712 auto part = vd.getDomainIteratorGPU();
714 real_number dt205 = dt*dt*0.5;
716 CUDA_LAUNCH(euler_int_gpu,part,vd.toKernel(),dt,dt205);
719 remove_marked<red>(vd);
724template<
typename vector_type,
typename NN_type>
727 real_number tot_ker = 0.0;
733 auto itg = NN.getNNIteratorBox(NN.getCell(xp));
736 auto q = itg.get_sort();
739 if (vd.template getProp<type>(q) != FLUID)
750 real_number r = sqrt(norm2(xp - xq));
752 real_number ker = Wab(r) * (MassFluid / rho_zero);
759 *press_tmp += vd.template getProp<Pressure>(q) * ker;
770 {*press_tmp = 1.0 / tot_ker * *press_tmp;}
773template<
typename Vector,
typename CellList>
774inline void sensor_pressure(
Vector & vd,
783 for (
size_t i = 0 ; i < probes.size() ; i++)
787 real_number press_tmp;
790 if (vd.getDecomposition().isLocal(probes.get(i)) ==
true)
792 CUDA_LAUNCH_DIM3(sensor_pressure_gpu,1,1,vd.toKernel_sorted(),NN.toKernel(),probes.
get(i),(real_number *)press_tmp_.toKernel());
794 vd.merge<Pressure>(NN);
797 press_tmp_.deviceToHost();
798 press_tmp = *(real_number *)press_tmp_.getPointer();
808 press_t.last().add(press_tmp);
812int main(
int argc,
char* argv[])
815 openfpm_init(&argc,&argv);
820#ifdef CUDIFY_USE_CUDA
821 cudaDeviceSetCacheConfig(cudaFuncCachePreferL1);
828 probes.add({0.8779,0.3,0.02});
829 probes.add({0.754,0.31,0.02});
833 size_t sz[3] = {413,179,133};
836 W_dap = 1.0/Wab(H/1.5);
839 size_t bc[3]={NON_PERIODIC,NON_PERIODIC,NON_PERIODIC};
844 particles vd(0,domain,bc,g,DEC_GRAN(128));
850 Box<3,real_number> fluid_box({dp/2.0,dp/2.0,dp/2.0},{0.4+dp/2.0,0.67-dp/2.0,0.3+dp/2.0});
856 max_fluid_height = fluid_it.getBoxMargins().
getHigh(2);
857 h_swl = fluid_it.getBoxMargins().
getHigh(2) - fluid_it.getBoxMargins().
getLow(2);
858 B = (coeff_sound)*(coeff_sound)*gravity*h_swl*rho_zero / gamma_;
859 cbar = coeff_sound * sqrt(gravity * h_swl);
862 while (fluid_it.isNext())
868 vd.getLastPos()[0] = fluid_it.get().get(0);
869 vd.getLastPos()[1] = fluid_it.get().get(1);
870 vd.getLastPos()[2] = fluid_it.get().get(2);
873 vd.template getLastProp<type>() = FLUID;
882 vd.template getLastProp<Pressure>() = rho_zero * gravity * (max_fluid_height - fluid_it.get().get(2));
884 vd.template getLastProp<rho>() = pow(vd.template getLastProp<Pressure>() / B + 1, 1.0/gamma_) * rho_zero;
885 vd.template getLastProp<rho_prev>() = vd.template getLastProp<rho>();
886 vd.template getLastProp<velocity>()[0] = 0.0;
887 vd.template getLastProp<velocity>()[1] = 0.0;
888 vd.template getLastProp<velocity>()[2] = 0.0;
890 vd.template getLastProp<velocity_prev>()[0] = 0.0;
891 vd.template getLastProp<velocity_prev>()[1] = 0.0;
892 vd.template getLastProp<velocity_prev>()[2] = 0.0;
903 Box<3,real_number> obstacle2({0.9+dp,0.24+dp/2.0,0.0},{1.02-dp/2.0,0.36-dp,0.45-dp/2.0});
907 holes.add(recipient2);
908 holes.add(obstacle1);
911 while (bound_box.isNext())
915 vd.getLastPos()[0] = bound_box.get().get(0);
916 vd.getLastPos()[1] = bound_box.get().get(1);
917 vd.getLastPos()[2] = bound_box.get().get(2);
919 vd.template getLastProp<type>() = BOUNDARY;
920 vd.template getLastProp<rho>() = rho_zero;
921 vd.template getLastProp<rho_prev>() = rho_zero;
922 vd.template getLastProp<velocity>()[0] = 0.0;
923 vd.template getLastProp<velocity>()[1] = 0.0;
924 vd.template getLastProp<velocity>()[2] = 0.0;
926 vd.template getLastProp<velocity_prev>()[0] = 0.0;
927 vd.template getLastProp<velocity_prev>()[1] = 0.0;
928 vd.template getLastProp<velocity_prev>()[2] = 0.0;
935 while (obstacle_box.isNext())
939 vd.getLastPos()[0] = obstacle_box.get().get(0);
940 vd.getLastPos()[1] = obstacle_box.get().get(1);
941 vd.getLastPos()[2] = obstacle_box.get().get(2);
943 vd.template getLastProp<type>() = BOUNDARY;
944 vd.template getLastProp<rho>() = rho_zero;
945 vd.template getLastProp<rho_prev>() = rho_zero;
946 vd.template getLastProp<velocity>()[0] = 0.0;
947 vd.template getLastProp<velocity>()[1] = 0.0;
948 vd.template getLastProp<velocity>()[2] = 0.0;
950 vd.template getLastProp<velocity_prev>()[0] = 0.0;
951 vd.template getLastProp<velocity_prev>()[1] = 0.0;
952 vd.template getLastProp<velocity_prev>()[2] = 0.0;
962 vd.addComputationCosts(md);
963 vd.getDecomposition().decompose();
970 vd.hostToDevicePos();
971 vd.template hostToDeviceProp<type,rho,rho_prev,Pressure,velocity>();
974 vd.ghost_get<type,rho,Pressure,velocity>(RUN_ON_DEVICE);
976 auto NN = vd.getCellListGPU(2*H / 2.0);
996 vd.map(RUN_ON_DEVICE);
999 vd.deviceToHostPos();
1000 vd.template deviceToHostProp<type>();
1004 vd.addComputationCosts(md);
1005 vd.getDecomposition().decompose();
1008 {std::cout <<
"REBALANCED " << it_reb << std::endl;}
1011 vd.map(RUN_ON_DEVICE);
1016 real_number max_visc = 0.0;
1018 vd.ghost_get<type,rho,Pressure,velocity>(RUN_ON_DEVICE);
1022 calc_forces(vd,NN,max_visc,cnt,fluid_ids,border_ids);
1029 real_number dt = calc_deltaT(vd,max_visc);
1047 vd.map(RUN_ON_DEVICE);
1048 vd.ghost_get<type,rho,Pressure,velocity>(RUN_ON_DEVICE);
1049 vd.updateCellList(NN);
1054 std::cout <<
"OUTPUT " << dt << std::endl;
1058 vd.deviceToHostPos();
1059 vd.deviceToHostProp<type,rho,rho_prev,Pressure,drho,force,velocity,velocity_prev,
red,red2>();
1064 auto ito = vd.getDomainIterator();
1072 vd_out.getLastPos()[0] = vd.getPos(p)[0];
1073 vd_out.getLastPos()[1] = vd.getPos(p)[1];
1074 vd_out.getLastPos()[2] = vd.getPos(p)[2];
1076 vd_out.template getLastProp<0>() = vd.template getProp<type>(p);
1078 vd_out.template getLastProp<1>()[0] = vd.template getProp<velocity>(p)[0];
1079 vd_out.template getLastProp<1>()[1] = vd.template getProp<velocity>(p)[1];
1080 vd_out.template getLastProp<1>()[2] = vd.template getProp<velocity>(p)[2];
1085 vd_out.write_frame(
"Particles",write,VTK_WRITER | FORMAT_BINARY);
1089 {std::cout <<
"TIME: " << t <<
" write " << it_time.
getwct() <<
" " << it_reb <<
" " << cnt <<
" Max visc: " << max_visc <<
" " << vd.size_local() << std::endl;}
1094 {std::cout <<
"TIME: " << t <<
" " << it_time.
getwct() <<
" " << it_reb <<
" " << cnt <<
" Max visc: " << max_visc <<
" " << vd.size_local() << std::endl;}
1099 std::cout <<
"Time to complete: " << tot_sim.
getwct() <<
" seconds" << std::endl;
1107int main(
int argc,
char* argv[])
This class represent an N-dimensional box.
__device__ __host__ T getLow(int i) const
get the i-coordinate of the low bound interval of the box
__device__ __host__ T getHigh(int i) const
get the high interval of the box
Class for FAST cell list implementation.
auto get(size_t cell, size_t ele) -> decltype(this->Mem_type::get(cell, ele))
Get an element in the cell.
This class define the domain decomposition interface.
static PointIteratorSkin< dim, T, typename vd_type::Decomposition_type > DrawSkin(vd_type &vd, size_t(&sz)[dim], Box< dim, T > &domain, Box< dim, T > &sub_A, Box< dim, T > &sub_B)
Draw particles in a box B excluding the area of a second box A (B - A)
static PointIterator< dim, T, typename vd_type::Decomposition_type > DrawBox(vd_type &vd, size_t(&sz)[dim], Box< dim, T > &domain, Box< dim, T > &sub)
Draw particles in a box.
This class implement the point shape in an N-dimensional space.
__device__ __host__ const T & get(unsigned int i) const
Get coordinate.
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.
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.
Implementation of 1-D std::vector like structure.
Class for cpu time benchmarking.
void stop()
Stop the timer.
void start()
Start the timer.
double getwct()
Return the elapsed real time.
size_t size_local() const
return the local size of the vector
void updateCellList(CellL &cell_list, bool no_se3=false, cl_construct_opt opt=cl_construct_opt::Full)
Update a cell list using the stored particles.
auto getPos(vect_dist_key_dx vec_key) -> decltype(v_pos.template get< 0 >(vec_key.getKey()))
Get the position of an element.
const vector_dist_prop & getPropVectorSort() const
return the property vector of all the particles
Vcluster< Memory > & getVC()
Get the Virtual Cluster machine.
Model for Dynamic load balancing.
aggregate of properties, from a list of object if create a struct that follow the OPENFPM native stru...
temporal buffer for reductions