49VoFPatchTransfer::VoFPatchTransfer
58 coeffDict_.getOrDefault<scalar>(
"deltaFactorToVoF", 1.0)
62 coeffDict_.getOrDefault<scalar>(
"deltaFactorToFilm", 0.5)
66 coeffDict_.getOrDefault<scalar>(
"alphaToVoF", 0.5)
70 coeffDict_.getOrDefault<scalar>(
"alphaToFilm", 0.1)
74 coeffDict_.getOrDefault<scalar>(
"transferRateCoeff", 0.1)
77 patchTransferredMasses_()
79 const polyBoundaryMesh&
pbm = film.regionMesh().boundaryMesh();
83 pbm.size() - film.regionMesh().globalData().processorPatches().size()
87 if (coeffDict_.readIfPresent(
"patches",
patchNames))
92 Info<<
" applying to " << patchIDs_.size() <<
" patches:" <<
nl;
94 for (
const label patchi : patchIDs_)
101 Info<<
" applying to all patches" <<
endl;
106 patchTransferredMasses_.resize(patchIDs_.size(), Zero);
108 if (patchIDs_.empty())
111 <<
"No patches selected"
119VoFPatchTransfer::~VoFPatchTransfer()
125void VoFPatchTransfer::correct
127 scalarField& availableMass,
128 scalarField& massToTransfer
135void VoFPatchTransfer::correct
137 scalarField& availableMass,
138 scalarField& massToTransfer,
139 scalarField& energyToTransfer
143 if (!patchIDs_.size())
return;
149 const polyBoundaryMesh&
pbm = film().regionMesh().boundaryMesh();
152 const twoPhaseMixtureThermo&
thermo
154 film().primaryMesh().lookupObject<twoPhaseMixtureThermo>
156 twoPhaseMixtureThermo::dictName
163 const tmp<volScalarField> trhoVoF(
thermo.thermo1().rho());
169 const label patchi = patchIDs_[pidi];
170 label primaryPatchi = -1;
172 forAll(film().intCoupledPatchIDs(), i)
174 const label filmPatchi = film().intCoupledPatchIDs()[i];
176 if (filmPatchi == patchi)
178 primaryPatchi = film().primaryPatchIDs()[i];
182 if (primaryPatchi != -1)
186 film().primaryMesh().
boundary()[primaryPatchi].deltaCoeffs()
188 film().toRegion(patchi, deltaCoeffs);
190 scalarField hp(heVoF.boundaryField()[primaryPatchi]);
191 film().toRegion(patchi, hp);
193 scalarField Tp(TVoF.boundaryField()[primaryPatchi]);
194 film().toRegion(patchi, Tp);
196 scalarField Cpp(CpVoF.boundaryField()[primaryPatchi]);
197 film().toRegion(patchi, Cpp);
199 scalarField rhop(rhoVoF.boundaryField()[primaryPatchi]);
200 film().toRegion(patchi, rhop);
202 scalarField alphap(alphaVoF.boundaryField()[primaryPatchi]);
203 film().toRegion(patchi, alphap);
207 film().primaryMesh().
boundary()[primaryPatchi]
208 .patchInternalField(film().primaryMesh().
V())
210 film().toRegion(patchi, Vp);
212 const polyPatch&
pp =
pbm[patchi];
216 scalar dMassPatch = 0;
220 const label celli = faceCells[facei];
226 delta[celli] > 2*deltaFactorToVoF_/deltaCoeffs[facei]
227 || alphap[facei] > alphaToVoF_
231 transferRateCoeff_*
delta[celli]*
rho[celli]*magSf[celli];
233 massToTransfer[celli] += dMass;
234 energyToTransfer[celli] += dMass*film().hs()[celli];
241 &&
delta[celli] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
242 && alphap[facei] < alphaToFilm_
246 -transferRateCoeff_*alphap[facei]*rhop[facei]*Vp[facei];
248 massToTransfer[celli] += dMass;
249 energyToTransfer[celli] += dMass*hp[facei];
252 availableMass[celli] -= dMass;
256 patchTransferredMasses_[pidi] += dMassPatch;
257 addToTransferredMass(dMassPatch);
261 transferModel::correct();
267 getModelProperty<scalarField>
269 "patchTransferredMasses",
274 scalarField patchTransferredMassTotals(patchTransferredMasses_);
275 Pstream::listCombineGather
277 patchTransferredMassTotals,
280 patchTransferredMasses0 += patchTransferredMassTotals;
282 setModelProperty<scalarField>
284 "patchTransferredMasses",
285 patchTransferredMasses0
288 patchTransferredMasses_ = 0;
293void VoFPatchTransfer::patchTransferredMassTotals
295 scalarField& patchMasses
299 if (!patchIDs_.size())
return;
303 getModelProperty<scalarField>
305 "patchTransferredMasses",
310 scalarField patchTransferredMassTotals(patchTransferredMasses_);
311 Pstream::listCombineGather(patchTransferredMassTotals, plusEqOp<scalar>());
315 const label patchi = patchIDs_[pidi];
316 patchMasses[patchi] +=
317 patchTransferredMasses[pidi] + patchTransferredMassTotals[pidi];
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Transfer mass between the film and the VoF in the continuous phase.
Base class for surface film models.
Base class for film transfer models, handling mass transfer between the film and the continuous phase...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const expr V(m.psi().mesh().V())
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
constexpr char nl
The newline '\n' character (0x0a).
wordList patchNames(nPatches)
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.