OpenFPM_pdata  4.1.0
Project that contain the implementation of distributed structures
 
Loading...
Searching...
No Matches
main_dbg.cpp
1
2//#define SE_CLASS1
3//#define STOP_ON_ERROR
4
5#include "Vector/vector_dist.hpp"
6#include <math.h>
7#include "Draw/DrawParticles.hpp"
8
9
10// A constant to indicate boundary particles
11#define BOUNDARY 0
12
13// A constant to indicate fluid particles
14#define FLUID 1
15
16size_t cnt = 0;
17
18// initial spacing between particles dp in the formulas
19const double dp = 0.0085;
20// Maximum height of the fluid water
21// is going to be calculated and filled later on
22double h_swl = 0.0;
23
24// c_s in the formulas (constant used to calculate the sound speed)
25const double coeff_sound = 20.0;
26
27// gamma in the formulas
28const double gamma_ = 7.0;
29
30// sqrt(3.0*dp*dp) support of the kernel
31const double H = 0.0147224318643;
32
33// Eta in the formulas
34const double Eta2 = 0.01 * H*H;
35
36// alpha in the formula
37const double visco = 0.1;
38
39// cbar in the formula (calculated later)
40double cbar = 0.0;
41
42// Mass of the fluid particles
43const double MassFluid = 0.000614125;
44
45// Mass of the boundary particles
46const double MassBound = 0.000614125;
47
48// End simulation time
49const double t_end = 1.5;
50
51// Gravity acceleration
52const double gravity = 9.81;
53
54// Reference densitu 1000Kg/m^3
55const double rho_zero = 1000.0;
56
57// Filled later require h_swl, it is b in the formulas
58double B = 0.0;
59
60// Constant used to define time integration
61const double CFLnumber = 0.2;
62
63// Minimum T
64const double DtMin = 0.00001;
65
66// Minimum Rho allowed
67const double RhoMin = 700.0;
68
69// Maximum Rho allowed
70const double RhoMax = 1300.0;
71
72// Filled in initialization
73double max_fluid_height = 0.0;
74
75// Properties
76
77// FLUID or BOUNDARY
78const size_t type = 0;
79
80// Density
81const int rho = 1;
82
83// Density at step n-1
84const int rho_prev = 2;
85
86// Pressure
87const int Pressure = 3;
88
89// Delta rho calculated in the force calculation
90const int drho = 4;
91
92// calculated force
93const int force = 5;
94
95// velocity
96const int velocity = 6;
97
98// velocity at previous step
99const int velocity_prev = 7;
100
101struct nb_f
102{
103 size_t id;
104 double fact;
105
106 // Used to reorder the neighborhood particles by id
107 bool operator<(const struct nb_f & pag) const
108 {
109 return (id < pag.id);
110 }
111};
112
113// Type of the vector containing particles
115// | | | | | | | |
116// | | | | | | | |
117// type density density Pressure delta force velocity velocity
118// at n-1 density at n - 1
119
120struct ModelCustom
121{
122 template<typename Decomposition, typename vector> inline void addComputation(Decomposition & dec, const vector & vd, size_t v, size_t p)
123 {
124 if (vd.template getProp<type>(p) == FLUID)
125 dec.addComputationCost(v,4);
126 else
127 dec.addComputationCost(v,3);
128 }
129
130 template<typename Decomposition> inline void applyModel(Decomposition & dec, size_t v)
131 {
132 dec.setSubSubDomainComputationCost(v, dec.getSubSubDomainComputationCost(v) * dec.getSubSubDomainComputationCost(v));
133 }
134};
135
136inline void EqState(particles & vd)
137{
138 auto it = vd.getDomainIterator();
139
140 while (it.isNext())
141 {
142 auto a = it.get();
143
144 double rho_a = vd.template getProp<rho>(a);
145 double rho_frac = rho_a / rho_zero;
146
147 vd.template getProp<Pressure>(a) = B*( rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac*rho_frac - 1.0);
148
149 ++it;
150 }
151}
152
153
154
155const double a2 = 1.0/M_PI/H/H/H;
156
157inline double Wab(double r)
158{
159 r /= H;
160
161 if (r < 1.0)
162 return (1.0 - 3.0/2.0*r*r + 3.0/4.0*r*r*r)*a2;
163 else if (r < 2.0)
164 return (1.0/4.0*(2.0 - r*r)*(2.0 - r*r)*(2.0 - r*r))*a2;
165 else
166 return 0.0;
167}
168
169
170
171const double c1 = -3.0/M_PI/H/H/H/H;
172const double d1 = 9.0/4.0/M_PI/H/H/H/H;
173const double c2 = -3.0/4.0/M_PI/H/H/H/H;
174const double a2_4 = 0.25*a2;
175// Filled later
176double W_dap = 0.0;
177
178inline void DWab(Point<3,double> & dx, Point<3,double> & DW, double r, bool print)
179{
180 const double qq=r/H;
181
182 if (qq < 1.0)
183 {
184 double qq2 = qq * qq;
185 double fac = (c1*qq + d1*qq2)/r;
186
187 DW.get(0) = fac*dx.get(0);
188 DW.get(1) = fac*dx.get(1);
189 DW.get(2) = fac*dx.get(2);
190 }
191 else if (qq < 2.0)
192 {
193 double wqq = (2.0 - qq);
194 double fac = c2 * wqq * wqq / r;
195
196 DW.get(0) = fac * dx.get(0);
197 DW.get(1) = fac * dx.get(1);
198 DW.get(2) = fac * dx.get(2);
199 }
200 else
201 {
202 DW.get(0) = 0.0;
203 DW.get(1) = 0.0;
204 DW.get(2) = 0.0;
205 }
206}
207
208
209// Tensile correction
210inline double Tensile(double r, double rhoa, double rhob, double prs1, double prs2)
211{
212 const double qq=r/H;
213 //-Cubic Spline kernel
214 double wab;
215 if(r>H)
216 {
217 double wqq1=2.0f-qq;
218 double wqq2=wqq1*wqq1;
219
220 wab=a2_4*(wqq2*wqq1);
221 }
222 else
223 {
224 double wqq2=qq*qq;
225 double wqq3=wqq2*qq;
226
227 wab=a2*(1.0f-1.5f*wqq2+0.75f*wqq3);
228 }
229
230 //-Tensile correction.
231 double fab=wab*W_dap;
232 fab*=fab; fab*=fab; //fab=fab^4
233 const double tensilp1=(prs1/(rhoa*rhoa))*(prs1>0? 0.01: -0.2);
234 const double tensilp2=(prs2/(rhob*rhob))*(prs2>0? 0.01: -0.2);
235
236 return (fab*(tensilp1+tensilp2));
237}
238
239
240inline double Pi(const Point<3,double> & dr, double rr2, Point<3,double> & dv, double rhoa, double rhob, double massb, double & visc)
241{
242 const double dot = dr.get(0)*dv.get(0) + dr.get(1)*dv.get(1) + dr.get(2)*dv.get(2);
243 const double dot_rr2 = dot/(rr2+Eta2);
244 visc=std::max(dot_rr2,visc);
245
246 if(dot < 0)
247 {
248 const float amubar=H*dot_rr2;
249 const float robar=(rhoa+rhob)*0.5f;
250 const float pi_visc=(-visco*cbar*amubar/robar);
251
252 return pi_visc;
253 }
254 else
255 return 0.0;
256}
257
258
259template<typename VerletList> inline double calc_forces(particles & vd, VerletList & NN, double & max_visc, double r_gskin)
260{
261 auto NN2 = vd.getCellList(r_gskin);
262
263 // Reset the ghost
264 auto itg = vd.getDomainAndGhostIterator();
265 while (itg.isNext())
266 {
267 auto p = itg.get();
268 // Reset force
269
270 if (p.getKey() < vd.size_local())
271 {
272 // Reset the force counter (- gravity on zeta direction)
273 vd.template getProp<force>(p)[0] = 0.0;
274 vd.template getProp<force>(p)[1] = 0.0;
275 vd.template getProp<force>(p)[2] = -gravity;
276 vd.template getProp<drho>(p) = 0.0;
277
278
279 // Reset the force counter (- gravity on zeta direction)
280 vd.template getProp<10>(p)[0] = 0.0;
281 vd.template getProp<10>(p)[1] = 0.0;
282 vd.template getProp<10>(p)[2] = -gravity;
283 vd.template getProp<11>(p) = 0.0;
284 }
285 else
286 {
287 vd.getProp<force>(p)[0] = 0.0;
288 vd.getProp<force>(p)[1] = 0.0;
289 vd.getProp<force>(p)[2] = 0.0;
290
291 vd.getProp<drho>(p) = 0.0;
292 }
293
294
295 vd.getProp<8>(p).clear();
296 vd.getProp<9>(p).clear();
297
298 ++itg;
299 }
300
301 // Get an iterator over particles
302 auto part = vd.getParticleIteratorCRS(NN.getInternalCellList());
303
304 double visc = 0;
305
306 // For each particle ...
307 while (part.isNext())
308 {
309 // ... a
310 auto a = part.get();
311
312 // Get the position xp of the particle
313 Point<3,double> xa = vd.getPos(a);
314
315 // Take the mass of the particle dependently if it is FLUID or BOUNDARY
316 double massa = (vd.getProp<type>(a) == FLUID)?MassFluid:MassBound;
317
318 // Get the density of the of the particle a
319 double rhoa = vd.getProp<rho>(a);
320
321 // Get the pressure of the particle a
322 double Pa = vd.getProp<Pressure>(a);
323
324 // Get the Velocity of the particle a
325 Point<3,double> va = vd.getProp<velocity>(a);
326
327 // We threat FLUID particle differently from BOUNDARY PARTICLES ...
328 if (vd.getProp<type>(a) != FLUID)
329 {
330 // If it is a boundary particle calculate the delta rho based on equation 2
331 // This require to run across the neighborhoods particles of a
332 auto Np = NN.template getNNIterator<NO_CHECK>(a);
333
334 // For each neighborhood particle
335 while (Np.isNext() == true)
336 {
337 // ... q
338 auto b = Np.get();
339
340 // Get the position xp of the particle
341 Point<3,double> xb = vd.getPos(b);
342
343 // if (p == q) skip this particle
344 if (a == b) {++Np; continue;};
345
346 // get the mass of the particle
347 double massb = (vd.getProp<type>(b) == FLUID)?MassFluid:MassBound;
348
349 // Get the velocity of the particle b
350 Point<3,double> vb = vd.getProp<velocity>(b);
351
352 // Get the pressure and density of particle b
353 double Pb = vd.getProp<Pressure>(b);
354 double rhob = vd.getProp<rho>(b);
355
356 // Get the distance between p and q
357 Point<3,double> dr = xa - xb;
358 // take the norm of this vector
359 double r2 = norm2(dr);
360
361 // If the particles interact ...
362 if (r2 < 4.0*H*H)
363 {
364 // ... calculate delta rho
365 double r = sqrt(r2);
366
367 Point<3,double> dv = va - vb;
368
370 DWab(dr,DW,r,false);
371
372 const double dot = dr.get(0)*dv.get(0) + dr.get(1)*dv.get(1) + dr.get(2)*dv.get(2);
373 const double dot_rr2 = dot/(r2+Eta2);
374 max_visc=std::max(dot_rr2,max_visc);
375
376 double scal = (dv.get(0)*DW.get(0)+dv.get(1)*DW.get(1)+dv.get(2)*DW.get(2));
377
378 vd.getProp<drho>(a) += massb*scal;
379 vd.getProp<drho>(b) += massa*scal;
380
381 struct nb_f nna;
382 nna.id = b;
383 nna.fact = 0.0;
384
385 struct nb_f nnb;
386
387 // If it is a fluid we have to calculate force
388 if (vd.getProp<type>(b) == FLUID)
389 {
390 Point<3,double> v_rel = va - vb;
391 double factor = - ((vd.getProp<Pressure>(a) + vd.getProp<Pressure>(b)) / (rhoa * rhob) + Tensile(r,rhoa,rhob,Pa,Pb) + Pi(dr,r2,v_rel,rhoa,rhob,massb,visc));
392
393 vd.getProp<force>(b)[0] -= massa * factor * DW.get(0);
394 vd.getProp<force>(b)[1] -= massa * factor * DW.get(1);
395 vd.getProp<force>(b)[2] -= massa * factor * DW.get(2);
396
397 nnb.id = vd.getProp<12>(a);
398 nnb.fact = -massa * factor * DW.get(2);
399 }
400 else
401 {
402 struct nb_f nnb;
403 nnb.id = vd.getProp<12>(a);
404 nnb.fact = 0.0;
405 }
406
407 vd.getProp<9>(a).add(nna);
408 vd.getProp<9>(b).add(nnb);
409
410 }
411
412 ++Np;
413 }
414
416
417 // Get an iterator over the neighborhood particles of p
418 auto Np2 = NN2.template getNNIterator<NO_CHECK>(NN2.getCell(vd.getPos(a)));
419
420 // For each neighborhood particle
421 while (Np2.isNext() == true)
422 {
423 // ... q
424 auto b = Np2.get();
425
426 // Get the position xp of the particle
427 Point<3,double> xb = vd.getPos(b);
428
429 // if (p == q) skip this particle
430 if (a == b) {++Np2; continue;};
431
432 // get the mass of the particle
433 double massb = (vd.getProp<type>(b) == FLUID)?MassFluid:MassBound;
434
435 // Get the velocity of the particle b
436 Point<3,double> vb = vd.getProp<velocity>(b);
437
438 // Get the pressure and density of particle b
439 double Pb = vd.getProp<Pressure>(b);
440 double rhob = vd.getProp<rho>(b);
441
442 // Get the distance between p and q
443 Point<3,double> dr = xa - xb;
444 // take the norm of this vector
445 double r2 = norm2(dr);
446
447 // If the particles interact ...
448 if (r2 < 4.0*H*H)
449 {
450 // ... calculate delta rho
451 double r = sqrt(r2);
452
453 Point<3,double> dv = va - vb;
454
456 DWab(dr,DW,r,false);
457
458 const double dot = dr.get(0)*dv.get(0) + dr.get(1)*dv.get(1) + dr.get(2)*dv.get(2);
459 const double dot_rr2 = dot/(r2+Eta2);
460
461 vd.getProp<11>(a) += massb*(dv.get(0)*DW.get(0)+dv.get(1)*DW.get(1)+dv.get(2)*DW.get(2));
462
463
464 struct nb_f nna;
465 nna.id = vd.getProp<12>(b);
466 nna.fact = 0.0;
467
468 vd.getProp<8>(a).add(nna);
469 }
470
471 ++Np2;
472 }
473
474 }
475 else
476 {
477 // If it is a fluid particle calculate based on equation 1 and 2
478
479 // Get an iterator over the neighborhood particles of p
480 auto Np = NN.template getNNIterator<NO_CHECK>(a);
481
482 // For each neighborhood particle
483 while (Np.isNext() == true)
484 {
485 // ... q
486 auto b = Np.get();
487
488 // Get the position xp of the particle
489 Point<3,double> xb = vd.getPos(b);
490
491 // if (p == q) skip this particle
492 if (a == b) {++Np; continue;};
493
494 double massb = (vd.getProp<type>(b) == FLUID)?MassFluid:MassBound;
495 Point<3,double> vb = vd.getProp<velocity>(b);
496 double Pb = vd.getProp<Pressure>(b);
497 double rhob = vd.getProp<rho>(b);
498
499 // Get the distance between p and q
500 Point<3,double> dr = xa - xb;
501 // take the norm of this vector
502 double r2 = norm2(dr);
503
504 // if they interact
505 if (r2 < 4.0*H*H)
506 {
507 double r = sqrt(r2);
508
509 Point<3,double> v_rel = va - vb;
510
512 DWab(dr,DW,r,false);
513
514 double factor = - ((vd.getProp<Pressure>(a) + vd.getProp<Pressure>(b)) / (rhoa * rhob) + Tensile(r,rhoa,rhob,Pa,Pb) + Pi(dr,r2,v_rel,rhoa,rhob,massb,visc));
515
516 vd.getProp<force>(a)[0] += massb * factor * DW.get(0);
517 vd.getProp<force>(a)[1] += massb * factor * DW.get(1);
518 vd.getProp<force>(a)[2] += massb * factor * DW.get(2);
519
520 vd.getProp<force>(b)[0] -= massa * factor * DW.get(0);
521 vd.getProp<force>(b)[1] -= massa * factor * DW.get(1);
522 vd.getProp<force>(b)[2] -= massa * factor * DW.get(2);
523
524 double scal = (v_rel.get(0)*DW.get(0)+v_rel.get(1)*DW.get(1)+v_rel.get(2)*DW.get(2));
525
526 vd.getProp<drho>(a) += massb*scal;
527 vd.getProp<drho>(b) += massa*scal;
528
529
530 if (vd.getProp<12>(a) == 15604 && vd.getProp<12>(b) == 15601)
531 {
532 std::cout << "DETECTED ERROR " << std::endl;
533 }
534
535 if (vd.getProp<12>(b) == 15604 && vd.getProp<12>(a) == 15601)
536 {
537 std::cout << "DETECTED ERROR " << create_vcluster().getProcessUnitID() << " a " << a << " b " << b << " " << vd.size_local_with_ghost() << std::endl;
538 }
539
540
541 struct nb_f nna;
542 nna.id = vd.getProp<12>(b);
543 nna.fact = massb * factor * DW.get(2);
544
545 struct nb_f nnb;
546 nnb.id = vd.getProp<12>(a);
547 nnb.fact = -massa * factor * DW.get(2);
548
549 vd.getProp<9>(a).add(nna);
550 vd.getProp<9>(b).add(nnb);
551 }
552
553 ++Np;
554 }
555
556
558
559
561
562 // Get an iterator over the neighborhood particles of p
563 auto Np2 = NN2.template getNNIterator<NO_CHECK>(NN2.getCell(vd.getPos(a)));
564
565 if (a == 0 && create_vcluster().getProcessUnitID() == 0)
566 {
567// std::cout << "Particle 0" << std::endl;
568 }
569
570 // For each neighborhood particle
571 while (Np2.isNext() == true)
572 {
573 // ... q
574 auto b = Np2.get();
575
576 // Get the position xp of the particle
577 Point<3,double> xb = vd.getPos(b);
578
579 // if (p == q) skip this particle
580 if (a == b) {++Np2; continue;};
581
582 double massb = (vd.getProp<type>(b) == FLUID)?MassFluid:MassBound;
583 Point<3,double> vb = vd.getProp<velocity>(b);
584 double Pb = vd.getProp<Pressure>(b);
585 double rhob = vd.getProp<rho>(b);
586
587 // Get the distance between p and q
588 Point<3,double> dr = xa - xb;
589 // take the norm of this vector
590 double r2 = norm2(dr);
591
592 // if they interact
593 if (r2 < 4.0*H*H)
594 {
595 double r = sqrt(r2);
596
597 Point<3,double> v_rel = va - vb;
598
600 DWab(dr,DW,r,false);
601
602 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,visc));
603
604 vd.getProp<10>(a)[0] += factor * DW.get(0);
605 vd.getProp<10>(a)[1] += factor * DW.get(1);
606 vd.getProp<10>(a)[2] += factor * DW.get(2);
607
608 vd.getProp<11>(a) += massb*(v_rel.get(0)*DW.get(0)+v_rel.get(1)*DW.get(1)+v_rel.get(2)*DW.get(2));
609
610 struct nb_f nna;
611 nna.id = vd.getProp<12>(b);
612 nna.fact = factor * DW.get(2);
613
614 vd.getProp<8>(a).add(nna);
615 }
616
617 ++Np2;
618 }
619 }
620
621 ++part;
622 }
623
624 vd.template ghost_put<add_,drho>();
625 vd.template ghost_put<add_,force>();
626 vd.template ghost_put<merge_,9>();
627 vd.template ghost_put<add_,10>();
628 vd.template ghost_put<add_,11>();
629
630 auto part3 = vd.getDomainIterator();
631
632 while (part3.isNext())
633 {
634 auto a = part3.get();
635
636 if (vd.getProp<force>(a)[0] != vd.getProp<10>(a)[0] ||
637 vd.getProp<force>(a)[1] != vd.getProp<10>(a)[1] ||
638 vd.getProp<force>(a)[2] != vd.getProp<10>(a)[2])
639 {
640 if (create_vcluster().getProcessUnitID() == 0)
641 {
642 std::cout << "LISTS " << vd.getProp<12>(a) << " " << vd.getProp<8>(a).size() << " " << vd.getProp<9>(a).size() << std::endl;
643 std::cout << "ERROR: " << vd.getProp<force>(a)[2] << " " << vd.getProp<10>(a)[2] << " part: " << a.getKey() << std::endl;
644
645 // Print factors and sum
646
647
648 vd.getProp<8>(a).sort();
649 vd.getProp<9>(a).sort();
650
651 double sum1 = 0.0;
652 double sum2 = 0.0;
653
654 for (size_t i = 0 ; i < vd.getProp<8>(a).size() ; i++)
655 {
656
657 std::cout << "FACT: " << vd.getProp<8>(a).get(i).fact << " " << vd.getProp<9>(a).get(i).fact << " " << vd.getProp<8>(a).get(i).id << " " << vd.getProp<9>(a).get(i).id << std::endl;
658
659 sum1 += vd.getProp<8>(a).get(i).fact;
660 sum2 += vd.getProp<9>(a).get(i).fact;
661 }
662
663 std::cout << "sum: " << sum1 << " " << sum2 << std::endl;
664 }
665
666 break;
667 }
668
669 vd.getProp<8>(a).sort();
670 vd.getProp<9>(a).sort();
671
672 if (create_vcluster().getProcessUnitID() == 0)
673 {
674 if (vd.getProp<8>(a).size() != vd.getProp<9>(a).size())
675 {
676 std::cout << "ERROR: " << vd.getProp<12>(a) << " " << vd.getProp<8>(a).size() << " " << vd.getProp<9>(a).size() << std::endl;
677
678 for (size_t i = 0 ; i < vd.getProp<8>(a).size() ; i++)
679 {
680 std::cout << "NNNN: " << vd.getProp<9>(a).get(i).id << " " << vd.getProp<8>(a).get(i).id << std::endl;
681 }
682
683 break;
684 }
685 }
686
687 ++part3;
688 }
689}
690
691
692void max_acceleration_and_velocity(particles & vd, double & max_acc, double & max_vel)
693{
694 // Calculate the maximum acceleration
695 auto part = vd.getDomainIterator();
696
697 while (part.isNext())
698 {
699 auto a = part.get();
700
701 Point<3,double> acc(vd.getProp<force>(a));
702 double acc2 = norm2(acc);
703
704 Point<3,double> vel(vd.getProp<velocity>(a));
705 double vel2 = norm2(vel);
706
707 if (vel2 >= max_vel)
708 max_vel = vel2;
709
710 if (acc2 >= max_acc)
711 max_acc = acc2;
712
713 ++part;
714 }
715 max_acc = sqrt(max_acc);
716 max_vel = sqrt(max_vel);
717
718 Vcluster & v_cl = create_vcluster();
719 v_cl.max(max_acc);
720 v_cl.max(max_vel);
721 v_cl.execute();
722}
723
724
725double calc_deltaT(particles & vd, double ViscDtMax)
726{
727 double Maxacc = 0.0;
728 double Maxvel = 0.0;
729 max_acceleration_and_velocity(vd,Maxacc,Maxvel);
730
731 //-dt1 depends on force per unit mass.
732 const double dt_f = (Maxacc)?sqrt(H/Maxacc):std::numeric_limits<int>::max();
733
734 //-dt2 combines the Courant and the viscous time-step controls.
735 const double dt_cv = H/(std::max(cbar,Maxvel*10.) + H*ViscDtMax);
736
737 //-dt new value of time step.
738 double dt=double(CFLnumber)*std::min(dt_f,dt_cv);
739 if(dt<double(DtMin))
740 dt=double(DtMin);
741
742 return dt;
743}
744
745
747
748
749void verlet_int(particles & vd, double dt, double & max_disp)
750{
751 // list of the particle to remove
752 to_remove.clear();
753
754 // particle iterator
755 auto part = vd.getDomainIterator();
756
757 double dt205 = dt*dt*0.5;
758 double dt2 = dt*2.0;
759
760 max_disp = 0;
761
762 // For each particle ...
763 while (part.isNext())
764 {
765 // ... a
766 auto a = part.get();
767
768 // if the particle is boundary
769 if (vd.template getProp<type>(a) == BOUNDARY)
770 {
771 // Update rho
772 double rhop = vd.template getProp<rho>(a);
773
774 // Update only the density
775 vd.template getProp<velocity>(a)[0] = 0.0;
776 vd.template getProp<velocity>(a)[1] = 0.0;
777 vd.template getProp<velocity>(a)[2] = 0.0;
778 vd.template getProp<rho>(a) = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
779
780 vd.template getProp<rho_prev>(a) = rhop;
781
782 ++part;
783 continue;
784 }
785
786 //-Calculate displacement and update position / Calcula desplazamiento y actualiza posicion.
787 double dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
788 double dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
789 double dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
790
791 vd.getPos(a)[0] += dx;
792 vd.getPos(a)[1] += dy;
793 vd.getPos(a)[2] += dz;
794
795 double d2 = dx*dx + dy*dy + dz*dz;
796
797 max_disp = (max_disp > d2)?max_disp:d2;
798
799 double velX = vd.template getProp<velocity>(a)[0];
800 double velY = vd.template getProp<velocity>(a)[1];
801 double velZ = vd.template getProp<velocity>(a)[2];
802 double rhop = vd.template getProp<rho>(a);
803
804 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity_prev>(a)[0] + vd.template getProp<force>(a)[0]*dt2;
805 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity_prev>(a)[1] + vd.template getProp<force>(a)[1]*dt2;
806 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity_prev>(a)[2] + vd.template getProp<force>(a)[2]*dt2;
807 vd.template getProp<rho>(a) = vd.template getProp<rho_prev>(a) + dt2*vd.template getProp<drho>(a);
808
809 // Check if the particle go out of range in space and in density
810 if (vd.getPos(a)[0] < 0.000263878 || vd.getPos(a)[1] < 0.000263878 || vd.getPos(a)[2] < 0.000263878 ||
811 vd.getPos(a)[0] > 0.000263878+1.59947 || vd.getPos(a)[1] > 0.000263878+0.672972 || vd.getPos(a)[2] > 0.000263878+0.903944 ||
812 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
813 {
814 to_remove.add(a.getKey());
815
816
817 // if we remove something the verlet are not anymore correct and must be reconstructed
818 max_disp = 100.0;
819 }
820
821 vd.template getProp<velocity_prev>(a)[0] = velX;
822 vd.template getProp<velocity_prev>(a)[1] = velY;
823 vd.template getProp<velocity_prev>(a)[2] = velZ;
824 vd.template getProp<rho_prev>(a) = rhop;
825
826 ++part;
827 }
828
829 Vcluster & v_cl = create_vcluster();
830 v_cl.max(max_disp);
831 v_cl.execute();
832
833 max_disp = sqrt(max_disp);
834
835 // remove the particles
836 vd.remove(to_remove,0);
837
838 // increment the iteration counter
839 cnt++;
840}
841
842void euler_int(particles & vd, double dt, double & max_disp)
843{
844 // list of the particle to remove
845 to_remove.clear();
846
847 // particle iterator
848 auto part = vd.getDomainIterator();
849
850 double dt205 = dt*dt*0.5;
851 double dt2 = dt*2.0;
852
853 max_disp = 0;
854
855 // For each particle ...
856 while (part.isNext())
857 {
858 // ... a
859 auto a = part.get();
860
861 // if the particle is boundary
862 if (vd.template getProp<type>(a) == BOUNDARY)
863 {
864 // Update rho
865 double rhop = vd.template getProp<rho>(a);
866
867 // Update only the density
868 vd.template getProp<velocity>(a)[0] = 0.0;
869 vd.template getProp<velocity>(a)[1] = 0.0;
870 vd.template getProp<velocity>(a)[2] = 0.0;
871 vd.template getProp<rho>(a) = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
872
873 vd.template getProp<rho_prev>(a) = rhop;
874
875 ++part;
876 continue;
877 }
878
879 //-Calculate displacement and update position / Calcula desplazamiento y actualiza posicion.
880 double dx = vd.template getProp<velocity>(a)[0]*dt + vd.template getProp<force>(a)[0]*dt205;
881 double dy = vd.template getProp<velocity>(a)[1]*dt + vd.template getProp<force>(a)[1]*dt205;
882 double dz = vd.template getProp<velocity>(a)[2]*dt + vd.template getProp<force>(a)[2]*dt205;
883
884 vd.getPos(a)[0] += dx;
885 vd.getPos(a)[1] += dy;
886 vd.getPos(a)[2] += dz;
887
888 double d2 = dx*dx + dy*dy + dz*dz;
889
890 max_disp = (max_disp > d2)?max_disp:d2;
891
892 double velX = vd.template getProp<velocity>(a)[0];
893 double velY = vd.template getProp<velocity>(a)[1];
894 double velZ = vd.template getProp<velocity>(a)[2];
895 double rhop = vd.template getProp<rho>(a);
896
897 vd.template getProp<velocity>(a)[0] = vd.template getProp<velocity>(a)[0] + vd.template getProp<force>(a)[0]*dt;
898 vd.template getProp<velocity>(a)[1] = vd.template getProp<velocity>(a)[1] + vd.template getProp<force>(a)[1]*dt;
899 vd.template getProp<velocity>(a)[2] = vd.template getProp<velocity>(a)[2] + vd.template getProp<force>(a)[2]*dt;
900 vd.template getProp<rho>(a) = vd.template getProp<rho>(a) + dt*vd.template getProp<drho>(a);
901
902 // Check if the particle go out of range in space and in density
903 if (vd.getPos(a)[0] < 0.000263878 || vd.getPos(a)[1] < 0.000263878 || vd.getPos(a)[2] < 0.000263878 ||
904 vd.getPos(a)[0] > 0.000263878+1.59947 || vd.getPos(a)[1] > 0.000263878+0.672972 || vd.getPos(a)[2] > 0.000263878+0.903944 ||
905 vd.template getProp<rho>(a) < RhoMin || vd.template getProp<rho>(a) > RhoMax)
906 {
907 to_remove.add(a.getKey());
908
909 // if we remove something the verlet are not anymore correct and must be reconstructed
910 max_disp = 100.0;
911 }
912
913 vd.template getProp<velocity_prev>(a)[0] = velX;
914 vd.template getProp<velocity_prev>(a)[1] = velY;
915 vd.template getProp<velocity_prev>(a)[2] = velZ;
916 vd.template getProp<rho_prev>(a) = rhop;
917
918 ++part;
919 }
920
921 // remove the particles
922 vd.remove(to_remove,0);
923
924 Vcluster & v_cl = create_vcluster();
925 v_cl.max(max_disp);
926 v_cl.execute();
927
928 max_disp = sqrt(max_disp);
929
930 // increment the iteration counter
931 cnt++;
932}
933
934
935int main(int argc, char* argv[])
936{
937 // initialize the library
938 openfpm_init(&argc,&argv);
939
940 // Here we define our domain a 2D box with internals from 0 to 1.0 for x and y
941 Box<3,double> domain({-0.05,-0.05,-0.05},{1.7010,0.7065,0.5025});
942 size_t sz[3] = {207,90,66};
943
944 // Fill W_dap
945 W_dap = 1.0/Wab(H/1.5);
946
947 // Here we define the boundary conditions of our problem
948 size_t bc[3]={NON_PERIODIC,NON_PERIODIC,NON_PERIODIC};
949
950 double skin = 0.25 * 2*H;
951 double r_gskin = 2*H + skin;
952
953 // extended boundary around the domain, and the processor domain
954 // by the support of the cubic kernel
955 Ghost<3,double> g(r_gskin);
956
957 particles vd(0,domain,bc,g,BIND_DEC_TO_GHOST);
958
959 // You can ignore all these dp/2.0 is a trick to reach the same initialization
960 // of Dual-SPH that use a different criteria to draw particles
961 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});
962
963 // return an iterator to the fluid particles to add to vd
964 auto fluid_it = DrawParticles::DrawBox(vd,sz,domain,fluid_box);
965
966 // here we fill some of the constants needed by the simulation
967 max_fluid_height = fluid_it.getBoxMargins().getHigh(2);
968 h_swl = fluid_it.getBoxMargins().getHigh(2) - fluid_it.getBoxMargins().getLow(2);
969 B = (coeff_sound)*(coeff_sound)*gravity*h_swl*rho_zero / gamma_;
970 cbar = coeff_sound * sqrt(gravity * h_swl);
971
972 // for each particle inside the fluid box ...
973 while (fluid_it.isNext())
974 {
975 // ... add a particle ...
976 vd.add();
977
978 // ... and set it position ...
979 vd.getLastPos()[0] = fluid_it.get().get(0);
980 vd.getLastPos()[1] = fluid_it.get().get(1);
981 vd.getLastPos()[2] = fluid_it.get().get(2);
982
983 // and its type.
984 vd.template getLastProp<type>() = FLUID;
985
986 // We also initialize the density of the particle and the hydro-static pressure given by
987 //
988 // rho_zero*g*h = P
989 //
990 // rho_p = (P/B + 1)^(1/Gamma) * rho_zero
991 //
992
993 vd.template getLastProp<Pressure>() = rho_zero * gravity * (max_fluid_height - fluid_it.get().get(2));
994
995 vd.template getLastProp<rho>() = pow(vd.template getLastProp<Pressure>() / B + 1, 1.0/gamma_) * rho_zero;
996 vd.template getLastProp<rho_prev>() = vd.template getLastProp<rho>();
997 vd.template getLastProp<velocity>()[0] = 0.0;
998 vd.template getLastProp<velocity>()[1] = 0.0;
999 vd.template getLastProp<velocity>()[2] = 0.0;
1000
1001 vd.template getLastProp<velocity_prev>()[0] = 0.0;
1002 vd.template getLastProp<velocity_prev>()[1] = 0.0;
1003 vd.template getLastProp<velocity_prev>()[2] = 0.0;
1004
1005 // next fluid particle
1006 ++fluid_it;
1007 }
1008
1009
1010 // Recipient
1011 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});
1012 Box<3,double> recipient2({dp,dp,dp},{1.6-dp/2.0,0.67-dp/2.0,0.4+dp/2.0});
1013
1014 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});
1015 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});
1016 Box<3,double> obstacle3({0.9+dp,0.24,0.0},{1.02,0.36,0.45});
1017
1019 holes.add(recipient2);
1020 holes.add(obstacle1);
1021 auto bound_box = DrawParticles::DrawSkin(vd,sz,domain,holes,recipient1);
1022
1023 while (bound_box.isNext())
1024 {
1025 vd.add();
1026
1027 vd.getLastPos()[0] = bound_box.get().get(0);
1028 vd.getLastPos()[1] = bound_box.get().get(1);
1029 vd.getLastPos()[2] = bound_box.get().get(2);
1030
1031 vd.template getLastProp<type>() = BOUNDARY;
1032 vd.template getLastProp<rho>() = rho_zero;
1033 vd.template getLastProp<rho_prev>() = rho_zero;
1034 vd.template getLastProp<velocity>()[0] = 0.0;
1035 vd.template getLastProp<velocity>()[1] = 0.0;
1036 vd.template getLastProp<velocity>()[2] = 0.0;
1037
1038 vd.template getLastProp<velocity_prev>()[0] = 0.0;
1039 vd.template getLastProp<velocity_prev>()[1] = 0.0;
1040 vd.template getLastProp<velocity_prev>()[2] = 0.0;
1041
1042 ++bound_box;
1043 }
1044
1045 auto obstacle_box = DrawParticles::DrawSkin(vd,sz,domain,obstacle2,obstacle1);
1046
1047 while (obstacle_box.isNext())
1048 {
1049 vd.add();
1050
1051 vd.getLastPos()[0] = obstacle_box.get().get(0);
1052 vd.getLastPos()[1] = obstacle_box.get().get(1);
1053 vd.getLastPos()[2] = obstacle_box.get().get(2);
1054
1055 vd.template getLastProp<type>() = BOUNDARY;
1056 vd.template getLastProp<rho>() = rho_zero;
1057 vd.template getLastProp<rho_prev>() = rho_zero;
1058 vd.template getLastProp<velocity>()[0] = 0.0;
1059 vd.template getLastProp<velocity>()[1] = 0.0;
1060 vd.template getLastProp<velocity>()[2] = 0.0;
1061
1062 vd.template getLastProp<velocity_prev>()[0] = 0.0;
1063 vd.template getLastProp<velocity_prev>()[1] = 0.0;
1064 vd.template getLastProp<velocity_prev>()[2] = 0.0;
1065
1066 ++obstacle_box;
1067 }
1068
1069 vd.map();
1070
1071 // Now that we fill the vector with particles
1072 ModelCustom md;
1073
1074 vd.addComputationCosts(md);
1075 vd.getDecomposition().decompose();
1076 vd.map();
1077
1078 size_t loc = vd.size_local();
1080
1081 Vcluster & v_cl = create_vcluster();
1082 v_cl.allGather(loc,gt);
1083 v_cl.execute();
1084
1085
1086 auto debug_it = vd.getDomainIterator();
1087
1088
1089 size_t tot = 0;
1090
1091
1092 for (size_t i ; i < create_vcluster().getProcessUnitID() ; i++)
1093 tot += gt.get(i);
1094
1095 while(debug_it.isNext())
1096 {
1097 auto a = debug_it.get();
1098
1099 vd.getProp<12>(a) = tot;
1100 tot++;
1101
1102 ++debug_it;
1103 }
1104
1105 vd.ghost_get<type,rho,Pressure,velocity,12>();
1106
1107 auto NN = vd.getVerletCrs(r_gskin);
1108
1109 // Evolve
1110
1111 size_t write = 0;
1112 size_t it = 0;
1113 size_t it_reb = 0;
1114 double t = 0.0;
1115 double tot_disp = 0.0;
1116 double max_disp;
1117 while (t <= t_end)
1118 {
1119 Vcluster & v_cl = create_vcluster();
1120 timer it_time;
1121
1122 it_reb++;
1123// if (2*tot_disp >= skin)
1124// {
1125 vd.map();
1126
1127 vd.getDecomposition().write("DecBEFORE");
1128
1129 if (it_reb > 5)
1130 {
1131 ModelCustom md;
1132 vd.addComputationCosts(md);
1133 vd.getDecomposition().decompose();
1134
1135 vd.map();
1136
1137 it_reb = 0;
1138
1139 if (v_cl.getProcessUnitID() == 0)
1140 std::cout << "REBALANCED " << std::endl;
1141
1142 }
1143
1144 vd.getDecomposition().write("DecAFTER");
1145
1146 // Calculate pressure from the density
1147 EqState(vd);
1148
1149 vd.ghost_get<type,rho,Pressure,velocity,12>();
1150
1151 //vd.updateVerlet(NN,r_gskin,VL_CRS_SYMMETRIC);
1152 NN = vd.getVerletCrs(r_gskin);
1153
1154 vd.write("Debug");
1155
1156 tot_disp = 0.0;
1157
1158 if (v_cl.getProcessUnitID() == 0)
1159 std::cout << "RECONSTRUCT Verlet " << std::endl;
1160// }
1161// else
1162// {
1163// // Calculate pressure from the density
1164// EqState(vd);
1165//
1166// vd.ghost_get<type,rho,Pressure,velocity>(SKIP_LABELLING);
1167// }
1168
1169 double max_visc = 0.0;
1170
1171 // Calc forces
1172 calc_forces(vd,NN,max_visc,r_gskin);
1173
1174 // Get the maximum viscosity term across processors
1175 v_cl.max(max_visc);
1176 v_cl.execute();
1177
1178 // Calculate delta t integration
1179 double dt = calc_deltaT(vd,max_visc);
1180
1181 // VerletStep or euler step
1182 it++;
1183 if (it < 40)
1184 verlet_int(vd,dt,max_disp);
1185 else
1186 {
1187 euler_int(vd,dt,max_disp);
1188 it = 0;
1189 }
1190
1191 tot_disp += max_disp;
1192 t += dt;
1193
1194 if (write < t*100)
1195 {
1196// vd.deleteGhost();
1197 vd.write("Geometry",write);
1198// vd.ghost_get<type,rho,Pressure,velocity>(SKIP_LABELLING);
1199 write++;
1200
1201 if (v_cl.getProcessUnitID() == 0)
1202 std::cout << "TIME: " << t << " write " << it_time.getwct() << " " << v_cl.getProcessUnitID() << " TOT disp: " << tot_disp << " " << cnt << std::endl;
1203 }
1204 else
1205 {
1206 if (v_cl.getProcessUnitID() == 0)
1207 std::cout << "TIME: " << t << " " << it_time.getwct() << " " << v_cl.getProcessUnitID() << " TOT disp: " << tot_disp << " " << cnt << std::endl;
1208 }
1209 }
1210
1211 openfpm_finalize();
1212}
1213
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
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.
size_t getProcessUnitID()
Get the process unit id.
bool allGather(T &send, openfpm::vector< T, Mem, gr > &v)
Gather the data from all processors.
void max(T &num)
Get the maximum number across all processors (or reduction with infinity norm)
Implementation of VCluster class.
Definition VCluster.hpp:59
Class for Verlet list implementation.
CellListImpl & getInternalCellList()
Get the internal cell-list used to construct the Verlet-list.
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
vect_dist_key_dx get()
Get the actual key.
Distributed vector.
size_t size_local() const
return the local size of the 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.
size_t size_local_with_ghost() const
return the local size of the vector
openfpm::vector_key_iterator_seq< typename vrl::Mem_type_type::local_index_type > getParticleIteratorCRS(vrl &NN)
Get a special particle iterator able to iterate across particles using symmetric crossing scheme.
auto getPos(vect_dist_key_dx vec_key) -> decltype(v_pos.template get< 0 >(vec_key.getKey()))
Get the position of an element.
CellL getCellList(St r_cut, bool no_se3=false)
Construct a cell list starting from the stored particles.
vector_dist_iterator getDomainIterator() const
Get an iterator that traverse the particles in the domain.
void ghost_get(size_t opt=WITH_POSITION)
It synchronize the properties and position of the ghost particles.
void map(size_t opt=NONE)
It move all the particles that does not belong to the local processor to the respective processor.
bool write(std::string out, int opt=VTK_WRITER)
Output particle position and properties.
auto getLastPos() -> decltype(v_pos.template get< 0 >(0))
Get the position of the last element.
VerletL getVerletCrs(St r_cut)
for each particle get the symmetric verlet list
vector_dist_iterator getDomainAndGhostIterator() const
Get an iterator that traverse the particles in the domain.
void add()
Add local particle.
void addComputationCosts(const self &vd, Model md=Model())
Add the computation cost on the decomposition coming from the particles.
Decomposition & getDecomposition()
Get the decomposition.
Model for Dynamic load balancing.
Definition main.cpp:222