OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
 
Loading...
Searching...
No Matches
main.cpp
1
60//#define SE_CLASS1
61//#define STOP_ON_ERROR
62
64#include "Vector/vector_dist.hpp"
65#include <math.h>
66#include "Draw/DrawParticles.hpp"
68
113// A constant to indicate boundary particles
114#define BOUNDARY 0
115
116// A constant to indicate fluid particles
117#define FLUID 1
118
119// initial spacing between particles dp in the formulas
120const double dp = 0.0085;
121// Maximum height of the fluid water
122// is going to be calculated and filled later on
123double h_swl = 0.0;
124
125// c_s in the formulas (constant used to calculate the sound speed)
126const double coeff_sound = 20.0;
127
128// gamma in the formulas
129const double gamma_ = 7.0;
130
131// sqrt(3.0*dp*dp) support of the kernel
132const double H = 0.0147224318643;
133
134// Eta in the formulas
135const double Eta2 = 0.01 * H*H;
136
137// alpha in the formula
138const double visco = 0.1;
139
140// cbar in the formula (calculated later)
141double cbar = 0.0;
142
143// Mass of the fluid particles
144const double MassFluid = 0.000614125;
145
146// Mass of the boundary particles
147const double MassBound = 0.000614125;
148
149// End simulation time
150#ifdef TEST_RUN
151const double t_end = 0.001;
152#else
153const double t_end = 1.5;
154#endif
155
156// Gravity acceleration
157const double gravity = 9.81;
158
159// Reference densitu 1000Kg/m^3
160const double rho_zero = 1000.0;
161
162// Filled later require h_swl, it is b in the formulas
163double B = 0.0;
164
165// Constant used to define time integration
166const double CFLnumber = 0.2;
167
168// Minimum T
169const double DtMin = 0.00001;
170
171// Minimum Rho allowed
172const double RhoMin = 700.0;
173
174// Maximum Rho allowed
175const double RhoMax = 1300.0;
176
177// Filled in initialization
178double max_fluid_height = 0.0;
179
180// Properties
181
182// FLUID or BOUNDARY
183const size_t type = 0;
184
185// Density
186const int rho = 1;
187
188// Density at step n-1
189const int rho_prev = 2;
190
191// Pressure
192const int Pressure = 3;
193
194// Delta rho calculated in the force calculation
195const int drho = 4;
196
197// calculated force
198const int force = 5;
199
200// velocity
201const int velocity = 6;
202
203// velocity at previous step
204const int velocity_prev = 7;
205
210// Type of the vector containing particles
212// | | | | | | | |
213// | | | | | | | |
214// type density density Pressure delta force velocity velocity
215// at n-1 density at n - 1
216
222{
223 template<typename Decomposition, typename vector> inline void addComputation(Decomposition & dec,
224 vector & vd,
225 size_t v,
226 size_t p)
227 {
228 if (vd.template getProp<type>(p) == FLUID)
229 dec.addComputationCost(v,4);
230 else
231 dec.addComputationCost(v,3);
232 }
233
234 template<typename Decomposition> inline void applyModel(Decomposition & dec, size_t v)
235 {
236 dec.setSubSubDomainComputationCost(v, dec.getSubSubDomainComputationCost(v) * dec.getSubSubDomainComputationCost(v));
237 }
238
239 double distributionTol()
240 {
241 return 1.01;
242 }
243};
244
262inline void EqState(particles & vd)
263{
264 auto it = vd.getDomainIterator();
265
266 while (it.isNext())
267 {
268 auto a = it.get();
269
270 double rho_a = vd.template getProp<rho>(a);
271 double rho_frac = rho_a / rho_zero;
272
273 vd.template getProp<Pressure>(a) = B*( rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac - 1.0);
274
275 ++it;
276 }
277}
278
297const double a2 = 1.0/M_PI/H/H/H;
298
299inline double Wab(double r)
300{
301 r /= H;
302
303 if (r < 1.0)
304 return (1.0 - 3.0/2.0*r*r + 3.0/4.0*r*r*r)*a2;
305 else if (r < 2.0)
306 return (1.0/4.0*(2.0 - r*r)*(2.0 - r*r)*(2.0 - r*r))*a2;
307 else
308 return 0.0;
309}
310
328const double c1 = -3.0/M_PI/H/H/H/H;
329const double d1 = 9.0/4.0/M_PI/H/H/H/H;
330const double c2 = -3.0/4.0/M_PI/H/H/H/H;
331const double a2_4 = 0.25*a2;
332// Filled later
333double W_dap = 0.0;
334
335inline void DWab(Point<3,double> & dx, Point<3,double> & DW, double r, bool print)
336{
337 const double qq=r/H;
338
339 double qq2 = qq * qq;
340 double fac1 = (c1*qq + d1*qq2)/r;
341 double b1 = (qq < 1.0)?1.0f:0.0f;
342
343 double wqq = (2.0 - qq);
344 double fac2 = c2 * wqq * wqq / r;
345 double b2 = (qq >= 1.0 && qq < 2.0)?1.0f:0.0f;
346
347 double factor = (b1*fac1 + b2*fac2);
348
349 DW.get(0) = factor * dx.get(0);
350 DW.get(1) = factor * dx.get(1);
351 DW.get(2) = factor * dx.get(2);
352}
353
374// Tensile correction
375inline double Tensile(double r, double rhoa, double rhob, double prs1, double prs2)
376{
377 const double qq=r/H;
378 //-Cubic Spline kernel
379 double wab;
380 if(r>H)
381 {
382 double wqq1=2.0f-qq;
383 double wqq2=wqq1*wqq1;
384
385 wab=a2_4*(wqq2*wqq1);
386 }
387 else
388 {
389 double wqq2=qq*qq;
390 double wqq3=wqq2*qq;
391
392 wab=a2*(1.0f-1.5f*wqq2+0.75f*wqq3);
393 }
394
395 //-Tensile correction.
396 double fab=wab*W_dap;
397 fab*=fab; fab*=fab; //fab=fab^4
398 const double tensilp1=(prs1/(rhoa*rhoa))*(prs1>0? 0.01: -0.2);
399 const double tensilp2=(prs2/(rhob*rhob))*(prs2>0? 0.01: -0.2);
400
401 return (fab*(tensilp1+tensilp2));
402}
403
422inline double Pi(const Point<3,double> & dr, double rr2, Point<3,double> & dv, double rhoa, double rhob, double massb, double & visc)
423{
424 const double dot = dr.get(0)*dv.get(0) + dr.get(1)*dv.get(1) + dr.get(2)*dv.get(2);
425 const double dot_rr2 = dot/(rr2+Eta2);
426 visc=std::max(dot_rr2,visc);
427
428 if(dot < 0)
429 {
430 const float amubar=H*dot_rr2;
431 const float robar=(rhoa+rhob)*0.5f;
432 const float pi_visc=(-visco*cbar*amubar/robar);
433
434 return pi_visc;
435 }
436 else
437 return 0.0;
438}
439
457template<typename CellList> inline void calc_forces(particles & vd, CellList & NN, double & max_visc)
458{
459 auto part = vd.getDomainIterator();
460
461 // Update the cell-list
462 vd.updateCellList(NN);
463
464 // For each particle ...
465 while (part.isNext())
466 {
467 // ... a
468 auto a = part.get();
469
470 // Get the position xp of the particle
471 Point<3,double> xa = vd.getPos(a);
472
473 // Take the mass of the particle dependently if it is FLUID or BOUNDARY
474 double massa = (vd.getProp<type>(a) == FLUID)?MassFluid:MassBound;
475
476 // Get the density of the of the particle a
477 double rhoa = vd.getProp<rho>(a);
478
479 // Get the pressure of the particle a
480 double Pa = vd.getProp<Pressure>(a);
481
482 // Get the Velocity of the particle a
483 Point<3,double> va = vd.getProp<velocity>(a);
484
485 // Reset the force counter (- gravity on zeta direction)
486 vd.template getProp<force>(a)[0] = 0.0;
487 vd.template getProp<force>(a)[1] = 0.0;
488 vd.template getProp<force>(a)[2] = -gravity;
489 vd.template getProp<drho>(a) = 0.0;
490
491 // We threat FLUID particle differently from BOUNDARY PARTICLES ...
492 if (vd.getProp<type>(a) != FLUID)
493 {
494 // If it is a boundary particle calculate the delta rho based on equation 2
495 // This require to run across the neighborhoods particles of a
496 auto Np = NN.template getNNIterator<NO_CHECK>(NN.getCell(vd.getPos(a)));
497
498 // For each neighborhood particle
499 while (Np.isNext() == true)
500 {
501 // ... q
502 auto b = Np.get();
503
504 // Get the position xp of the particle
505 Point<3,double> xb = vd.getPos(b);
506
507 // if (p == q) skip this particle
508 if (a.getKey() == b) {++Np; continue;};
509
510 // get the mass of the particle
511 double massb = (vd.getProp<type>(b) == FLUID)?MassFluid:MassBound;
512
513 // Get the velocity of the particle b
514 Point<3,double> vb = vd.getProp<velocity>(b);
515
516 // Get the pressure and density of particle b
517 double Pb = vd.getProp<Pressure>(b);
518 double rhob = vd.getProp<rho>(b);
519
520 // Get the distance between p and q
521 Point<3,double> dr = xa - xb;
522 // take the norm of this vector
523 double r2 = norm2(dr);
524
525 // If the particles interact ...
526 if (r2 < 4.0*H*H)
527 {
528 // ... calculate delta rho
529 double r = sqrt(r2);
530
531 Point<3,double> dv = va - vb;
532
534 DWab(dr,DW,r,false);
535
536 const double dot = dr.get(0)*dv.get(0) + dr.get(1)*dv.get(1) + dr.get(2)*dv.get(2);
537 const double dot_rr2 = dot/(r2+Eta2);
538 max_visc=std::max(dot_rr2,max_visc);
539
540 vd.getProp<drho>(a) += massb*(dv.get(0)*DW.get(0)+dv.get(1)*DW.get(1)+dv.get(2)*DW.get(2));
541 }
542
543 ++Np;
544 }
545 }
546 else
547 {
548 // If it is a fluid particle calculate based on equation 1 and 2
549
550 // Get an iterator over the neighborhood particles of p
551 auto Np = NN.template getNNIterator<NO_CHECK>(NN.getCell(vd.getPos(a)));
552
553 // For each neighborhood particle
554 while (Np.isNext() == true)
555 {
556 // ... q
557 auto b = Np.get();
558
559 // Get the position xp of the particle
560 Point<3,double> xb = vd.getPos(b);
561
562 // if (p == q) skip this particle
563 if (a.getKey() == b) {++Np; continue;};
564
565 double massb = (vd.getProp<type>(b) == FLUID)?MassFluid:MassBound;
566 Point<3,double> vb = vd.getProp<velocity>(b);
567 double Pb = vd.getProp<Pressure>(b);
568 double rhob = vd.getProp<rho>(b);
569
570 // Get the distance between p and q
571 Point<3,double> dr = xa - xb;
572 // take the norm of this vector
573 double r2 = norm2(dr);
574
575 // if they interact
576 if (r2 < 4.0*H*H)
577 {
578 double r = sqrt(r2);
579
580 Point<3,double> v_rel = va - vb;
581
583 DWab(dr,DW,r,false);
584
585 double factor = - massb*((vd.getProp<Pressure>(a) + vd.getProp<Pressure>(b)) / (rhoa * rhob) + Tensile(r,rhoa,rhob,Pa,Pb) + Pi(dr,r2,v_rel,rhoa,rhob,massb,max_visc));
586
587 vd.getProp<force>(a)[0] += factor * DW.get(0);
588 vd.getProp<force>(a)[1] += factor * DW.get(1);
589 vd.getProp<force>(a)[2] += factor * DW.get(2);
590
591 vd.getProp<drho>(a) += massb*(v_rel.get(0)*DW.get(0)+v_rel.get(1)*DW.get(1)+v_rel.get(2)*DW.get(2));
592 }
593
594 ++Np;
595 }
596 }
597
598 ++part;
599 }
600}
601
620void max_acceleration_and_velocity(particles & vd, double & max_acc, double & max_vel)
621{
622 // Calculate the maximum acceleration
623 auto part = vd.getDomainIterator();
624
625 while (part.isNext())
626 {
627 auto a = part.get();
628
629 Point<3,double> acc(vd.getProp<force>(a));
630 double acc2 = norm2(acc);
631
632 Point<3,double> vel(vd.getProp<velocity>(a));
633 double vel2 = norm2(vel);
634
635 if (vel2 >= max_vel)
636 max_vel = vel2;
637
638 if (acc2 >= max_acc)
639 max_acc = acc2;
640
641 ++part;
642 }
643 max_acc = sqrt(max_acc);
644 max_vel = sqrt(max_vel);
645
646 Vcluster<> & v_cl = create_vcluster();
647 v_cl.max(max_acc);
648 v_cl.max(max_vel);
649 v_cl.execute();
650}
651
677double calc_deltaT(particles & vd, double ViscDtMax)
678{
679 double Maxacc = 0.0;
680 double Maxvel = 0.0;
681 max_acceleration_and_velocity(vd,Maxacc,Maxvel);
682
683 //-dt1 depends on force per unit mass.
684 const double dt_f = (Maxacc)?sqrt(H/Maxacc):std::numeric_limits<int>::max();
685
686 //-dt2 combines the Courant and the viscous time-step controls.
687 const double dt_cv = H/(std::max(cbar,Maxvel*10.) + H*ViscDtMax);
688
689 //-dt new value of time step.
690 double dt=double(CFLnumber)*std::min(dt_f,dt_cv);
691 if(dt<double(DtMin))
692 dt=double(DtMin);
693
694 return dt;
695}
696
731
732size_t cnt = 0;
733
734void verlet_int(particles & vd, double dt)
735{
736 // list of the particle to remove
737 to_remove.clear();
738
739 // particle iterator
740 auto part = vd.getDomainIterator();
741
742 double dt205 = dt*dt*0.5;
743 double dt2 = dt*2.0;
744
745 // For each particle ...
746 while (part.isNext())
747 {
748 // ... a
749 auto a = part.get();
750
751 // if the particle is boundary
752 if (vd.template getProp<type>(a) == BOUNDARY)
753 {
754 // Update rho
755 double rhop = vd.template getProp<rho>(a);
756
757 // Update only the density
758 vd.template getProp<velocity>(a)[0] = 0.0;
759 vd.template getProp<velocity>(a)[1] = 0.0;
760 vd.template getProp<velocity>(a)[2] = 0.0;
761 double rhonew = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
762 vd.template getProp<rho>(a) = (rhonew < rho_zero)?rho_zero:rhonew;
763
764 vd.template getProp<rho_prev>(a) = rhop;
765
766 ++part;
767 continue;
768 }
769
770 //-Calculate displacement and update position / Calcula desplazamiento y actualiza posicion.
771 double dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
772 double dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
773 double dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
774
775 vd.getPos(a)[0] += dx;
776 vd.getPos(a)[1] += dy;
777 vd.getPos(a)[2] += dz;
778
779 double velX = vd.template getProp<velocity>(a)[0];
780 double velY = vd.template getProp<velocity>(a)[1];
781 double velZ = vd.template getProp<velocity>(a)[2];
782 double rhop = vd.template getProp<rho>(a);
783
784 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity_prev>(a)[0] + vd.template getProp<force>(a)[0]*dt2;
785 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity_prev>(a)[1] + vd.template getProp<force>(a)[1]*dt2;
786 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity_prev>(a)[2] + vd.template getProp<force>(a)[2]*dt2;
787 vd.template getProp<rho>(a) = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
788
789 // Check if the particle go out of range in space and in density
790 if (vd.getPos(a)[0] < 0.000263878 || vd.getPos(a)[1] < 0.000263878 || vd.getPos(a)[2] < 0.000263878 ||
791 vd.getPos(a)[0] > 0.000263878+1.59947 || vd.getPos(a)[1] > 0.000263878+0.672972 || vd.getPos(a)[2] > 0.000263878+0.903944 ||
792 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
793 {
794 to_remove.add(a.getKey());
795 }
796
797 vd.template getProp<velocity_prev>(a)[0] = velX;
798 vd.template getProp<velocity_prev>(a)[1] = velY;
799 vd.template getProp<velocity_prev>(a)[2] = velZ;
800 vd.template getProp<rho_prev>(a) = rhop;
801
802 ++part;
803 }
804
805 // remove the particles
806 vd.remove(to_remove,0);
807
808 // increment the iteration counter
809 cnt++;
810}
811
812void euler_int(particles & vd, double dt)
813{
814 // list of the particle to remove
815 to_remove.clear();
816
817 // particle iterator
818 auto part = vd.getDomainIterator();
819
820 double dt205 = dt*dt*0.5;
821 double dt2 = dt*2.0;
822
823 // For each particle ...
824 while (part.isNext())
825 {
826 // ... a
827 auto a = part.get();
828
829 // if the particle is boundary
830 if (vd.template getProp<type>(a) == BOUNDARY)
831 {
832 // Update rho
833 double rhop = vd.template getProp<rho>(a);
834
835 // Update only the density
836 vd.template getProp<velocity>(a)[0] = 0.0;
837 vd.template getProp<velocity>(a)[1] = 0.0;
838 vd.template getProp<velocity>(a)[2] = 0.0;
839 double rhonew = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
840 vd.template getProp<rho>(a) = (rhonew < rho_zero)?rho_zero:rhonew;
841
842 vd.template getProp<rho_prev>(a) = rhop;
843
844 ++part;
845 continue;
846 }
847
848 //-Calculate displacement and update position / Calcula desplazamiento y actualiza posicion.
849 double dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
850 double dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
851 double dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
852
853 vd.getPos(a)[0] += dx;
854 vd.getPos(a)[1] += dy;
855 vd.getPos(a)[2] += dz;
856
857 double velX = vd.template getProp<velocity>(a)[0];
858 double velY = vd.template getProp<velocity>(a)[1];
859 double velZ = vd.template getProp<velocity>(a)[2];
860 double rhop = vd.template getProp<rho>(a);
861
862 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity>(a)[0] + vd.template getProp<force>(a)[0]*dt;
863 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity>(a)[1] + vd.template getProp<force>(a)[1]*dt;
864 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity>(a)[2] + vd.template getProp<force>(a)[2]*dt;
865 vd.template getProp<rho>(a) = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
866
867 // Check if the particle go out of range in space and in density
868 if (vd.getPos(a)[0] < 0.000263878 || vd.getPos(a)[1] < 0.000263878 || vd.getPos(a)[2] < 0.000263878 ||
869 vd.getPos(a)[0] > 0.000263878+1.59947 || vd.getPos(a)[1] > 0.000263878+0.672972 || vd.getPos(a)[2] > 0.000263878+0.903944 ||
870 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
871 {
872 to_remove.add(a.getKey());
873 }
874
875 vd.template getProp<velocity_prev>(a)[0] = velX;
876 vd.template getProp<velocity_prev>(a)[1] = velY;
877 vd.template getProp<velocity_prev>(a)[2] = velZ;
878 vd.template getProp<rho_prev>(a) = rhop;
879
880 ++part;
881 }
882
883 // remove the particles
884 vd.remove(to_remove,0);
885
886 // increment the iteration counter
887 cnt++;
888}
889
917template<typename Vector, typename CellList>
918inline void sensor_pressure(Vector & vd,
919 CellList & NN,
922{
923 Vcluster<> & v_cl = create_vcluster();
924
925 press_t.add();
926
927 for (size_t i = 0 ; i < probes.size() ; i++)
928 {
929 float press_tmp = 0.0f;
930 float tot_ker = 0.0;
931
932 // if the probe is inside the processor domain
933 if (vd.getDecomposition().isLocal(probes.get(i)) == true)
934 {
935 // Get the position of the probe i
936 Point<3,double> xp = probes.get(i);
937
938 // get the iterator over the neighbohood particles of the probes position
939 auto itg = NN.template getNNIterator<NO_CHECK>(NN.getCell(probes.get(i)));
940 while (itg.isNext())
941 {
942 auto q = itg.get();
943
944 // Only the fluid particles are importants
945 if (vd.template getProp<type>(q) != FLUID)
946 {
947 ++itg;
948 continue;
949 }
950
951 // Get the position of the neighborhood particle q
952 Point<3,double> xq = vd.getPos(q);
953
954 // Calculate the contribution of the particle to the pressure
955 // of the probe
956 double r = sqrt(norm2(xp - xq));
957
958 double ker = Wab(r) * (MassFluid / rho_zero);
959
960 // Also keep track of the calculation of the summed
961 // kernel
962 tot_ker += ker;
963
964 // Add the total pressure contribution
965 press_tmp += vd.template getProp<Pressure>(q) * ker;
966
967 // next neighborhood particle
968 ++itg;
969 }
970
971 // We calculate the pressure normalizing the
972 // sum over all kernels
973 if (tot_ker == 0.0)
974 press_tmp = 0.0;
975 else
976 press_tmp = 1.0 / tot_ker * press_tmp;
977
978 }
979
980 // This is not necessary in principle, but if you
981 // want to make all processor aware of the history of the calculated
982 // pressure we have to execute this
983 v_cl.sum(press_tmp);
984 v_cl.execute();
985
986 // We add the calculated pressure into the history
987 press_t.last().add(press_tmp);
988 }
989}
990
993int main(int argc, char* argv[])
994{
1011
1012 // initialize the library
1013 openfpm_init(&argc,&argv);
1014
1015 // It contain for each time-step the value detected by the probes
1018
1019 probes.add({0.8779,0.3,0.02});
1020 probes.add({0.754,0.31,0.02});
1021
1022 // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
1023 Box<3,double> domain({-0.05,-0.05,-0.05},{1.7010,0.7065,0.5025});
1024 size_t sz[3] = {207,90,66};
1025
1026 // Fill W_dap
1027 W_dap = 1.0/Wab(H/1.5);
1028
1029 // Here we define the boundary conditions of our problem
1030 size_t bc[3]={NON_PERIODIC,NON_PERIODIC,NON_PERIODIC};
1031
1032 // extended boundary around the domain, and the processor domain
1033 Ghost<3,double> g(2*H);
1034
1036
1073
1074 particles vd(0,domain,bc,g,DEC_GRAN(512));
1075
1077
1113
1114 // You can ignore all these dp/2.0 is a trick to reach the same initialization
1115 // of Dual-SPH that use a different criteria to draw particles
1116 Box<3,double> 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});
1117
1118 // return an iterator to the fluid particles to add to vd
1119 auto fluid_it = DrawParticles::DrawBox(vd,sz,domain,fluid_box);
1120
1121 // here we fill some of the constants needed by the simulation
1122 max_fluid_height = fluid_it.getBoxMargins().getHigh(2);
1123 h_swl = fluid_it.getBoxMargins().getHigh(2) - fluid_it.getBoxMargins().getLow(2);
1124 B = (coeff_sound)*(coeff_sound)*gravity*h_swl*rho_zero / gamma_;
1125 cbar = coeff_sound * sqrt(gravity * h_swl);
1126
1127 // for each particle inside the fluid box ...
1128 while (fluid_it.isNext())
1129 {
1130 // ... add a particle ...
1131 vd.add();
1132
1133 // ... and set it position ...
1134 vd.getLastPos()[0] = fluid_it.get().get(0);
1135 vd.getLastPos()[1] = fluid_it.get().get(1);
1136 vd.getLastPos()[2] = fluid_it.get().get(2);
1137
1138 // and its type.
1139 vd.template getLastProp<type>() = FLUID;
1140
1141 // We also initialize the density of the particle and the hydro-static pressure given by
1142 //
1143 // rho_zero*g*h = P
1144 //
1145 // rho_p = (P/B + 1)^(1/Gamma) * rho_zero
1146 //
1147
1148 vd.template getLastProp<Pressure>() = rho_zero * gravity * (max_fluid_height - fluid_it.get().get(2));
1149
1150 vd.template getLastProp<rho>() = pow(vd.template getLastProp<Pressure>() / B + 1, 1.0/gamma_) * rho_zero;
1151 vd.template getLastProp<rho_prev>() = vd.template getLastProp<rho>();
1152 vd.template getLastProp<velocity>()[0] = 0.0;
1153 vd.template getLastProp<velocity>()[1] = 0.0;
1154 vd.template getLastProp<velocity>()[2] = 0.0;
1155
1156 vd.template getLastProp<velocity_prev>()[0] = 0.0;
1157 vd.template getLastProp<velocity_prev>()[1] = 0.0;
1158 vd.template getLastProp<velocity_prev>()[2] = 0.0;
1159
1160 // next fluid particle
1161 ++fluid_it;
1162 }
1163
1165
1190
1191 // Recipient
1192 Box<3,double> recipient1({0.0,0.0,0.0},{1.6+dp/2.0,0.67+dp/2.0,0.4+dp/2.0});
1193 Box<3,double> recipient2({dp,dp,dp},{1.6-dp/2.0,0.67-dp/2.0,0.4+dp/2.0});
1194
1195 Box<3,double> obstacle1({0.9,0.24-dp/2.0,0.0},{1.02+dp/2.0,0.36,0.45+dp/2.0});
1196 Box<3,double> obstacle2({0.9+dp,0.24+dp/2.0,0.0},{1.02-dp/2.0,0.36-dp,0.45-dp/2.0});
1197 Box<3,double> obstacle3({0.9+dp,0.24,0.0},{1.02,0.36,0.45});
1198
1200 holes.add(recipient2);
1201 holes.add(obstacle1);
1202 auto bound_box = DrawParticles::DrawSkin(vd,sz,domain,holes,recipient1);
1203
1204 while (bound_box.isNext())
1205 {
1206 vd.add();
1207
1208 vd.getLastPos()[0] = bound_box.get().get(0);
1209 vd.getLastPos()[1] = bound_box.get().get(1);
1210 vd.getLastPos()[2] = bound_box.get().get(2);
1211
1212 vd.template getLastProp<type>() = BOUNDARY;
1213 vd.template getLastProp<rho>() = rho_zero;
1214 vd.template getLastProp<rho_prev>() = rho_zero;
1215 vd.template getLastProp<velocity>()[0] = 0.0;
1216 vd.template getLastProp<velocity>()[1] = 0.0;
1217 vd.template getLastProp<velocity>()[2] = 0.0;
1218
1219 vd.template getLastProp<velocity_prev>()[0] = 0.0;
1220 vd.template getLastProp<velocity_prev>()[1] = 0.0;
1221 vd.template getLastProp<velocity_prev>()[2] = 0.0;
1222
1223 ++bound_box;
1224 }
1225
1227
1245
1246 auto obstacle_box = DrawParticles::DrawSkin(vd,sz,domain,obstacle2,obstacle1);
1247
1248 while (obstacle_box.isNext())
1249 {
1250 vd.add();
1251
1252 vd.getLastPos()[0] = obstacle_box.get().get(0);
1253 vd.getLastPos()[1] = obstacle_box.get().get(1);
1254 vd.getLastPos()[2] = obstacle_box.get().get(2);
1255
1256 vd.template getLastProp<type>() = BOUNDARY;
1257 vd.template getLastProp<rho>() = rho_zero;
1258 vd.template getLastProp<rho_prev>() = rho_zero;
1259 vd.template getLastProp<velocity>()[0] = 0.0;
1260 vd.template getLastProp<velocity>()[1] = 0.0;
1261 vd.template getLastProp<velocity>()[2] = 0.0;
1262
1263 vd.template getLastProp<velocity_prev>()[0] = 0.0;
1264 vd.template getLastProp<velocity_prev>()[1] = 0.0;
1265 vd.template getLastProp<velocity_prev>()[2] = 0.0;
1266
1267 ++obstacle_box;
1268 }
1269
1270 vd.map();
1271
1273
1356
1357 // Now that we fill the vector with particles
1358 ModelCustom md;
1359
1360 vd.addComputationCosts(md);
1361 vd.getDecomposition().decompose();
1362 vd.map();
1363
1365
1366 vd.ghost_get<type,rho,Pressure,velocity>();
1367
1368 auto NN = vd.getCellList(2*H);
1369
1370 // Evolve
1371
1387
1388 size_t write = 0;
1389 size_t it = 0;
1390 size_t it_reb = 0;
1391 double t = 0.0;
1392 while (t <= t_end)
1393 {
1394 Vcluster<> & v_cl = create_vcluster();
1395 timer it_time;
1396
1398 it_reb++;
1399 if (it_reb == 200)
1400 {
1401 vd.map();
1402
1403 it_reb = 0;
1404 ModelCustom md;
1405 vd.addComputationCosts(md);
1406 vd.getDecomposition().decompose();
1407
1408 if (v_cl.getProcessUnitID() == 0)
1409 std::cout << "REBALANCED " << std::endl;
1410 }
1411
1412 vd.map();
1413
1414 // Calculate pressure from the density
1415 EqState(vd);
1416
1417 double max_visc = 0.0;
1418
1419 vd.ghost_get<type,rho,Pressure,velocity>();
1420
1421 // Calc forces
1422 calc_forces(vd,NN,max_visc);
1423
1424 // Get the maximum viscosity term across processors
1425 v_cl.max(max_visc);
1426 v_cl.execute();
1427
1428 // Calculate delta t integration
1429 double dt = calc_deltaT(vd,max_visc);
1430
1431 // VerletStep or euler step
1432 it++;
1433 if (it < 40)
1434 verlet_int(vd,dt);
1435 else
1436 {
1437 euler_int(vd,dt);
1438 it = 0;
1439 }
1440
1441 t += dt;
1442
1443 if (write < t*100)
1444 {
1445 // sensor_pressure calculation require ghost and update cell-list
1446 vd.map();
1447 vd.ghost_get<type,rho,Pressure,velocity>();
1448 vd.updateCellList(NN);
1449
1450 // calculate the pressure at the sensor points
1451 sensor_pressure(vd,NN,press_t,probes);
1452
1453 vd.write_frame("Geometry",write);
1454 write++;
1455
1456 if (v_cl.getProcessUnitID() == 0)
1457 {std::cout << "TIME: " << t << " write " << it_time.getwct() << " " << v_cl.getProcessUnitID() << " " << cnt << " Max visc: " << max_visc << std::endl;}
1458 }
1459 else
1460 {
1461 if (v_cl.getProcessUnitID() == 0)
1462 {std::cout << "TIME: " << t << " " << it_time.getwct() << " " << v_cl.getProcessUnitID() << " " << cnt << " Max visc: " << max_visc << std::endl;}
1463 }
1464 }
1465
1467
1482
1483 openfpm_finalize();
1484
1486
1495}
1496
This class represent an N-dimensional box.
Definition Box.hpp:61
__device__ __host__ T getLow(int i) const
get the i-coordinate of the low bound interval of the box
Definition Box.hpp:556
__device__ __host__ T getHigh(int i) const
get the high interval of the box
Definition Box.hpp:567
Class for FAST cell list implementation.
Definition CellList.hpp:357
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.
Definition Point.hpp:28
__device__ __host__ const T & get(unsigned int i) const
Get coordinate.
Definition Point.hpp:172
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.
Definition VCluster.hpp:59
Sparse Matrix implementation stub object when OpenFPM is compiled with no linear algebra support.
Definition Vector.hpp:40
Implementation of 1-D std::vector like structure.
Class for cpu time benchmarking.
Definition timer.hpp:28
double getwct()
Return the elapsed real time.
Definition timer.hpp:130
Distributed vector.
void remove(openfpm::vector< size_t > &keys, size_t start=0)
Remove a set of elements from the distributed vector.
auto getProp(vect_dist_key_dx vec_key) -> decltype(v_prp.template get< id >(vec_key.getKey()))
Get the property of an element.
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.
vector_dist_iterator getDomainIterator() const
Get an iterator that traverse the particles in the domain.
Model for Dynamic load balancing.
Definition main.cpp:222