KinematicSurfaceFilm.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2021-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "KinematicSurfaceFilm.H"
29 #include "surfaceFilmRegionModel.H"
30 #include "liquidFilmModel.H"
32 #include "unitConversion.H"
33 #include "Pstream.H"
34 
35 using namespace Foam::constant::mathematical;
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 template<class CloudType>
40 const Foam::Enum
41 <
43 >
45 ({
46  { interactionType::absorb, "absorb" },
47  { interactionType::bounce, "bounce" },
48  { interactionType::splashBai, "splashBai" },
49 });
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
54 template<class CloudType>
56 (
57  const vector& v
58 ) const
59 {
60  vector tangent(Zero);
61  scalar magTangent = 0.0;
62 
63  while (magTangent < SMALL)
64  {
65  const vector vTest(rndGen_.sample01<vector>());
66  tangent = vTest - (vTest & v)*v;
67  magTangent = mag(tangent);
68  }
69 
70  return tangent/magTangent;
71 }
72 
73 
74 template<class CloudType>
76 (
77  const vector& tanVec1,
78  const vector& tanVec2,
79  const vector& nf
80 ) const
81 {
82  // Azimuthal angle [rad]
83  const scalar phiSi = twoPi*rndGen_.sample01<scalar>();
84 
85  // Ejection angle [rad]
86  const scalar thetaSi = degToRad(rndGen_.sample01<scalar>()*(50 - 5) + 5);
87 
88  // Direction vector of new parcel
89  const scalar alpha = sin(thetaSi);
90  const scalar dcorr = cos(thetaSi);
91  const vector normal(alpha*(tanVec1*cos(phiSi) + tanVec2*sin(phiSi)));
92  vector dirVec(dcorr*nf);
93  dirVec += normal;
94 
95  return dirVec/mag(dirVec);
96 }
97 
98 
99 template<class CloudType>
101 {
102  const fvMesh& mesh = this->owner().mesh();
103 
104  // Set up region film
105  if (!filmModel_)
106  {
107  filmModel_ =
108  mesh.time().objectRegistry::template getObjectPtr<regionFilm>
109  (
110  "surfaceFilmProperties"
111  );
112  }
113 
114  // Set up area films
115  if (areaFilms_.empty())
116  {
117  for
118  (
119  const areaFilm& regionFa
120  : mesh.time().objectRegistry::template csorted<areaFilm>()
121  )
122  {
123  areaFilms_.push_back(const_cast<areaFilm*>(&regionFa));
124  }
125  }
126 }
127 
128 
129 template<class CloudType>
130 void Foam::KinematicSurfaceFilm<CloudType>::init(bool binitThermo)
131 {
132  if (binitThermo)
133  {
134  this->coeffDict().readEntry("pRef", pRef_);
135  this->coeffDict().readEntry("TRef", TRef_);
136  thermo_ = new liquidMixtureProperties(this->coeffDict().subDict("thermo"));
137  }
138 }
139 
140 
141 template<class CloudType>
142 template<class filmType>
144 (
145  filmType& film,
146  const parcelType& p,
147  const polyPatch& pp,
148  const label facei,
149  const scalar mass,
150  bool& keepParticle
151 )
152 {
153  DebugInfo<< "Parcel " << p.origId() << " absorbInteraction" << endl;
154 
155  // Patch face normal
156  const vector& nf = pp.faceNormals()[facei];
157 
158  // Patch velocity
159  const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
160 
161  // Relative parcel velocity
162  const vector Urel(p.U() - Up);
163 
164  // Parcel normal velocity
165  const vector Un(nf*(Urel & nf));
166 
167  // Parcel tangential velocity
168  const vector Ut(Urel - Un);
169 
170  film.addSources
171  (
172  pp.index(),
173  facei,
174  mass, // mass
175  mass*Ut, // tangential momentum
176  mass*mag(Un), // impingement pressure
177  0 // energy
178  );
179 
180  this->nParcelsTransferred()++;
181 
182  this->totalMassTransferred() += mass;
184  keepParticle = false;
185 }
186 
187 
188 template<class CloudType>
190 (
191  parcelType& p,
192  const polyPatch& pp,
193  const label facei,
194  bool& keepParticle
195 ) const
196 {
197  DebugInfo<< "Parcel " << p.origId() << " bounceInteraction" << endl;
198 
199  // Patch face normal
200  const vector& nf = pp.faceNormals()[facei];
201 
202  // Patch velocity
203  const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
204 
205  // Relative parcel velocity
206  const vector Urel(p.U() - Up);
207 
208  // Flip parcel normal velocity component
209  p.U() -= 2.0*nf*(Urel & nf);
210 
211  keepParticle = true;
212 }
213 
214 
215 template<class CloudType>
216 template<class filmType>
218 (
219  filmType& filmModel,
220  const scalar sigma,
221  const scalar mu,
222  const parcelType& p,
223  const polyPatch& pp,
224  const label facei,
225  bool& keepParticle
226 )
227 {
228  DebugInfo<< "Parcel " << p.origId() << " drySplashInteraction" << endl;
229 
230  // Patch face velocity and normal
231  const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
232  const vector& nf = pp.faceNormals()[facei];
233 
234  // Local pressure
235  //const scalar pc = thermo_.thermo().p()[p.cell()];
236 
237  // Retrieve parcel properties
238  const scalar m = p.mass()*p.nParticle();
239  const scalar rho = p.rho();
240  const scalar d = p.d();
241  const vector Urel(p.U() - Up);
242  const vector Un(nf*(Urel & nf));
243 
244  // Laplace number
245  const scalar La = rho*sigma*d/sqr(mu);
246 
247  // Weber number
248  const scalar We = rho*magSqr(Un)*d/sigma;
249 
250  // Critical Weber number
251  const scalar Wec = Adry_*pow(La, -0.183);
252 
253  if (We < Wec) // Adhesion - assume absorb
254  {
255  absorbInteraction<filmType>
256  (filmModel, p, pp, facei, m, keepParticle);
257  }
258  else // Splash
259  {
260  // Ratio of incident mass to splashing mass
261  const scalar mRatio = 0.2 + 0.6*rndGen_.sample01<scalar>();
262  splashInteraction<filmType>
263  (filmModel, p, pp, facei, mRatio, We, Wec, sigma, keepParticle);
264  }
265 }
266 
267 
268 template<class CloudType>
269 template<class filmType>
271 (
272  filmType& filmModel,
273  const scalar sigma,
274  const scalar mu,
275  parcelType& p,
276  const polyPatch& pp,
277  const label facei,
278  bool& keepParticle
279 )
280 {
281  DebugInfo<< "Parcel " << p.origId() << " wetSplashInteraction" << endl;
282 
283  // Patch face velocity and normal
284  const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
285  const vector& nf = pp.faceNormals()[facei];
286 
287  // Retrieve parcel properties
288  const scalar m = p.mass()*p.nParticle();
289  const scalar rho = p.rho();
290  const scalar d = p.d();
291  vector& U = p.U();
292  const vector Urel(p.U() - Up);
293  const vector Un(nf*(Urel & nf));
294  const vector Ut(Urel - Un);
295 
296  // Laplace number
297  const scalar La = rho*sigma*d/sqr(mu);
298 
299  // Weber number
300  const scalar We = rho*magSqr(Un)*d/sigma;
301 
302  // Critical Weber number
303  const scalar Wec = Awet_*pow(La, -0.183);
304 
305  if (We < 2) // Adhesion - assume absorb
306  {
307  absorbInteraction<filmType>
308  (filmModel, p, pp, facei, m, keepParticle);
309  }
310  else if ((We >= 2) && (We < 20)) // Bounce
311  {
312  // Incident angle of impingement
313  const scalar theta = piByTwo - acos(U/mag(U) & nf);
314 
315  // Restitution coefficient
316  const scalar epsilon = 0.993 - theta*(1.76 - theta*(1.56 - theta*0.49));
317 
318  // Update parcel velocity
319  U = -epsilon*(Un) + 5.0/7.0*(Ut);
320 
321  keepParticle = true;
322  return;
323  }
324  else if ((We >= 20) && (We < Wec)) // Spread - assume absorb
325  {
326  absorbInteraction<filmType>
327  (filmModel, p, pp, facei, m, keepParticle);
328  }
329  else // Splash
330  {
331  // Ratio of incident mass to splashing mass
332  // splash mass can be > incident mass due to film entrainment
333  const scalar mRatio = 0.2 + 0.9*rndGen_.sample01<scalar>();
334  splashInteraction<filmType>
335  (filmModel, p, pp, facei, mRatio, We, Wec, sigma, keepParticle);
336  }
337 }
338 
339 
340 template<class CloudType>
341 template<class filmType>
343 (
344  filmType& filmModel,
345  const parcelType& p,
346  const polyPatch& pp,
347  const label facei,
348  const scalar mRatio,
349  const scalar We,
350  const scalar Wec,
351  const scalar sigma,
352  bool& keepParticle
353 )
354 {
355  // Patch face velocity and normal
356  const fvMesh& mesh = this->owner().mesh();
357  const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
358  const vector& nf = pp.faceNormals()[facei];
359 
360  // Determine direction vectors tangential to patch normal
361  const vector tanVec1(tangentVector(nf));
362  const vector tanVec2(nf^tanVec1);
363 
364  // Retrieve parcel properties
365  const scalar np = p.nParticle();
366  const scalar m = p.mass()*np;
367  const scalar d = p.d();
368  const vector Urel(p.U() - Up);
369  const vector Un(nf*(Urel & nf));
370  const vector Ut(Urel - Un);
371  const vector& posC = mesh.C()[p.cell()];
372  const vector& posCf = mesh.Cf().boundaryField()[pp.index()][facei];
373 
374  // Total mass of (all) splashed parcels
375  const scalar mSplash = m*mRatio;
376 
377  // Number of splashed particles per incoming particle
378  const scalar Ns = 5.0*(We/Wec - 1.0);
379 
380  // Average diameter of splashed particles
381  const scalar dBarSplash = 1/cbrt(6.0)*cbrt(mRatio/Ns)*d + ROOTVSMALL;
382 
383  // Cumulative diameter splash distribution
384  const scalar dMax = dMaxSplash_ > 0 ? dMaxSplash_ : 0.9*cbrt(mRatio)*d;
385  const scalar dMin = dMinSplash_ > 0 ? dMinSplash_ : 0.1*dMax;
386  const scalar K = exp(-dMin/dBarSplash) - exp(-dMax/dBarSplash);
387 
388  // Surface energy of secondary parcels [J]
389  scalar ESigmaSec = 0;
390 
391  // Sample splash distribution to determine secondary parcel diameters
392  scalarList dNew(parcelsPerSplash_);
393  scalarList npNew(parcelsPerSplash_);
394  forAll(dNew, i)
395  {
396  const scalar y = rndGen_.sample01<scalar>();
397  dNew[i] = -dBarSplash*log(exp(-dMin/dBarSplash) - y*K);
398  npNew[i] = mRatio*np*pow3(d)/pow3(dNew[i])/parcelsPerSplash_;
399  ESigmaSec += npNew[i]*sigma*p.areaS(dNew[i]);
400  }
401 
402  // Incident kinetic energy [J]
403  const scalar EKIn = 0.5*m*magSqr(Un);
404 
405  // Incident surface energy [J]
406  const scalar ESigmaIn = np*sigma*p.areaS(d);
407 
408  // Dissipative energy
409  const scalar Ed = max(0.8*EKIn, np*Wec/12*pi*sigma*sqr(d));
410 
411  // Total energy [J]
412  const scalar EKs = EKIn + ESigmaIn - ESigmaSec - Ed;
413 
414  // Switch to absorb if insufficient energy for splash
415  if (EKs <= 0)
416  {
417  absorbInteraction<filmType>
418  (filmModel, p, pp, facei, m, keepParticle);
419  return;
420  }
421 
422  // Helper variables to calculate magUns0
423  const scalar logD = log(d);
424  const scalar coeff2 = log(dNew[0]) - logD + ROOTVSMALL;
425  scalar coeff1 = 0.0;
426 
427  // Note: loop from i = 1 to (p-1)
428  for (int i = 1; i < parcelsPerSplash_; ++i)
429  {
430  coeff1 += sqr(log(dNew[i]) - logD);
431  }
432 
433  // Magnitude of the normal velocity of the first splashed parcel
434  const scalar magUns0 =
435  sqrt(2.0*parcelsPerSplash_*EKs/mSplash/(1.0 + coeff1/sqr(coeff2)));
436 
437  // Set splashed parcel properties
438  forAll(dNew, i)
439  {
440  const vector dirVec = splashDirection(tanVec1, tanVec2, -nf);
441 
442  // Create a new parcel by copying source parcel
443  parcelType* pPtr = new parcelType(p);
444 
445  pPtr->origId() = pPtr->getNewParticleID();
446 
447  pPtr->origProc() = Pstream::myProcNo();
448 
449  if (splashParcelType_ >= 0)
450  {
451  pPtr->typeId() = splashParcelType_;
452  }
453 
454  // Perturb new parcels towards the owner cell centre
455  pPtr->track(0.5*rndGen_.sample01<scalar>()*(posC - posCf), 0);
456 
457  pPtr->nParticle() = npNew[i];
458 
459  pPtr->d() = dNew[i];
460 
461  pPtr->U() = dirVec*(mag(Cf_*Ut) + magUns0*(log(dNew[i]) - logD)/coeff2);
462 
463  // Apply correction to velocity for 2-D cases
465 
466  // Add the new parcel
467  this->owner().addParticle(pPtr);
468 
469  nParcelsSplashed_++;
470  }
471 
472  // Transfer remaining part of parcel to film 0 - splashMass can be -ve
473  // if entraining from the film
474  const scalar mDash = m - mSplash;
475  absorbInteraction<filmType>
476  (filmModel, p, pp, facei, mDash, keepParticle);
477 }
478 
479 
480 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
481 
482 template<class CloudType>
484 (
485  const dictionary& dict,
486  CloudType& owner,
487  const word& type,
488  bool initThermo
489 )
490 :
492  rndGen_(owner.rndGen()),
493  thermo_(nullptr),
494  filmModel_(nullptr),
495  areaFilms_(),
496  interactionType_
497  (
498  interactionTypeNames.get("interactionType", this->coeffDict())
499  ),
500  parcelTypes_(this->coeffDict().getOrDefault("parcelTypes", labelList())),
501  deltaWet_(Zero),
502  splashParcelType_(0),
503  parcelsPerSplash_(0),
504  dMaxSplash_(-1),
505  dMinSplash_(-1),
506  Adry_(Zero),
507  Awet_(Zero),
508  Cf_(Zero),
509  nParcelsSplashed_(0)
510 {
511  Info<< " Applying " << interactionTypeNames[interactionType_]
512  << " interaction model" << endl;
513 
515  {
516  this->coeffDict().readEntry("deltaWet", deltaWet_);
518  this->coeffDict().getOrDefault("splashParcelType", -1);
519  parcelsPerSplash_ =
520  this->coeffDict().getOrDefault("parcelsPerSplash", 2);
521  dMinSplash_ = this->coeffDict().getOrDefault("dMinSplash", -1);
522  dMaxSplash_ = this->coeffDict().getOrDefault("dMaxSplash", -1);
523 
524  this->coeffDict().readEntry("Adry", Adry_);
525  this->coeffDict().readEntry("Awet", Awet_);
526  this->coeffDict().readEntry("Cf", Cf_);
527  init(initThermo);
528  }
529 }
530 
531 
532 template<class CloudType>
534 (
536  bool initThermo
537 )
538 :
540  rndGen_(sfm.rndGen_),
541  thermo_(nullptr),
542  filmModel_(nullptr),
543  areaFilms_(),
544  interactionType_(sfm.interactionType_),
545  parcelTypes_(sfm.parcelTypes_),
546  deltaWet_(sfm.deltaWet_),
547  splashParcelType_(sfm.splashParcelType_),
548  parcelsPerSplash_(sfm.parcelsPerSplash_),
549  dMaxSplash_(sfm.dMaxSplash_),
550  dMinSplash_(sfm.dMinSplash_),
551  Adry_(sfm.Adry_),
552  Awet_(sfm.Awet_),
553  Cf_(sfm.Cf_),
554  nParcelsSplashed_(sfm.nParcelsSplashed_)
555 {
557  {
558  init(initThermo);
559  }
560 }
561 
562 
563 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
564 
565 template<class CloudType>
567 (
568  parcelType& p,
569  const polyPatch& pp,
570  bool& keepParticle
571 )
572 {
573  if (parcelTypes_.size() && parcelTypes_.find(p.typeId()) == -1)
574  {
575  if (debug)
576  {
577  Pout<< "transferParcel: ignoring particle with typeId="
578  << p.typeId()
579  << endl;
580  }
581 
582  return false;
583  }
584 
585  const label patchi = pp.index();
586  const label meshFacei = p.face();
587  const label facei = pp.whichFace(meshFacei);
588 
589  initFilmModels();
590 
591  // Check the singleLayer film models
592  if (this->filmModel_ && this->filmModel_->isRegionPatch(patchi))
593  {
594  auto& film = *filmModel_;
595 
596  switch (interactionType_)
597  {
598  case interactionType::bounce:
599  {
600  bounceInteraction(p, pp, facei, keepParticle);
601 
602  break;
603  }
604 
605  case interactionType::absorb:
606  {
607  const scalar m = p.nParticle()*p.mass();
608 
609  absorbInteraction<regionFilm>
610  (film, p, pp, facei, m, keepParticle);
611 
612  break;
613  }
614 
615  case interactionType::splashBai:
616  {
617  const scalarField X(thermo_->size(), 1);
618  const scalar mu = thermo_->mu(pRef_, TRef_, X);
619  const scalar sigma = thermo_->sigma(pRef_, TRef_, X);
620 
621  const bool dry
622  (
623  this->deltaFilmPatch_[patchi][facei] < deltaWet_
624  );
625 
626  if (dry)
627  {
628  drySplashInteraction<regionFilm>
629  (film, sigma, mu, p, pp, facei, keepParticle);
630  }
631  else
632  {
633  wetSplashInteraction<regionFilm>
634  (film, sigma, mu, p, pp, facei, keepParticle);
635  }
636 
637  break;
638  }
639 
640  default:
641  {
643  << "Unknown interaction type enumeration"
644  << abort(FatalError);
645  }
646  }
647 
648  // Transfer parcel/parcel interactions complete
649  return true;
650  }
651 
652 
653  // Check the area film models
654  for (areaFilm& film : this->areaFilms_)
655  {
656  const label filmFacei
657  (
658  film.isRegionPatch(patchi)
659  ? film.regionMesh().whichFace(meshFacei)
660  : -1
661  );
662 
663  if (filmFacei < 0)
664  {
665  // Film model does not include this patch face
666  continue;
667  }
668 
669  switch (interactionType_)
670  {
671  case interactionType::bounce:
672  {
673  bounceInteraction(p, pp, facei, keepParticle);
674 
675  break;
676  }
677 
678  case interactionType::absorb:
679  {
680  const scalar m = p.nParticle()*p.mass();
681 
682  absorbInteraction<areaFilm>
683  (
684  film, p, pp, facei, m, keepParticle
685  );
686  break;
687  }
688 
689  case interactionType::splashBai:
690  {
691  auto& liqFilm =
692  refCast
694  (
695  film
696  );
697 
698  const scalarField X(liqFilm.thermo().size(), 1);
699  const scalar pRef = film.pRef();
700  const scalar TRef = liqFilm.Tref();
701 
702  const scalar mu = liqFilm.thermo().mu(pRef, TRef, X);
703  const scalar sigma = liqFilm.thermo().sigma(pRef, TRef, X);
704 
705  const bool dry = film.h()[filmFacei] < this->deltaWet_;
706 
707  if (dry)
708  {
709  drySplashInteraction<areaFilm>
710  (film, sigma, mu, p, pp, facei, keepParticle);
711  }
712  else
713  {
714  wetSplashInteraction<areaFilm>
715  (film, sigma, mu, p, pp, facei, keepParticle);
716  }
717 
718  break;
719  }
720 
721  default:
722  {
724  << "Unknown interaction type enumeration"
725  << abort(FatalError);
726  }
727  }
728 
729  // Transfer parcel/parcel interactions complete
730  return true;
731  }
732 
733  // Parcel not interacting with film
734  return false;
735 }
736 
737 
738 template<class CloudType>
740 (
741  const label filmPatchi,
742  const label primaryPatchi,
744 )
745 {
747  (
748  filmPatchi,
749  primaryPatchi,
750  filmModel
751  );
752 }
753 
754 
755 template<class CloudType>
757 (
758  const areaFilm& film
759 )
760 {
762 }
763 
764 
765 template<class CloudType>
767 (
768  parcelType& p,
769  const label filmFacei
770 ) const
771 {
773 }
774 
775 
776 template<class CloudType>
778 {
780 
781  label nSplash0 = this->template getModelProperty<label>("nParcelsSplashed");
782  label nSplashTotal =
783  nSplash0 + returnReduce(nParcelsSplashed_, sumOp<label>());
784 
785  Log_<< " - new splash parcels = " << nSplashTotal << endl;
786 
787  if (this->writeTime())
788  {
789  this->setModelProperty("nParcelsSplashed", nSplashTotal);
790  nParcelsSplashed_ = 0;
791  }
792 }
793 
794 
795 // ************************************************************************* //
interactionType interactionType_
Interaction type enumeration.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
interactionType
Options for the interaction types.
KinematicSurfaceFilm(const dictionary &dict, CloudType &owner, const word &type=typeName, bool initThermo=true)
Construct from components.
dictionary dict
dimensionedScalar acos(const dimensionedScalar &ds)
bool isRegionPatch(const label patchi) const
True if patchi on the primary region is coupled to this region.
void wetSplashInteraction(filmType &, const scalar sigma, const scalar mu, parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle)
Parcel interaction with wetted surface.
dimensionedScalar log(const dimensionedScalar &ds)
scalar Adry_
Dry surface roughness coefficient.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const areaScalarField & h() const
Access const reference h.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
scalar Cf_
Skin friction typically in the range 0.6 < Cf < 0.8.
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:882
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
Random rndGen
Definition: createFields.H:23
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
label whichFace(const label meshFacei) const
The area-face corresponding to the mesh-face, -1 if not found.
Definition: faMeshI.H:129
void bounceInteraction(parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle) const
Bounce parcel (flip parcel normal velocity)
static const Enum< interactionType > interactionTypeNames
Names for interactionType.
void splashInteraction(filmType &, const parcelType &p, const polyPatch &pp, const label facei, const scalar mRatio, const scalar We, const scalar Wec, const scalar sigma, bool &keepParticle)
Bai parcel splash interaction model.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
CGAL::Exact_predicates_exact_constructions_kernel K
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
scalar deltaWet_
Film thickness beyond which patch is assumed to be wet.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void init(bool binitThermo)
Initialise thermo.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
#define Log_
Report write to Foam::Info if the class log switch is true.
dimensionedScalar exp(const dimensionedScalar &ds)
constexpr scalar twoPi(2 *M_PI)
Mathematical constants.
A class for handling words, derived from Foam::string.
Definition: word.H:63
vector splashDirection(const vector &tanVec1, const vector &tanVec2, const vector &nf) const
Return splashed parcel direction.
virtual void cacheFilmFields(const areaFilm &film)
Cache the film fields in preparation for injection.
dimensionedScalar cbrt(const dimensionedScalar &ds)
const faMesh & regionMesh() const
Return the region mesh database.
void absorbInteraction(filmType &, const parcelType &p, const polyPatch &pp, const label facei, const scalar mass, bool &keepParticle)
Absorb parcel into film.
constexpr scalar pi(M_PI)
constexpr scalar piByTwo(0.5 *M_PI)
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:122
Urel
Definition: pEqn.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
vector tangentVector(const vector &v) const
Return a vector tangential to input vector, v.
dimensionedScalar sin(const dimensionedScalar &ds)
Kinematic parcel surface film model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
virtual void info()
Write surface film info.
virtual bool transferParcel(parcelType &p, const polyPatch &pp, bool &keepParticle)
Transfer parcel from cloud to surface film.
void initFilmModels()
Initialise pointers of films.
scalar dMinSplash_
Minimum splash particle diameter for Chi-square distribution.
Templated wall surface film model class.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
scalar Awet_
Wet surface roughness coefficient.
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionFilm &)
Cache the film fields in preparation for injection.
U
Definition: pEqn.H:72
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
scalar dMaxSplash_
Maximum splash particle diameter for Chi-square distribution.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
scalar epsilon
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar TRef("TRef", dimTemperature, laminarTransport)
void drySplashInteraction(filmType &, const scalar sigma, const scalar mu, const parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle)
Parcel interaction with dry surface.
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:680
const volVectorField & C() const
Return cell centres as volVectorField.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
label splashParcelType_
Splash parcel type label - id assigned to identify parcel for.
virtual void info()
Write surface film info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127