36 template<
class CloudType>
41 "(absorb bounce splashBai)" 48 template<
class CloudType>
52 forAll(interactionTypeNames_, i)
54 if (interactionTypeNames_[i] == it)
61 <<
"Unknown interaction type " << it
62 <<
". Valid interaction types include: " << interactionTypeNames_
69 template<
class CloudType>
75 if (it >= interactionTypeNames_.size())
81 return interactionTypeNames_[it];
87 template<
class CloudType>
96 const scalar phiSi =
twoPi*rndGen_.sample01<scalar>();
99 const scalar thetaSi =
pi/180.0*(rndGen_.sample01<scalar>()*(50 - 5) + 5);
103 const scalar dcorr =
cos(thetaSi);
104 const vector normal = alpha*(tanVec1*
cos(phiSi) + tanVec2*
sin(phiSi));
108 return dirVec/
mag(dirVec);
112 template<
class CloudType>
125 Info<<
"Parcel " << p.origId() <<
" absorbInteraction" <<
endl;
139 const vector& Up = this->owner().U().boundaryField()[pp.
index()][facei];
145 const vector Un = nf*(Urel & nf);
148 const vector Ut = Urel - Un;
153 const scalar pc = carrierThermo.
p()[p.cell()];
162 mass*liq.
Hs(pc, p.T())
165 this->nParcelsTransferred()++;
167 keepParticle =
false;
171 template<
class CloudType>
182 Info<<
"Parcel " << p.origId() <<
" bounceInteraction" <<
endl;
189 const vector& Up = this->owner().U().boundaryField()[pp.
index()][facei];
195 p.U() -= 2.0*nf*(Urel & nf);
201 template<
class CloudType>
213 Info<<
"Parcel " << p.origId() <<
" drySplashInteraction" <<
endl;
226 const vector& Up = this->owner().U().boundaryField()[pp.
index()][facei];
230 const scalar pc = carrierThermo.
p()[p.cell()];
233 const scalar m = p.mass()*p.nParticle();
234 const scalar
rho = p.rho();
235 const scalar d = p.d();
236 const scalar sigma = liq.
sigma(pc, p.T());
237 const scalar
mu = liq.
mu(pc, p.T());
239 const vector Un = nf*(Urel & nf);
242 const scalar La = rho*sigma*d/
sqr(mu);
248 const scalar Wec = Adry_*
pow(La, -0.183);
252 absorbInteraction(filmModel, p, pp, facei, m, keepParticle);
257 const scalar mRatio = 0.2 + 0.6*rndGen_.sample01<scalar>();
259 (filmModel, p, pp, facei, mRatio, We, Wec, sigma, keepParticle);
264 template<
class CloudType>
276 Info<<
"Parcel " << p.origId() <<
" wetSplashInteraction" <<
endl;
289 const vector& Up = this->owner().U().boundaryField()[pp.
index()][facei];
293 const scalar pc = carrierThermo.
p()[p.cell()];
296 const scalar m = p.mass()*p.nParticle();
297 const scalar
rho = p.rho();
298 const scalar d = p.d();
300 const scalar sigma = liq.
sigma(pc, p.T());
301 const scalar
mu = liq.
mu(pc, p.T());
303 const vector Un = nf*(Urel & nf);
304 const vector Ut = Urel - Un;
307 const scalar La = rho*sigma*d/
sqr(mu);
313 const scalar Wec = Awet_*
pow(La, -0.183);
317 absorbInteraction(filmModel, p, pp, facei, m, keepParticle);
319 else if ((We >= 2) && (We < 20))
322 const scalar theta =
pi/2 -
acos(U/
mag(U) & nf);
325 const scalar
epsilon = 0.993 - theta*(1.76 - theta*(1.56 - theta*0.49));
328 U = -epsilon*(Un) + 5.0/7.0*(Ut);
333 else if ((We >= 20) && (We < Wec))
335 absorbInteraction(filmModel, p, pp, facei, m, keepParticle);
341 const scalar mRatio = 0.2 + 0.9*rndGen_.sample01<scalar>();
343 (filmModel, p, pp, facei, mRatio, We, Wec, sigma, keepParticle);
348 template<
class CloudType>
364 const vector& Up = this->owner().U().boundaryField()[pp.
index()][facei];
369 const vector tanVec2 = nf^tanVec1;
372 const scalar np = p.nParticle();
373 const scalar m = p.mass()*np;
374 const scalar d = p.d();
376 const vector Un = nf*(Urel & nf);
377 const vector Ut = Urel - Un;
378 const vector& posC = mesh.
C()[p.cell()];
382 const scalar mSplash = m*mRatio;
385 const scalar Ns = 5.0*(We/Wec - 1.0);
388 const scalar dBarSplash = 1/
cbrt(6.0)*
cbrt(mRatio/Ns)*d + rootVSmall;
391 const scalar dMax = 0.9*
cbrt(mRatio)*d;
392 const scalar dMin = 0.1*dMax;
393 const scalar
K =
exp(-dMin/dBarSplash) -
exp(-dMax/dBarSplash);
396 scalar ESigmaSec = 0;
403 const scalar
y = rndGen_.sample01<scalar>();
404 dNew[i] = -dBarSplash*
log(
exp(-dMin/dBarSplash) - y*K);
405 npNew[i] = mRatio*np*
pow3(d)/
pow3(dNew[i])/parcelsPerSplash_;
406 ESigmaSec += npNew[i]*sigma*p.areaS(dNew[i]);
410 const scalar EKIn = 0.5*m*
magSqr(Un);
413 const scalar ESigmaIn = np*sigma*p.areaS(d);
416 const scalar Ed =
max(0.8*EKIn, np*Wec/12*
pi*sigma*
sqr(d));
419 const scalar EKs = EKIn + ESigmaIn - ESigmaSec - Ed;
424 absorbInteraction(filmModel, p, pp, facei, m, keepParticle);
429 const scalar logD =
log(d);
430 const scalar coeff2 =
log(dNew[0]) - logD + rootVSmall;
434 coeff1 +=
sqr(
log(dNew[i]) - logD);
438 const scalar magUns0 =
439 sqrt(2.0*parcelsPerSplash_*EKs/mSplash/(1.0 + coeff1/
sqr(coeff2)));
444 const vector dirVec = splashDirection(tanVec1, tanVec2, -nf);
449 pPtr->origId() = pPtr->getNewParticleID();
453 if (splashParcelType_ >= 0)
455 pPtr->typeId() = splashParcelType_;
459 pPtr->track(0.5*rndGen_.sample01<scalar>()*(posC - posCf), 0);
461 pPtr->nParticle() = npNew[i];
465 pPtr->U() = dirVec*(
mag(Cf_*Ut) + magUns0*(
log(dNew[i]) - logD)/coeff2);
471 this->owner().addParticle(pPtr);
478 const scalar mDash = m - mSplash;
479 absorbInteraction(filmModel, p, pp, facei, mDash, keepParticle);
485 template<
class CloudType>
493 rndGen_(owner.rndGen()),
498 interactionTypeEnum(this->coeffDict().
lookup(
"interactionType"))
501 splashParcelType_(0),
502 parcelsPerSplash_(0),
508 Info<<
" Applying " << interactionTypeStr(interactionType_)
509 <<
" interaction model" <<
endl;
511 if (interactionType_ == itSplashBai)
513 this->coeffDict().lookup(
"deltaWet") >> deltaWet_;
515 this->coeffDict().lookupOrDefault(
"splashParcelType", -1);
517 this->coeffDict().lookupOrDefault(
"parcelsPerSplash", 2);
518 this->coeffDict().lookup(
"Adry") >> Adry_;
519 this->coeffDict().lookup(
"Awet") >> Awet_;
520 this->coeffDict().lookup(
"Cf") >> Cf_;
525 template<
class CloudType>
548 template<
class CloudType>
555 template<
class CloudType>
567 this->owner().db().time().objectRegistry::template
571 "surfaceFilmProperties" 581 switch (interactionType_)
585 bounceInteraction(p, pp, facei, keepParticle);
591 const scalar m = p.nParticle()*p.mass();
592 absorbInteraction(filmModel, p, pp, facei, m, keepParticle);
598 bool dry = this->deltaFilmPatch_[
patchi][facei] < deltaWet_;
602 drySplashInteraction(filmModel, p, pp, facei, keepParticle);
606 wetSplashInteraction(filmModel, p, pp, facei, keepParticle);
614 <<
"Unknown interaction type enumeration" 628 template<
class CloudType>
631 const label filmPatchi,
632 const label primaryPatchi,
644 refCast<const regionModels::surfaceFilmModels::thermoSingleLayer>
650 filmModel.
toPrimary(filmPatchi, TFilmPatch_);
653 thermalFilmModel.
thermo().
Cpv()().boundaryField()[filmPatchi];
654 filmModel.
toPrimary(filmPatchi, CpFilmPatch_);
658 template<
class CloudType>
662 const label filmFacei
668 p.T() = TFilmPatch_[filmFacei];
669 p.Cp() = CpFilmPatch_[filmFacei];
673 template<
class CloudType>
678 label nSplash0 = this->
template getModelProperty<label>(
"nParcelsSplashed");
682 os <<
" New film splash parcels = " << nSplashTotal <<
endl;
684 if (this->writeTime())
686 this->setModelProperty(
"nParcelsSplashed", nSplashTotal);
687 nParcelsSplashed_ = 0;
scalar Awet_
Wet surface roughness coefficient.
virtual void info(Ostream &os)
Write surface film info to stream.
Thermo parcel surface film model.
virtual scalar Hs(scalar p, scalar T) const =0
Liquid sensible enthalpy [J/kg].
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource)=0
External hook to add sources to the film.
scalar Cf_
Skin friction typically in the range 0.6 < Cf < 0.8.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
dimensionedScalar log(const dimensionedScalar &ds)
Base class for surface film models.
A list of keyword definitions, which are a keyword followed by any number of values (e...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
bool isRegionPatch(const label primaryPatchi) const
Return true if patchi on the primary region is a coupled.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
ThermoSurfaceFilm(const dictionary &dict, CloudType &owner)
Construct from components.
scalar Adry_
Dry surface roughness coefficient.
virtual ~ThermoSurfaceFilm()
Destructor.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
scalarList TFilmPatch_
Film temperature / patch face.
CGAL::Exact_predicates_exact_constructions_kernel K
virtual volScalarField & p()=0
Pressure [Pa].
virtual tmp< volScalarField > Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
void drySplashInteraction(regionModels::surfaceFilmModels::surfaceFilmRegionModel &, const parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle)
Parcel interaction with dry surface.
Macros for easy insertion into run-time selection tables.
Templated base class for thermodynamic cloud.
virtual void info(Ostream &os)
Write surface film info to stream.
virtual scalar sigma(scalar p, scalar T) const =0
Surface tension [N/m].
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
void toPrimary(const label regionPatchi, List< Type > ®ionField) const
Convert a local region field to the primary region.
interactionType interactionType_
Interaction type enumeration.
void absorbInteraction(regionModels::surfaceFilmModels::surfaceFilmRegionModel &, const parcelType &p, const polyPatch &pp, const label facei, const scalar mass, bool &keepParticle)
Absorb parcel into film.
stressControl lookup("compactNormalStress") >> compactNormalStress
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
interactionType interactionTypeEnum(const word &it) const
static wordList interactionTypeNames_
Word descriptions of interaction type names.
A class for handling words, derived from string.
Base-class for fluid thermodynamic properties.
dimensionedScalar cbrt(const dimensionedScalar &ds)
The thermophysical properties of a liquid.
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionModels::surfaceFilmModels::surfaceFilmRegionModel &filmModel)
Cache the film fields in preparation for injection.
const Field< PointType > & faceNormals() const
Return face normals for patch.
errorManip< error > abort(error &err)
const scalar twoPi(2 *pi)
vector splashDirection(const vector &tanVec1, const vector &tanVec2, const vector &nf) const
Return splashed parcel direction.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionedScalar mu
Atomic mass unit.
const rhoThermo & thermo() const
Film thermo.
virtual void setParcelProperties(parcelType &p, const label filmFacei) const
Set the individual parcel properties.
void bounceInteraction(parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle) const
Bounce parcel (flip parcel normal velocity)
Vector< Cmpt > perpendicular(const Vector< Cmpt > &v)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const liquidMixtureProperties & liquids() const
Return reference to the global (additional) liquids.
virtual const volScalarField & T() const =0
Temperature [K].
CloudType::parcelType parcelType
Convenience typedef to the cloud's parcel type.
void splashInteraction(regionModels::surfaceFilmModels::surfaceFilmRegionModel &, 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.
Mesh data needed to do the Finite Volume discretisation.
label index() const
Return the index of this patch in the boundaryMesh.
label parcelsPerSplash_
Number of new parcels resulting from splash event.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void cacheFilmFields(const label filmPatchi, const label primaryPatchi, const regionModels::surfaceFilmModels::surfaceFilmRegionModel &)
Cache the film fields in preparation for injection.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
label nParcelsSplashed_
Counter for number of new splash parcels.
Random & rndGen_
Reference to the cloud random number generator.
const volVectorField & C() const
Return cell centres as volVectorField.
const PtrList< liquidProperties > & properties() const
Return the liquid properties.
A patch is a list of labels that address the faces in the global face list.
word interactionTypeStr(const interactionType &it) const
virtual scalar mu(scalar p, scalar T) const =0
Liquid viscosity [Pa s].
virtual bool transferParcel(parcelType &p, const polyPatch &pp, bool &keepParticle)
Transfer parcel from cloud to surface film.
void wetSplashInteraction(regionModels::surfaceFilmModels::surfaceFilmRegionModel &, parcelType &p, const polyPatch &pp, const label facei, bool &keepParticle)
Parcel interaction with wetted surface.
scalar deltaWet_
Film thickness beyond which patch is assumed to be wet.
Templated base class for dsmc cloud.
label splashParcelType_
Splash parcel type label - id assigned to identify parcel for.
scalarList CpFilmPatch_
Film specific heat capacity / patch face.
Thermodynamic form of single-cell layer surface film model.
label whichFace(const label l) const
Return label of face in patch from global face label.