27 #include "alphaContactAngleFvPatchScalarField.H" 32 #include "surfaceInterpolate.H" 41 void Foam::multiphaseMixture::calcAlphas()
46 forAllIter(PtrDictionary<phase>, phases_, iter)
48 alphas_ += level*iter();
74 phases_(
lookup(
"phases"), phase::iNew(U, phi)),
119 sigmas_(
lookup(
"sigmas")),
120 dimSigma_(1, 0, -2, 0, 0),
139 tmp<volScalarField> trho = iter()*iter().rho();
142 for (++iter; iter != phases_.end(); ++iter)
144 rho += iter()*iter().rho();
156 tmp<scalarField> trho = iter().boundaryField()[
patchi]*iter().rho().value();
159 for (++iter; iter != phases_.end(); ++iter)
161 rho += iter().boundaryField()[
patchi]*iter().rho().value();
173 tmp<volScalarField> tmu = iter()*iter().rho()*iter().nu();
176 for (++iter; iter != phases_.end(); ++iter)
178 mu += iter()*iter().rho()*iter().nu();
190 tmp<scalarField> tmu =
191 iter().boundaryField()[
patchi]
192 *iter().rho().value()
196 for (++iter; iter != phases_.end(); ++iter)
199 iter().boundaryField()[
patchi]
200 *iter().rho().value()
213 tmp<surfaceScalarField> tmuf =
217 for (++iter; iter != phases_.end(); ++iter)
251 tmp<surfaceScalarField> tstf
255 "surfaceTensionForce",
265 const phase& alpha1 = iter1();
270 for (; iter2 != phases_.end(); ++iter2)
272 const phase& alpha2 = iter2();
275 sigmas_.
find(interfacePair(alpha1, alpha2));
277 if (sigma == sigmas_.
end())
280 <<
"Cannot find interface " << interfacePair(alpha1, alpha2)
281 <<
" in list of sigma values" 302 const Time& runTime = mesh_.
time();
308 scalar cAlpha(alphaControls.lookup<scalar>(
"cAlpha"));
329 !(++alphaSubCycle).
end();
333 rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
350 forAllIter(PtrDictionary<phase>, phases_, iter)
379 return gradAlphaf/(
mag(gradAlphaf) + deltaN_);
390 return nHatfv(alpha1, alpha2) & mesh_.
Sf();
400 tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
408 tnHatfv.ref().boundaryFieldRef()
419 tmp<volScalarField> tnearInt
432 max(tnearInt(),
pos0(iter() - 0.01)*
pos0(0.99 - iter()));
439 void Foam::multiphaseMixture::solveAlphas
444 static label nSolves=-1;
453 UPtrList<const volScalarField> alphas(phases_.size());
454 PtrList<surfaceScalarField> alphaPhis(phases_.size());
457 forAllIter(PtrDictionary<phase>, phases_, iter)
459 const phase& alpha = iter();
461 alphas.set(phasei, &alpha);
468 "phi" + alpha.name() +
"Corr",
480 forAllIter(PtrDictionary<phase>, phases_, iter2)
482 phase& alpha2 = iter2();
484 if (&alpha2 == &alpha)
continue;
532 forAllIter(PtrDictionary<phase>, phases_, iter)
534 phase& alpha = iter();
544 rhoPhi_ += alphaPhi*alpha.rho();
546 Info<< alpha.name() <<
" volume fraction, min, max = " 547 << alpha.weightedAverage(mesh_.
V()).value()
548 <<
' ' <<
min(alpha).value()
549 <<
' ' <<
max(alpha).value()
557 Info<<
"Phase-sum volume fraction, min, max = " 558 << sumAlpha.weightedAverage(mesh_.
V()).value()
559 <<
' ' <<
min(sumAlpha).value()
560 <<
' ' <<
max(sumAlpha).value()
565 forAllIter(PtrDictionary<phase>, phases_, iter)
567 phase& alpha = iter();
568 alpha += alpha*sumCorr;
579 lookup(
"sigmas") >> sigmas_;
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const surfaceVectorField & Sf() const
Return cell face area vectors.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
virtual bool read()
Read object.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctContactAngle(const volScalarField &alpha1, const volScalarField &alpha2, const volVectorField::Boundary &Ubf, const dimensionedScalar &deltaN, surfaceVectorField::Boundary &nHatbf)
Correct the contact angle for the two volume fraction fields.
const word alphaScheme(mesh.schemes().div(divAlphaName)[1].wordToken())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
IOobject(const word &name, const fileName &instance, const objectRegistry ®istry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
tmp< volScalarField > trho
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
surfaceScalarField phir(fvc::flux(mixture.Udm()))
Calculate the snGrad of the given volField.
const fvSolution & solution() const
Return the fvSchemes.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
GeometricField< vector, fvPatchField, volMesh > volVectorField
CGAL::Exact_predicates_exact_constructions_kernel K
multiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > rho() const
Return the mixture density.
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
dictionary()
Construct top-level dictionary null.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
bool read()
Read base phaseProperties dictionary.
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Calculate the gradient of the given field.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void correct()
Correct the mixture properties.
Calculate the face-flux of the given field.
friend class const_iterator
Declare friendship with the const_iterator.
const Type & value() const
Return const reference to value.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
IOdictionary(const IOobject &io, const word &wantedType)
Construct given an IOobject, supply wanted typeName.
Calculate the divergence of the given field.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
void solve()
Solve for the mixture phase-fractions.
const dimensionSet dimMass
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
word alpharScheme("div(phirb,alpha)")
dimensioned< scalar > mag(const dimensioned< Type > &)
MULES: Multidimensional universal limiter for explicit solution.
const doubleScalar e
Elementary charge.
A class for managing temporary objects.
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return the local objectRegistry.
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
dimensionedScalar deltaT() const
Return time step.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
friend class const_iterator
tmp< surfaceScalarField > surfaceTensionForce() const
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.