39template<
class CloudType>
46 { interactionType::absorb,
"absorb" },
47 { interactionType::bounce,
"bounce" },
54template<
class CloudType>
61 scalar magTangent = 0.0;
63 while (magTangent < SMALL)
66 tangent = vTest - (vTest & v)*v;
67 magTangent =
mag(tangent);
70 return tangent/magTangent;
74template<
class CloudType>
83 const scalar phiSi =
twoPi*rndGen_.sample01<scalar>();
86 const scalar thetaSi =
degToRad(rndGen_.sample01<scalar>()*(50 - 5) + 5);
90 const scalar dcorr =
cos(thetaSi);
95 return dirVec/
mag(dirVec);
99template<
class CloudType>
108 mesh.time().objectRegistry::template getObjectPtr<regionFilm>
110 "surfaceFilmProperties"
116 if (areaFilms_.empty())
123template<
class CloudType>
128 this->coeffDict().readEntry(
"pRef", pRef_);
129 this->coeffDict().readEntry(
"TRef", TRef_);
135template<
class CloudType>
136template<
class filmType>
147 DebugInfo<<
"Parcel " <<
p.origId() <<
" absorbInteraction" <<
endl;
150 const vector& nf =
pp.faceNormals()[facei];
153 const vector& Up = this->owner().U().boundaryField()[
pp.index()][facei];
174 this->nParcelsTransferred()++;
178 keepParticle =
false;
182template<
class CloudType>
191 DebugInfo<<
"Parcel " <<
p.origId() <<
" bounceInteraction" <<
endl;
194 const vector& nf =
pp.faceNormals()[facei];
197 const vector& Up = this->owner().U().boundaryField()[
pp.index()][facei];
203 p.U() -= 2.0*nf*(
Urel & nf);
209template<
class CloudType>
210template<
class filmType>
222 DebugInfo<<
"Parcel " <<
p.origId() <<
" drySplashInteraction" <<
endl;
225 const vector& Up = this->owner().U().boundaryField()[
pp.index()][facei];
226 const vector& nf =
pp.faceNormals()[facei];
232 const scalar m =
p.mass()*
p.nParticle();
233 const scalar
rho =
p.rho();
234 const scalar d =
p.d();
245 const scalar Wec = Adry_*
pow(La, -0.183);
249 absorbInteraction<filmType>
250 (filmModel,
p,
pp, facei, m, keepParticle);
255 const scalar mRatio = 0.2 + 0.6*rndGen_.sample01<scalar>();
256 splashInteraction<filmType>
257 (filmModel,
p,
pp, facei, mRatio, We, Wec,
sigma, keepParticle);
262template<
class CloudType>
263template<
class filmType>
275 DebugInfo<<
"Parcel " <<
p.origId() <<
" wetSplashInteraction" <<
endl;
278 const vector& Up = this->owner().U().boundaryField()[
pp.index()][facei];
279 const vector& nf =
pp.faceNormals()[facei];
282 const scalar m =
p.mass()*
p.nParticle();
283 const scalar
rho =
p.rho();
284 const scalar d =
p.d();
297 const scalar Wec = Awet_*
pow(La, -0.183);
301 absorbInteraction<filmType>
302 (filmModel,
p,
pp, facei, m, keepParticle);
304 else if ((We >= 2) && (We < 20))
310 const scalar
epsilon = 0.993 - theta*(1.76 - theta*(1.56 - theta*0.49));
318 else if ((We >= 20) && (We < Wec))
320 absorbInteraction<filmType>
321 (filmModel,
p,
pp, facei, m, keepParticle);
327 const scalar mRatio = 0.2 + 0.9*rndGen_.sample01<scalar>();
328 splashInteraction<filmType>
329 (filmModel,
p,
pp, facei, mRatio, We, Wec,
sigma, keepParticle);
334template<
class CloudType>
335template<
class filmType>
351 const vector& Up = this->owner().U().boundaryField()[
pp.index()][facei];
352 const vector& nf =
pp.faceNormals()[facei];
355 const vector tanVec1(tangentVector(nf));
356 const vector tanVec2(nf^tanVec1);
359 const scalar np =
p.nParticle();
360 const scalar m =
p.mass()*np;
361 const scalar d =
p.d();
366 const vector& posCf =
mesh.Cf().boundaryField()[
pp.index()][facei];
369 const scalar mSplash = m*mRatio;
372 const scalar Ns = 5.0*(We/Wec - 1.0);
375 const scalar dBarSplash = 1/
cbrt(6.0)*
cbrt(mRatio/Ns)*d + ROOTVSMALL;
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);
383 scalar ESigmaSec = 0;
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]);
397 const scalar EKIn = 0.5*m*
magSqr(Un);
400 const scalar ESigmaIn = np*
sigma*
p.areaS(d);
406 const scalar EKs = EKIn + ESigmaIn - ESigmaSec - Ed;
411 absorbInteraction<filmType>
412 (filmModel,
p,
pp, facei, m, keepParticle);
417 const scalar logD =
log(d);
418 const scalar coeff2 =
log(dNew[0]) - logD + ROOTVSMALL;
422 for (
int i = 1; i < parcelsPerSplash_; ++i)
424 coeff1 +=
sqr(
log(dNew[i]) - logD);
428 const scalar magUns0 =
429 sqrt(2.0*parcelsPerSplash_*EKs/mSplash/(1.0 + coeff1/
sqr(coeff2)));
434 const vector dirVec = splashDirection(tanVec1, tanVec2, -nf);
437 parcelType* pPtr =
new parcelType(
p);
439 pPtr->origId() = pPtr->getNewParticleID();
443 if (splashParcelType_ >= 0)
445 pPtr->typeId() = splashParcelType_;
449 pPtr->track(0.5*rndGen_.sample01<scalar>()*(posC - posCf), 0);
451 pPtr->nParticle() = npNew[i];
455 pPtr->U() = dirVec*(
mag(Cf_*Ut) + magUns0*(
log(dNew[i]) - logD)/coeff2);
461 this->owner().addParticle(pPtr);
468 const scalar mDash = m - mSplash;
469 absorbInteraction<filmType>
470 (filmModel,
p,
pp, facei, mDash, keepParticle);
476template<
class CloudType>
492 interactionTypeNames.get(
"interactionType", this->coeffDict())
494 parcelTypes_(this->coeffDict().getOrDefault(
"parcelTypes",
labelList())),
496 splashParcelType_(0),
497 parcelsPerSplash_(0),
506 <<
" interaction model" <<
endl;
526template<
class CloudType>
534 rndGen_(sfm.rndGen_),
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_),
548 nParcelsSplashed_(sfm.nParcelsSplashed_)
559template<
class CloudType>
567 if (parcelTypes_.size() && parcelTypes_.find(
p.typeId()) == -1)
571 Pout<<
"transferParcel: ignoring particle with typeId="
579 const label patchi =
pp.index();
580 const label meshFacei =
p.face();
581 const label facei =
pp.whichFace(meshFacei);
586 if (this->filmModel_ && this->filmModel_->isRegionPatch(patchi))
588 auto& film = *filmModel_;
590 switch (interactionType_)
592 case interactionType::bounce:
594 bounceInteraction(
p,
pp, facei, keepParticle);
599 case interactionType::absorb:
601 const scalar m =
p.nParticle()*
p.mass();
603 absorbInteraction<regionFilm>
604 (film,
p,
pp, facei, m, keepParticle);
609 case interactionType::splashBai:
612 const scalar
mu = thermo_->mu(pRef_, TRef_, X);
613 const scalar
sigma = thermo_->sigma(pRef_, TRef_, X);
617 this->deltaFilmPatch_[patchi][facei] < deltaWet_
622 drySplashInteraction<regionFilm>
623 (film,
sigma,
mu,
p,
pp, facei, keepParticle);
627 wetSplashInteraction<regionFilm>
628 (film,
sigma,
mu,
p,
pp, facei, keepParticle);
637 <<
"Unknown interaction type enumeration"
648 for (areaFilm& film : this->areaFilms_)
650 const label filmFacei
652 film.isRegionPatch(patchi)
653 ? film.regionMesh().whichFace(meshFacei)
663 switch (interactionType_)
665 case interactionType::bounce:
667 bounceInteraction(
p,
pp, facei, keepParticle);
672 case interactionType::absorb:
674 const scalar m =
p.nParticle()*
p.mass();
676 absorbInteraction<areaFilm>
678 film,
p,
pp, facei, m, keepParticle
683 case interactionType::splashBai:
693 const scalar pRef = film.pRef();
694 const scalar
TRef = liqFilm.Tref();
696 const scalar
mu = liqFilm.thermo().mu(pRef,
TRef, X);
697 const scalar
sigma = liqFilm.thermo().sigma(pRef,
TRef, X);
699 const bool dry = film.h()[filmFacei] < this->deltaWet_;
703 drySplashInteraction<areaFilm>
704 (film,
sigma,
mu,
p,
pp, facei, keepParticle);
708 wetSplashInteraction<areaFilm>
709 (film,
sigma,
mu,
p,
pp, facei, keepParticle);
718 <<
"Unknown interaction type enumeration"
732template<
class CloudType>
735 const label filmPatchi,
736 const label primaryPatchi,
749template<
class CloudType>
759template<
class CloudType>
763 const label filmFacei
770template<
class CloudType>
775 label nSplash0 = this->
template getModelProperty<label>(
"nParcelsSplashed");
779 Log_<<
" - new splash parcels = " << nSplashTotal <<
endl;
781 if (this->writeTime())
783 this->setModelProperty(
"nParcelsSplashed", nSplashTotal);
784 nParcelsSplashed_ = 0;
CGAL::Exact_predicates_exact_constructions_kernel K
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...
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.
@ splashBai
Bai splash model.
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...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
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.
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
label getNewParticleID() const
Get unique particle creation id.
label origId() const noexcept
Return the particle ID on the originating processor.
label origProc() const noexcept
Return the originating processor ID.
A patch is a list of labels that address the faces in the global face list.
Base class for surface film models.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define Log_
Report write to Foam::Info if the class log switch is true.
#define DebugInfo
Report an information message using Foam::Info.
constexpr scalar pi(M_PI)
constexpr scalar piByTwo(0.5 *M_PI)
constexpr scalar twoPi(2 *M_PI)
Namespace for handling debugging switches.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
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.
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)
static constexpr const zero Zero
Global zero (0).
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.
List< scalar > scalarList
List of scalar.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
dimensionedScalar TRef("TRef", dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Unit conversion functions.