53 #define PRINT_STACKTRACE
58 #include "Vector/vector_dist.hpp"
60 #include "Draw/DrawParticles.hpp"
64 typedef float real_number;
73 const real_number dp = 0.0085;
76 real_number h_swl = 0.0;
79 const real_number coeff_sound = 20.0;
82 const real_number gamma_ = 7.0;
85 const real_number H = 0.0147224318643;
88 const real_number Eta2 = 0.01 * H*H;
90 const real_number FourH2 = 4.0 * H*H;
93 const real_number visco = 0.1;
96 real_number cbar = 0.0;
99 const real_number MassFluid = 0.000614125;
102 const real_number MassBound = 0.000614125;
108 const real_number t_end = 0.001;
110 const real_number t_end = 1.5;
114 const real_number gravity = 9.81;
117 const real_number rho_zero = 1000.0;
123 const real_number CFLnumber = 0.2;
126 const real_number DtMin = 0.00001;
129 const real_number RhoMin = 700.0;
132 const real_number RhoMax = 1300.0;
135 real_number max_fluid_height = 0.0;
140 const size_t type = 0;
146 const int rho_prev = 2;
149 const int Pressure = 3;
158 const int velocity = 6;
161 const int velocity_prev = 7;
168 typedef 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;
177 template<
typename Decomposition,
typename vector>
inline void addComputation(
183 if (vd.template getProp<type>(p) == FLUID)
184 dec.addComputationCost(v,4);
186 dec.addComputationCost(v,3);
189 template<
typename Decomposition>
inline void applyModel(
Decomposition & dec,
size_t v)
191 dec.setSubSubDomainComputationCost(v, dec.getSubSubDomainComputationCost(v) * dec.getSubSubDomainComputationCost(v));
194 real_number distributionTol()
200 template<
typename vd_type>
201 __global__
void EqState_gpu(vd_type vd, real_number B)
203 auto a = GET_PARTICLE(vd);
205 real_number rho_a = vd.template getProp<rho>(a);
206 real_number rho_frac = rho_a / rho_zero;
208 vd.template getProp<Pressure>(a) = B*( rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac - 1.0);
213 auto it = vd.getDomainIteratorGPU();
215 CUDA_LAUNCH(EqState_gpu,it,vd.toKernel(),B);
219 const real_number a2 = 1.0/M_PI/H/H/H;
221 inline __device__ __host__ real_number Wab(real_number r)
226 return (1.0 - 3.0/2.0*r*r + 3.0/4.0*r*r*r)*a2;
228 return (1.0/4.0*(2.0 - r)*(2.0 - r)*(2.0 - r))*a2;
234 const real_number c1 = -3.0/M_PI/H/H/H/H;
235 const real_number d1 = 9.0/4.0/M_PI/H/H/H/H;
236 const real_number c2 = -3.0/4.0/M_PI/H/H/H/H;
237 const real_number a2_4 = 0.25*a2;
239 real_number W_dap = 0.0;
243 const real_number qq=r/H;
245 real_number qq2 = qq * qq;
246 real_number fac1 = (c1*qq + d1*qq2)/r;
247 real_number b1 = (qq < 1.0f)?1.0f:0.0f;
249 real_number wqq = (2.0f - qq);
250 real_number fac2 = c2 * wqq * wqq / r;
251 real_number b2 = (qq >= 1.0f && qq < 2.0f)?1.0f:0.0f;
253 real_number factor = (b1*fac1 + b2*fac2);
255 DW.
get(0) = factor * dx.
get(0);
256 DW.
get(1) = factor * dx.
get(1);
257 DW.
get(2) = factor * dx.
get(2);
261 inline __device__ __host__ real_number Tensile(real_number r, real_number rhoa, real_number rhob, real_number prs1, real_number prs2, real_number W_dap)
263 const real_number qq=r/H;
268 real_number wqq1=2.0f-qq;
269 real_number wqq2=wqq1*wqq1;
271 wab=a2_4*(wqq2*wqq1);
275 real_number wqq2=qq*qq;
276 real_number wqq3=wqq2*qq;
278 wab=a2*(1.0f-1.5f*wqq2+0.75f*wqq3);
282 real_number fab=wab*W_dap;
284 const real_number tensilp1=(prs1/(rhoa*rhoa))*(prs1>0.0f? 0.01f: -0.2f);
285 const real_number tensilp2=(prs2/(rhob*rhob))*(prs2>0.0f? 0.01f: -0.2f);
287 return (fab*(tensilp1+tensilp2));
291 inline __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)
293 const real_number dot = dr.
get(0)*dv.
get(0) + dr.
get(1)*dv.
get(1) + dr.
get(2)*dv.
get(2);
294 const real_number dot_rr2 = dot/(rr2+Eta2);
295 visc=(dot_rr2 < visc)?visc:dot_rr2;
299 const float amubar=H*dot_rr2;
300 const float robar=(rhoa+rhob)*0.5f;
301 const float pi_visc=(-visco*cbar*amubar/robar);
309 template<
typename particles_type,
typename flu
id_
ids_type,
typename NN_type>
310 __global__
void calc_forces_fluid_gpu(
particles_type vd, fluid_ids_type fids, NN_type NN, real_number W_dap, real_number cbar)
315 GET_PARTICLE_BY_ID(a,fids);
317 real_number max_visc = 0.0f;
323 unsigned int typea = vd.template getProp<type>(a);
326 real_number rhoa = vd.template getProp<rho>(a);
329 real_number Pa = vd.template getProp<Pressure>(a);
335 force_.
get(0) = 0.0f;
336 force_.
get(1) = 0.0f;
337 force_.
get(2) = -gravity;
338 real_number drho_ = 0.0f;
341 auto Np = NN.getNNIteratorBox(NN.getCell(xa));
344 while (Np.isNext() ==
true)
347 auto b = Np.get_sort();
355 unsigned int typeb = vd.template getProp<type>(b);
357 real_number massb = (typeb == FLUID)?MassFluid:MassBound;
359 real_number Pb = vd.template getProp<Pressure>(b);
360 real_number rhob = vd.template getProp<rho>(b);
365 real_number r2 = norm2(dr);
368 if (r2 < FourH2 && r2 >= 1e-16)
370 real_number r = sqrtf(r2);
377 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));
380 factor = (typea == BOUNDARY && typeb == BOUNDARY)?0.0f:factor;
382 force_.
get(0) += factor * DW.
get(0);
383 force_.
get(1) += factor * DW.
get(1);
384 force_.
get(2) += factor * DW.
get(2);
386 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));
387 scal = (typea == BOUNDARY && typeb == BOUNDARY)?0.0f:scal;
395 vd.template getProp<red>(a) = max_visc;
397 vd.template getProp<force>(a)[0] = force_.
get(0);
398 vd.template getProp<force>(a)[1] = force_.
get(1);
399 vd.template getProp<force>(a)[2] = force_.
get(2);
400 vd.template getProp<drho>(a) = drho_;
403 template<
typename particles_type,
typename flu
id_
ids_type,
typename NN_type>
404 __global__
void calc_forces_border_gpu(
particles_type vd, fluid_ids_type fbord, NN_type NN, real_number W_dap, real_number cbar)
409 GET_PARTICLE_BY_ID(a,fbord);
411 real_number max_visc = 0.0f;
417 unsigned int typea = vd.template getProp<type>(a);
422 real_number drho_ = 0.0f;
425 auto Np = NN.getNNIteratorBox(NN.getCell(xa));
428 while (Np.isNext() ==
true)
431 auto b = Np.get_sort();
439 unsigned int typeb = vd.template getProp<type>(b);
441 real_number massb = (typeb == FLUID)?MassFluid:MassBound;
448 real_number r2 = norm2(dr);
451 if (r2 < FourH2 && r2 >= 1e-16)
453 real_number r = sqrtf(r2);
458 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));
459 scal = (typea == BOUNDARY && typeb == BOUNDARY)?0.0f:scal;
467 vd.template getProp<red>(a) = max_visc;
469 vd.template getProp<drho>(a) = drho_;
474 __device__
static bool check(
int c)
480 struct type_is_border
482 __device__
static bool check(
int c)
484 return c == BOUNDARY;
491 vd.template updateCellListGPU<type,rho,Pressure,velocity>(NN);
501 auto part = fluid_ids.getGPUIterator(96);
502 CUDA_LAUNCH(calc_forces_fluid_gpu,part,vd.toKernel(),fluid_ids.toKernel(),NN.toKernel(),W_dap,cbar);
504 part = border_ids.getGPUIterator(96);
505 CUDA_LAUNCH(calc_forces_border_gpu,part,vd.toKernel(),border_ids.toKernel(),NN.toKernel(),W_dap,cbar);
509 vd.template restoreOrder<drho,force,red>(NN);
511 max_visc = reduce_local<red,_max_>(vd);
514 template<
typename vector_type>
515 __global__
void max_acceleration_and_velocity_gpu(
vector_type vd)
517 auto a = GET_PARTICLE(vd);
520 vd.template getProp<red>(a) = norm(acc);
523 vd.template getProp<red2>(a) = norm(vel);
526 void max_acceleration_and_velocity(
particles & vd, real_number & max_acc, real_number & max_vel)
529 auto part = vd.getDomainIteratorGPU();
531 CUDA_LAUNCH(max_acceleration_and_velocity_gpu,part,vd.toKernel());
533 max_acc = reduce_local<red,_max_>(vd);
534 max_vel = reduce_local<red2,_max_>(vd);
543 real_number calc_deltaT(
particles & vd, real_number ViscDtMax)
545 real_number Maxacc = 0.0;
546 real_number Maxvel = 0.0;
547 max_acceleration_and_velocity(vd,Maxacc,Maxvel);
550 const real_number dt_f = (Maxacc)?sqrt(H/Maxacc):std::numeric_limits<float>::max();
553 const real_number dt_cv = H/(std::max(cbar,Maxvel*10.f) + H*ViscDtMax);
556 real_number dt=real_number(CFLnumber)*std::min(dt_f,dt_cv);
557 if(dt<real_number(DtMin))
558 {dt=real_number(DtMin);}
563 template<
typename vector_dist_type>
564 __global__
void verlet_int_gpu(
vector_dist_type vd, real_number dt, real_number dt2, real_number dt205)
567 auto a = GET_PARTICLE(vd);
570 if (vd.template getProp<type>(a) == BOUNDARY)
573 real_number rhop = vd.template getProp<rho>(a);
576 vd.template getProp<velocity>(a)[0] = 0.0;
577 vd.template getProp<velocity>(a)[1] = 0.0;
578 vd.template getProp<velocity>(a)[2] = 0.0;
579 real_number rhonew = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
580 vd.template getProp<rho>(a) = (rhonew < rho_zero)?rho_zero:rhonew;
582 vd.template getProp<rho_prev>(a) = rhop;
584 vd.template getProp<red>(a) = 0;
590 real_number dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
591 real_number dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
592 real_number dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
598 real_number velX = vd.template getProp<velocity>(a)[0];
599 real_number velY = vd.template getProp<velocity>(a)[1];
600 real_number velZ = vd.template getProp<velocity>(a)[2];
602 real_number rhop = vd.template getProp<rho>(a);
604 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity_prev>(a)[0] + vd.template getProp<force>(a)[0]*dt2;
605 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity_prev>(a)[1] + vd.template getProp<force>(a)[1]*dt2;
606 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity_prev>(a)[2] + vd.template getProp<force>(a)[2]*dt2;
607 vd.template getProp<rho>(a) = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
612 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
614 vd.template getProp<red>(a) = 1;
619 vd.template getProp<red>(a) = 0;
623 vd.template getProp<velocity_prev>(a)[0] = velX;
624 vd.template getProp<velocity_prev>(a)[1] = velY;
625 vd.template getProp<velocity_prev>(a)[2] = velZ;
626 vd.template getProp<rho_prev>(a) = rhop;
631 void verlet_int(
particles & vd, real_number dt)
634 auto part = vd.getDomainIteratorGPU();
636 real_number dt205 = dt*dt*0.5;
637 real_number dt2 = dt*2.0;
639 CUDA_LAUNCH(verlet_int_gpu,part,vd.toKernel(),dt,dt2,dt205);
642 remove_marked<red>(vd);
648 template<
typename vector_type>
649 __global__
void euler_int_gpu(
vector_type vd,real_number dt, real_number dt205)
652 auto a = GET_PARTICLE(vd);
655 if (vd.template getProp<type>(a) == BOUNDARY)
658 real_number rhop = vd.template getProp<rho>(a);
661 vd.template getProp<velocity>(a)[0] = 0.0;
662 vd.template getProp<velocity>(a)[1] = 0.0;
663 vd.template getProp<velocity>(a)[2] = 0.0;
664 real_number rhonew = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
665 vd.template getProp<rho>(a) = (rhonew < rho_zero)?rho_zero:rhonew;
667 vd.template getProp<rho_prev>(a) = rhop;
669 vd.template getProp<red>(a) = 0;
675 real_number dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
676 real_number dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
677 real_number dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
683 real_number velX = vd.template getProp<velocity>(a)[0];
684 real_number velY = vd.template getProp<velocity>(a)[1];
685 real_number velZ = vd.template getProp<velocity>(a)[2];
686 real_number rhop = vd.template getProp<rho>(a);
688 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity>(a)[0] + vd.template getProp<force>(a)[0]*dt;
689 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity>(a)[1] + vd.template getProp<force>(a)[1]*dt;
690 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity>(a)[2] + vd.template getProp<force>(a)[2]*dt;
691 vd.template getProp<rho>(a) = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
696 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
697 {vd.template getProp<red>(a) = 1;}
699 {vd.template getProp<red>(a) = 0;}
701 vd.template getProp<velocity_prev>(a)[0] = velX;
702 vd.template getProp<velocity_prev>(a)[1] = velY;
703 vd.template getProp<velocity_prev>(a)[2] = velZ;
704 vd.template getProp<rho_prev>(a) = rhop;
707 void euler_int(
particles & vd, real_number dt)
711 auto part = vd.getDomainIteratorGPU();
713 real_number dt205 = dt*dt*0.5;
715 CUDA_LAUNCH(euler_int_gpu,part,vd.toKernel(),dt,dt205);
718 remove_marked<red>(vd);
723 template<
typename vector_type,
typename NN_type>
726 real_number tot_ker = 0.0;
732 auto itg = NN.getNNIteratorBox(NN.getCell(xp));
735 auto q = itg.get_sort();
738 if (vd.template getProp<type>(q) != FLUID)
749 real_number r = sqrt(norm2(xp - xq));
751 real_number ker = Wab(r) * (MassFluid / rho_zero);
758 *press_tmp += vd.template getProp<Pressure>(q) * ker;
769 {*press_tmp = 1.0 / tot_ker * *press_tmp;}
772 template<
typename Vector,
typename CellList>
773 inline void sensor_pressure(
Vector & vd,
782 for (
size_t i = 0 ; i < probes.size() ; i++)
786 real_number press_tmp;
789 if (vd.getDecomposition().isLocal(probes.get(i)) ==
true)
791 vd.template updateCellListGPU<type,Pressure>(NN);
794 CUDA_LAUNCH_DIM3(sensor_pressure_gpu,1,1,vd.toKernel(),NN.toKernel(),probe,(real_number *)press_tmp_.toKernel());
796 vd.template restoreOrder<>(NN);
799 press_tmp_.deviceToHost();
800 press_tmp = *(real_number *)press_tmp_.getPointer();
810 press_t.last().add(press_tmp);
814 int main(
int argc,
char* argv[])
817 openfpm_init(&argc,&argv);
822 #ifdef CUDIFY_USE_CUDA
823 cudaDeviceSetCacheConfig(cudaFuncCachePreferL1);
830 probes.add({0.8779f,0.3f,0.02f});
831 probes.add({0.754f,0.31f,0.02f});
835 size_t sz[3] = {207,90,66};
838 W_dap = 1.0/Wab(H/1.5);
841 size_t bc[3]={NON_PERIODIC,NON_PERIODIC,NON_PERIODIC};
846 particles vd(0,domain,bc,g,DEC_GRAN(128));
852 Box<3,real_number> fluid_box({dp/2.0f,dp/2.0f,dp/2.0f},{0.4f+dp/2.0f,0.67f-dp/2.0f,0.3f+dp/2.0f});
858 max_fluid_height = fluid_it.getBoxMargins().
getHigh(2);
859 h_swl = fluid_it.getBoxMargins().
getHigh(2) - fluid_it.getBoxMargins().
getLow(2);
860 B = (coeff_sound)*(coeff_sound)*gravity*h_swl*rho_zero / gamma_;
861 cbar = coeff_sound * sqrt(gravity * h_swl);
864 while (fluid_it.isNext())
870 vd.getLastPos()[0] = fluid_it.get().get(0);
871 vd.getLastPos()[1] = fluid_it.get().get(1);
872 vd.getLastPos()[2] = fluid_it.get().get(2);
875 vd.template getLastProp<type>() = FLUID;
884 vd.template getLastProp<Pressure>() = rho_zero * gravity * (max_fluid_height - fluid_it.get().get(2));
886 vd.template getLastProp<rho>() = pow(vd.template getLastProp<Pressure>() / B + 1, 1.0/gamma_) * rho_zero;
887 vd.template getLastProp<rho_prev>() = vd.template getLastProp<rho>();
888 vd.template getLastProp<velocity>()[0] = 0.0;
889 vd.template getLastProp<velocity>()[1] = 0.0;
890 vd.template getLastProp<velocity>()[2] = 0.0;
892 vd.template getLastProp<velocity_prev>()[0] = 0.0;
893 vd.template getLastProp<velocity_prev>()[1] = 0.0;
894 vd.template getLastProp<velocity_prev>()[2] = 0.0;
901 Box<3,real_number> recipient1({0.0f,0.0f,0.0f},{1.6f+dp/2.0f,0.67f+dp/2.0f,0.4f+dp/2.0f});
904 Box<3,real_number> obstacle1({0.9f,0.24f-dp/2.0f,0.0f},{1.02f+dp/2.0f,0.36f,0.45f+dp/2.0f});
905 Box<3,real_number> obstacle2({0.9f+dp,0.24f+dp/2.0f,0.0f},{1.02f-dp/2.0f,0.36f-dp,0.45f-dp/2.0f});
909 holes.add(recipient2);
910 holes.add(obstacle1);
913 while (bound_box.isNext())
917 vd.getLastPos()[0] = bound_box.get().get(0);
918 vd.getLastPos()[1] = bound_box.get().get(1);
919 vd.getLastPos()[2] = bound_box.get().get(2);
921 vd.template getLastProp<type>() = BOUNDARY;
922 vd.template getLastProp<rho>() = rho_zero;
923 vd.template getLastProp<rho_prev>() = rho_zero;
924 vd.template getLastProp<velocity>()[0] = 0.0;
925 vd.template getLastProp<velocity>()[1] = 0.0;
926 vd.template getLastProp<velocity>()[2] = 0.0;
928 vd.template getLastProp<velocity_prev>()[0] = 0.0;
929 vd.template getLastProp<velocity_prev>()[1] = 0.0;
930 vd.template getLastProp<velocity_prev>()[2] = 0.0;
937 while (obstacle_box.isNext())
941 vd.getLastPos()[0] = obstacle_box.get().get(0);
942 vd.getLastPos()[1] = obstacle_box.get().get(1);
943 vd.getLastPos()[2] = obstacle_box.get().get(2);
945 vd.template getLastProp<type>() = BOUNDARY;
946 vd.template getLastProp<rho>() = rho_zero;
947 vd.template getLastProp<rho_prev>() = rho_zero;
948 vd.template getLastProp<velocity>()[0] = 0.0;
949 vd.template getLastProp<velocity>()[1] = 0.0;
950 vd.template getLastProp<velocity>()[2] = 0.0;
952 vd.template getLastProp<velocity_prev>()[0] = 0.0;
953 vd.template getLastProp<velocity_prev>()[1] = 0.0;
954 vd.template getLastProp<velocity_prev>()[2] = 0.0;
964 vd.addComputationCosts(md);
965 vd.getDecomposition().decompose();
972 vd.hostToDevicePos();
973 vd.template hostToDeviceProp<type,rho,rho_prev,Pressure,velocity>();
976 vd.ghost_get<type,rho,Pressure,velocity>(RUN_ON_DEVICE);
978 auto NN = vd.getCellListGPU(2*H / 2.0, CL_NON_SYMMETRIC | CL_GPU_REORDER, 2);
997 vd.map(RUN_ON_DEVICE);
1000 vd.deviceToHostPos();
1001 vd.template deviceToHostProp<type>();
1005 vd.addComputationCosts(md);
1006 vd.getDecomposition().decompose();
1009 {std::cout <<
"REBALANCED " << it_reb << std::endl;}
1012 vd.map(RUN_ON_DEVICE);
1017 real_number max_visc = 0.0;
1019 vd.ghost_get<type,rho,Pressure,velocity>(RUN_ON_DEVICE);
1023 calc_forces(vd,NN,max_visc,cnt,fluid_ids,border_ids);
1030 real_number dt = calc_deltaT(vd,max_visc);
1048 vd.map(RUN_ON_DEVICE);
1049 vd.ghost_get<type,rho,Pressure,velocity>(RUN_ON_DEVICE);
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;
1107 int 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.
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.
Vcluster< Memory > & getVC()
Get the Virtual Cluster machine.
size_t size_local() const
return the local size of the vector
const vector_dist_prop & getPropVector() const
return the property vector of all the particles
auto getPos(vect_dist_key_dx vec_key) -> decltype(vPos.template get< 0 >(vec_key.getKey()))
Get the position of an element.
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