Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
30#include "liquidFilmModel.H"
32#include "unitConversion.H"
33#include "Pstream.H"
34
35using namespace Foam::constant::mathematical;
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39template<class CloudType>
40const Foam::Enum
41<
43>
45({
46 { interactionType::absorb, "absorb" },
47 { interactionType::bounce, "bounce" },
49});
50
51
52// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53
54template<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
74template<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
99template<class CloudType>
101{
102 const fvMesh& mesh = this->owner().mesh();
103
104 // Set up region film. This is attached to Time
105 if (!filmModel_)
106 {
107 filmModel_ =
108 mesh.time().objectRegistry::template getObjectPtr<regionFilm>
109 (
110 "surfaceFilmProperties"
111 );
112 }
113
114 // Check finite-area films
115 // - allow non-const access to the area films
116 if (areaFilms_.empty())
119 }
120}
121
122
123template<class CloudType>
125{
126 if (binitThermo)
127 {
128 this->coeffDict().readEntry("pRef", pRef_);
129 this->coeffDict().readEntry("TRef", TRef_);
130 thermo_ = new liquidMixtureProperties(this->coeffDict().subDict("thermo"));
131 }
132}
133
134
135template<class CloudType>
136template<class filmType>
138(
139 filmType& film,
140 const parcelType& p,
141 const polyPatch& pp,
142 const label facei,
143 const scalar mass,
144 bool& keepParticle
145)
146{
147 DebugInfo<< "Parcel " << p.origId() << " absorbInteraction" << endl;
148
149 // Patch face normal
150 const vector& nf = pp.faceNormals()[facei];
151
152 // Patch velocity
153 const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
154
155 // Relative parcel velocity
156 const vector Urel(p.U() - Up);
157
158 // Parcel normal velocity
159 const vector Un(nf*(Urel & nf));
160
161 // Parcel tangential velocity
162 const vector Ut(Urel - Un);
163
164 film.addSources
165 (
166 pp.index(),
167 facei,
168 mass, // mass
169 mass*Ut, // tangential momentum
170 mass*mag(Un), // impingement pressure
171 0 // energy
172 );
173
174 this->nParcelsTransferred()++;
175
176 this->totalMassTransferred() += mass;
177
178 keepParticle = false;
179}
180
181
182template<class CloudType>
184(
185 parcelType& p,
186 const polyPatch& pp,
187 const label facei,
188 bool& keepParticle
189) const
190{
191 DebugInfo<< "Parcel " << p.origId() << " bounceInteraction" << endl;
192
193 // Patch face normal
194 const vector& nf = pp.faceNormals()[facei];
195
196 // Patch velocity
197 const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
198
199 // Relative parcel velocity
200 const vector Urel(p.U() - Up);
201
202 // Flip parcel normal velocity component
203 p.U() -= 2.0*nf*(Urel & nf);
205 keepParticle = true;
206}
207
208
209template<class CloudType>
210template<class filmType>
212(
213 filmType& filmModel,
214 const scalar sigma,
215 const scalar mu,
216 const parcelType& p,
217 const polyPatch& pp,
218 const label facei,
219 bool& keepParticle
220)
221{
222 DebugInfo<< "Parcel " << p.origId() << " drySplashInteraction" << endl;
223
224 // Patch face velocity and normal
225 const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
226 const vector& nf = pp.faceNormals()[facei];
227
228 // Local pressure
229 //const scalar pc = thermo_.thermo().p()[p.cell()];
230
231 // Retrieve parcel properties
232 const scalar m = p.mass()*p.nParticle();
233 const scalar rho = p.rho();
234 const scalar d = p.d();
235 const vector Urel(p.U() - Up);
236 const vector Un(nf*(Urel & nf));
237
238 // Laplace number
239 const scalar La = rho*sigma*d/sqr(mu);
240
241 // Weber number
242 const scalar We = rho*magSqr(Un)*d/sigma;
243
244 // Critical Weber number
245 const scalar Wec = Adry_*pow(La, -0.183);
246
247 if (We < Wec) // Adhesion - assume absorb
248 {
249 absorbInteraction<filmType>
250 (filmModel, p, pp, facei, m, keepParticle);
251 }
252 else // Splash
253 {
254 // Ratio of incident mass to splashing mass
255 const scalar mRatio = 0.2 + 0.6*rndGen_.sample01<scalar>();
256 splashInteraction<filmType>
257 (filmModel, p, pp, facei, mRatio, We, Wec, sigma, keepParticle);
258 }
259}
260
261
262template<class CloudType>
263template<class filmType>
265(
266 filmType& filmModel,
267 const scalar sigma,
268 const scalar mu,
269 parcelType& p,
270 const polyPatch& pp,
271 const label facei,
272 bool& keepParticle
273)
274{
275 DebugInfo<< "Parcel " << p.origId() << " wetSplashInteraction" << endl;
276
277 // Patch face velocity and normal
278 const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
279 const vector& nf = pp.faceNormals()[facei];
280
281 // Retrieve parcel properties
282 const scalar m = p.mass()*p.nParticle();
283 const scalar rho = p.rho();
284 const scalar d = p.d();
285 vector& U = p.U();
286 const vector Urel(p.U() - Up);
287 const vector Un(nf*(Urel & nf));
288 const vector Ut(Urel - Un);
289
290 // Laplace number
291 const scalar La = rho*sigma*d/sqr(mu);
292
293 // Weber number
294 const scalar We = rho*magSqr(Un)*d/sigma;
295
296 // Critical Weber number
297 const scalar Wec = Awet_*pow(La, -0.183);
298
299 if (We < 2) // Adhesion - assume absorb
300 {
301 absorbInteraction<filmType>
302 (filmModel, p, pp, facei, m, keepParticle);
303 }
304 else if ((We >= 2) && (We < 20)) // Bounce
305 {
306 // Incident angle of impingement
307 const scalar theta = piByTwo - acos(U/mag(U) & nf);
308
309 // Restitution coefficient
310 const scalar epsilon = 0.993 - theta*(1.76 - theta*(1.56 - theta*0.49));
311
312 // Update parcel velocity
313 U = -epsilon*(Un) + 5.0/7.0*(Ut);
314
315 keepParticle = true;
316 return;
317 }
318 else if ((We >= 20) && (We < Wec)) // Spread - assume absorb
319 {
320 absorbInteraction<filmType>
321 (filmModel, p, pp, facei, m, keepParticle);
322 }
323 else // Splash
324 {
325 // Ratio of incident mass to splashing mass
326 // splash mass can be > incident mass due to film entrainment
327 const scalar mRatio = 0.2 + 0.9*rndGen_.sample01<scalar>();
328 splashInteraction<filmType>
329 (filmModel, p, pp, facei, mRatio, We, Wec, sigma, keepParticle);
330 }
331}
332
333
334template<class CloudType>
335template<class filmType>
337(
338 filmType& filmModel,
339 const parcelType& p,
340 const polyPatch& pp,
341 const label facei,
342 const scalar mRatio,
343 const scalar We,
344 const scalar Wec,
345 const scalar sigma,
346 bool& keepParticle
347)
348{
349 // Patch face velocity and normal
350 const fvMesh& mesh = this->owner().mesh();
351 const vector& Up = this->owner().U().boundaryField()[pp.index()][facei];
352 const vector& nf = pp.faceNormals()[facei];
353
354 // Determine direction vectors tangential to patch normal
355 const vector tanVec1(tangentVector(nf));
356 const vector tanVec2(nf^tanVec1);
357
358 // Retrieve parcel properties
359 const scalar np = p.nParticle();
360 const scalar m = p.mass()*np;
361 const scalar d = p.d();
362 const vector Urel(p.U() - Up);
363 const vector Un(nf*(Urel & nf));
364 const vector Ut(Urel - Un);
365 const vector& posC = mesh.C()[p.cell()];
366 const vector& posCf = mesh.Cf().boundaryField()[pp.index()][facei];
367
368 // Total mass of (all) splashed parcels
369 const scalar mSplash = m*mRatio;
370
371 // Number of splashed particles per incoming particle
372 const scalar Ns = 5.0*(We/Wec - 1.0);
373
374 // Average diameter of splashed particles
375 const scalar dBarSplash = 1/cbrt(6.0)*cbrt(mRatio/Ns)*d + ROOTVSMALL;
376
377 // Cumulative diameter splash distribution
378 const scalar dMax = dMaxSplash_ > 0 ? dMaxSplash_ : 0.9*cbrt(mRatio)*d;
379 const scalar dMin = dMinSplash_ > 0 ? dMinSplash_ : 0.1*dMax;
380 const scalar K = exp(-dMin/dBarSplash) - exp(-dMax/dBarSplash);
381
382 // Surface energy of secondary parcels [J]
383 scalar ESigmaSec = 0;
384
385 // Sample splash distribution to determine secondary parcel diameters
386 scalarList dNew(parcelsPerSplash_);
387 scalarList npNew(parcelsPerSplash_);
388 forAll(dNew, i)
389 {
390 const scalar y = rndGen_.sample01<scalar>();
391 dNew[i] = -dBarSplash*log(exp(-dMin/dBarSplash) - y*K);
392 npNew[i] = mRatio*np*pow3(d)/pow3(dNew[i])/parcelsPerSplash_;
393 ESigmaSec += npNew[i]*sigma*p.areaS(dNew[i]);
394 }
395
396 // Incident kinetic energy [J]
397 const scalar EKIn = 0.5*m*magSqr(Un);
398
399 // Incident surface energy [J]
400 const scalar ESigmaIn = np*sigma*p.areaS(d);
401
402 // Dissipative energy
403 const scalar Ed = max(0.8*EKIn, np*Wec/12*pi*sigma*sqr(d));
404
405 // Total energy [J]
406 const scalar EKs = EKIn + ESigmaIn - ESigmaSec - Ed;
407
408 // Switch to absorb if insufficient energy for splash
409 if (EKs <= 0)
410 {
411 absorbInteraction<filmType>
412 (filmModel, p, pp, facei, m, keepParticle);
413 return;
414 }
415
416 // Helper variables to calculate magUns0
417 const scalar logD = log(d);
418 const scalar coeff2 = log(dNew[0]) - logD + ROOTVSMALL;
419 scalar coeff1 = 0.0;
420
421 // Note: loop from i = 1 to (p-1)
422 for (int i = 1; i < parcelsPerSplash_; ++i)
423 {
424 coeff1 += sqr(log(dNew[i]) - logD);
425 }
426
427 // Magnitude of the normal velocity of the first splashed parcel
428 const scalar magUns0 =
429 sqrt(2.0*parcelsPerSplash_*EKs/mSplash/(1.0 + coeff1/sqr(coeff2)));
430
431 // Set splashed parcel properties
432 forAll(dNew, i)
433 {
434 const vector dirVec = splashDirection(tanVec1, tanVec2, -nf);
435
436 // Create a new parcel by copying source parcel
437 parcelType* pPtr = new parcelType(p);
438
439 pPtr->origId() = pPtr->getNewParticleID();
440
441 pPtr->origProc() = Pstream::myProcNo();
442
443 if (splashParcelType_ >= 0)
444 {
445 pPtr->typeId() = splashParcelType_;
446 }
447
448 // Perturb new parcels towards the owner cell centre
449 pPtr->track(0.5*rndGen_.sample01<scalar>()*(posC - posCf), 0);
450
451 pPtr->nParticle() = npNew[i];
452
453 pPtr->d() = dNew[i];
454
455 pPtr->U() = dirVec*(mag(Cf_*Ut) + magUns0*(log(dNew[i]) - logD)/coeff2);
456
457 // Apply correction to velocity for 2-D cases
458 meshTools::constrainDirection(mesh, mesh.solutionD(), pPtr->U());
459
460 // Add the new parcel
461 this->owner().addParticle(pPtr);
462
463 nParcelsSplashed_++;
464 }
465
466 // Transfer remaining part of parcel to film 0 - splashMass can be -ve
467 // if entraining from the film
468 const scalar mDash = m - mSplash;
469 absorbInteraction<filmType>
470 (filmModel, p, pp, facei, mDash, keepParticle);
471}
472
473
474// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
475
476template<class CloudType>
478(
479 const dictionary& dict,
480 CloudType& owner,
481 const word& type,
482 bool initThermo
483)
484:
486 rndGen_(owner.rndGen()),
487 thermo_(nullptr),
488 filmModel_(nullptr),
489 areaFilms_(),
490 interactionType_
491 (
492 interactionTypeNames.get("interactionType", this->coeffDict())
493 ),
494 parcelTypes_(this->coeffDict().getOrDefault("parcelTypes", labelList())),
495 deltaWet_(Zero),
496 splashParcelType_(0),
497 parcelsPerSplash_(0),
498 dMaxSplash_(-1),
499 dMinSplash_(-1),
500 Adry_(Zero),
501 Awet_(Zero),
502 Cf_(Zero),
503 nParcelsSplashed_(0)
504{
506 << " interaction model" << endl;
507
509 {
510 this->coeffDict().readEntry("deltaWet", deltaWet_);
512 this->coeffDict().getOrDefault("splashParcelType", -1);
513 parcelsPerSplash_ =
514 this->coeffDict().getOrDefault("parcelsPerSplash", 2);
515 dMinSplash_ = this->coeffDict().getOrDefault("dMinSplash", -1);
516 dMaxSplash_ = this->coeffDict().getOrDefault("dMaxSplash", -1);
517
518 this->coeffDict().readEntry("Adry", Adry_);
519 this->coeffDict().readEntry("Awet", Awet_);
520 this->coeffDict().readEntry("Cf", Cf_);
521 init(initThermo);
522 }
523}
524
525
526template<class CloudType>
528(
530 bool initThermo
531)
532:
534 rndGen_(sfm.rndGen_),
535 thermo_(nullptr),
536 filmModel_(nullptr),
537 areaFilms_(),
538 interactionType_(sfm.interactionType_),
539 parcelTypes_(sfm.parcelTypes_),
540 deltaWet_(sfm.deltaWet_),
541 splashParcelType_(sfm.splashParcelType_),
542 parcelsPerSplash_(sfm.parcelsPerSplash_),
543 dMaxSplash_(sfm.dMaxSplash_),
544 dMinSplash_(sfm.dMinSplash_),
545 Adry_(sfm.Adry_),
546 Awet_(sfm.Awet_),
547 Cf_(sfm.Cf_),
548 nParcelsSplashed_(sfm.nParcelsSplashed_)
549{
551 {
552 init(initThermo);
554}
555
556
557// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
558
559template<class CloudType>
561(
562 parcelType& p,
563 const polyPatch& pp,
564 bool& keepParticle
565)
566{
567 if (parcelTypes_.size() && parcelTypes_.find(p.typeId()) == -1)
568 {
569 if (debug)
570 {
571 Pout<< "transferParcel: ignoring particle with typeId="
572 << p.typeId()
573 << endl;
574 }
575
576 return false;
577 }
578
579 const label patchi = pp.index();
580 const label meshFacei = p.face();
581 const label facei = pp.whichFace(meshFacei);
582
583 initFilmModels();
584
585 // Check the singleLayer film models
586 if (this->filmModel_ && this->filmModel_->isRegionPatch(patchi))
587 {
588 auto& film = *filmModel_;
589
590 switch (interactionType_)
591 {
592 case interactionType::bounce:
593 {
594 bounceInteraction(p, pp, facei, keepParticle);
595
596 break;
597 }
598
599 case interactionType::absorb:
600 {
601 const scalar m = p.nParticle()*p.mass();
602
603 absorbInteraction<regionFilm>
604 (film, p, pp, facei, m, keepParticle);
605
606 break;
607 }
608
609 case interactionType::splashBai:
610 {
611 const scalarField X(thermo_->size(), 1);
612 const scalar mu = thermo_->mu(pRef_, TRef_, X);
613 const scalar sigma = thermo_->sigma(pRef_, TRef_, X);
614
615 const bool dry
616 (
617 this->deltaFilmPatch_[patchi][facei] < deltaWet_
618 );
619
620 if (dry)
621 {
622 drySplashInteraction<regionFilm>
623 (film, sigma, mu, p, pp, facei, keepParticle);
624 }
625 else
626 {
627 wetSplashInteraction<regionFilm>
628 (film, sigma, mu, p, pp, facei, keepParticle);
629 }
630
631 break;
632 }
633
634 default:
635 {
637 << "Unknown interaction type enumeration"
638 << abort(FatalError);
639 }
640 }
641
642 // Transfer parcel/parcel interactions complete
643 return true;
644 }
645
646
647 // Check the area film models
648 for (areaFilm& film : this->areaFilms_)
649 {
650 const label filmFacei
651 (
652 film.isRegionPatch(patchi)
653 ? film.regionMesh().whichFace(meshFacei)
654 : -1
655 );
656
657 if (filmFacei < 0)
658 {
659 // Film model does not include this patch face
660 continue;
661 }
662
663 switch (interactionType_)
664 {
665 case interactionType::bounce:
666 {
667 bounceInteraction(p, pp, facei, keepParticle);
668
669 break;
670 }
671
672 case interactionType::absorb:
673 {
674 const scalar m = p.nParticle()*p.mass();
675
676 absorbInteraction<areaFilm>
677 (
678 film, p, pp, facei, m, keepParticle
679 );
680 break;
681 }
682
683 case interactionType::splashBai:
684 {
685 auto& liqFilm =
686 refCast
688 (
689 film
690 );
691
692 const scalarField X(liqFilm.thermo().size(), 1);
693 const scalar pRef = film.pRef();
694 const scalar TRef = liqFilm.Tref();
695
696 const scalar mu = liqFilm.thermo().mu(pRef, TRef, X);
697 const scalar sigma = liqFilm.thermo().sigma(pRef, TRef, X);
698
699 const bool dry = film.h()[filmFacei] < this->deltaWet_;
700
701 if (dry)
702 {
703 drySplashInteraction<areaFilm>
704 (film, sigma, mu, p, pp, facei, keepParticle);
705 }
706 else
707 {
708 wetSplashInteraction<areaFilm>
709 (film, sigma, mu, p, pp, facei, keepParticle);
710 }
711
712 break;
713 }
714
715 default:
716 {
718 << "Unknown interaction type enumeration"
719 << abort(FatalError);
720 }
721 }
722
723 // Transfer parcel/parcel interactions complete
724 return true;
725 }
727 // Parcel not interacting with film
728 return false;
729}
730
731
732template<class CloudType>
734(
735 const label filmPatchi,
736 const label primaryPatchi,
738)
739{
741 (
742 filmPatchi,
743 primaryPatchi,
744 filmModel
745 );
746}
747
748
749template<class CloudType>
751(
752 const areaFilm& film
761(
762 parcelType& p,
763 const label filmFacei
764) const
765{
767}
768
769
770template<class CloudType>
772{
774
775 label nSplash0 = this->template getModelProperty<label>("nParcelsSplashed");
776 label nSplashTotal =
777 nSplash0 + returnReduce(nParcelsSplashed_, sumOp<label>());
778
779 Log_<< " - new splash parcels = " << nSplashTotal << endl;
780
781 if (this->writeTime())
782 {
783 this->setModelProperty("nParcelsSplashed", nSplashTotal);
784 nParcelsSplashed_ = 0;
785 }
786}
787
788
789// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
scalar y
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
label typeId() const
Return type id.
const vector & U() const
Return const access to velocity.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Kinematic parcel surface film model.
virtual void cacheFilmFields(const areaFilm &film)
Cache the film fields in preparation for injection.
void wetSplashInteraction(filmType &, const scalar sigma, const scalar mu, parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle)
Parcel interaction with wetted surface.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
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.
label splashParcelType_
Splash parcel type label - id assigned to identify parcel for.
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.
KinematicSurfaceFilm(const dictionary &dict, CloudType &owner, const word &type=typeName, bool initThermo=true)
Construct from components.
scalar Cf_
Skin friction typically in the range 0.6 < Cf < 0.8.
liquidMixtureProperties * thermo_
Region Film liquid thermo.
Random & rndGen_
Reference to the cloud random number generator.
scalar dMaxSplash_
Maximum splash particle diameter for Chi-square distribution.
scalar Awet_
Wet surface roughness coefficient.
void absorbInteraction(filmType &, const parcelType &p, const polyPatch &pp, const label facei, const scalar mass, bool &keepParticle)
Absorb parcel into film.
regionModels::areaSurfaceFilmModels::liquidFilmBase areaFilm
vector splashDirection(const vector &tanVec1, const vector &tanVec2, const vector &nf) const
Return splashed parcel direction.
void bounceInteraction(parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle) const
Bounce parcel (flip parcel normal velocity).
labelList parcelTypes_
Particle type IDs that can interact with the film.
static const Enum< interactionType > interactionTypeNames
Names for interactionType.
interactionType
Options for the interaction types.
scalar Adry_
Dry surface roughness coefficient.
vector tangentVector(const vector &v) const
Return a vector tangential to input vector, v.
scalar deltaWet_
Film thickness beyond which patch is assumed to be wet.
interactionType interactionType_
Interaction type enumeration.
void init(bool binitThermo)
Initialise thermo.
virtual void info()
Write surface film info.
scalar dMinSplash_
Minimum splash particle diameter for Chi-square distribution.
UPtrList< areaFilm > areaFilms_
UPointers to area films.
regionFilm * filmModel_
Pointer to filmModel.
void initFilmModels()
Initialise pointers of films.
label parcelsPerSplash_
Number of new parcels resulting from splash event.
CloudType::parcelType parcelType
Convenience typedef to the cloud's parcel type.
scalar pRef_
Region Film reference pressure.
scalar TRef_
Region Film reference temperature.
virtual bool transferParcel(parcelType &p, const polyPatch &pp, bool &keepParticle)
Transfer parcel from cloud to surface film.
label nParcelsSplashed_
Counter for number of new splash parcels.
Templated wall surface film model class.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
label nParcelsTransferred() const noexcept
The number of parcels transferred to the film model.
scalar totalMassTransferred() const noexcept
The total mass transferred.
Field< scalarField > deltaFilmPatch_
Film height of all film patches / patch face.
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionFilm &)
Cache the film fields in preparation for injection.
static UPtrList< areaFilm > sorted_areaFilms(const polyMesh &)
Return a sorted list of area-film objects that are registered on the faMeshesRegistry.
virtual void info()
Write surface film info.
SurfaceFilmModel(CloudType &owner)
Construct null from owner.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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,...
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...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition particle.C:620
label getNewParticleID() const
Get unique particle creation id.
Definition particleI.H:96
label origId() const noexcept
Return the particle ID on the originating processor.
Definition particleI.H:194
label origProc() const noexcept
Return the originating processor ID.
Definition particleI.H:182
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
scalar epsilon
Urel
Definition pEqn.H:56
#define Log_
Report write to Foam::Info if the class log switch is true.
#define DebugInfo
Report an information message using Foam::Info.
Mathematical constants.
constexpr scalar pi(M_PI)
constexpr scalar piByTwo(0.5 *M_PI)
constexpr scalar twoPi(2 *M_PI)
Namespace for handling debugging switches.
Definition debug.C:45
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition meshTools.C:680
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
volScalarField & alpha
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
dimensionedScalar TRef("TRef", dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.
Random rndGen