33#include "surfaceInterpolate.H"
39#include "alphaContactAngleFvPatchScalarField.H"
43void Foam::multiphaseMixture::calcAlphas()
48 for (
const phase& ph : phases_)
50 alphas_ += level * ph;
68 "transportProperties",
121 sigmas_(
lookup(
"sigmas")),
122 dimSigma_(1, 0, -2, 0, 0),
129 rhoPhi_.setOriented();
141 auto iter = phases_.cbegin();
143 tmp<volScalarField>
trho = iter()*iter().rho();
146 for (++iter; iter != phases_.cend(); ++iter)
148 rho += iter()*iter().rho();
155Foam::tmp<Foam::scalarField>
158 auto iter = phases_.cbegin();
163 for (++iter; iter != phases_.cend(); ++iter)
165 rho += iter().boundaryField()[patchi]*iter().rho().value();
172Foam::tmp<Foam::volScalarField>
175 auto iter = phases_.cbegin();
180 for (++iter; iter != phases_.cend(); ++iter)
182 mu += iter()*iter().rho()*iter().nu();
189Foam::tmp<Foam::scalarField>
192 auto iter = phases_.cbegin();
196 iter().boundaryField()[patchi]
197 *iter().rho().value()
203 for (++iter; iter != phases_.cend(); ++iter)
207 iter().boundaryField()[patchi]
208 *iter().rho().value()
217Foam::tmp<Foam::surfaceScalarField>
220 auto iter = phases_.cbegin();
226 for (++iter; iter != phases_.cend(); ++iter)
236Foam::tmp<Foam::volScalarField>
243Foam::tmp<Foam::scalarField>
246 return nu_.boundaryField()[patchi];
250Foam::tmp<Foam::surfaceScalarField>
257Foam::tmp<Foam::surfaceScalarField>
262 "surfaceTensionForce",
277 for (++iter2; iter2 != phases_.cend(); ++iter2)
286 <<
"Cannot find interface " << interfacePair(
alpha1,
alpha2)
287 <<
" in list of sigma values"
320 mesh_.newIOobject(
"rhoPhiSum"),
330 !(++alphaSubCycle).
end();
351 for (
phase& ph : phases_)
358Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixture::nHatfv
380 return gradAlphaf/(
mag(gradAlphaf) + deltaN_);
384Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixture::nHatf
401void Foam::multiphaseMixture::correctContactAngle
420 correctBoundaryContactAngle(acap, patchi,
alpha1,
alpha2, nHatb);
427 correctBoundaryContactAngle(acap, patchi,
alpha2,
alpha1, nHatb);
433void Foam::multiphaseMixture::correctBoundaryContactAngle
448 mesh_.Sf().boundaryField()[patchi]
449 /mesh_.magSf().boundaryField()[patchi]
452 const auto tp = acap.thetaProps().cfind(interfacePair(
alpha1,
alpha2));
457 <<
"Cannot find interface " << interfacePair(
alpha1,
alpha2)
458 <<
"\n in table of theta properties for patch "
459 << acap.patch().name()
463 const bool matched = (tp.key().first() ==
alpha1.
name());
465 const scalar theta0 =
degToRad(tp().theta0(matched));
468 const scalar uTheta = tp().uTheta();
473 const scalar thetaA =
degToRad(tp().thetaA(matched));
474 const scalar thetaR =
degToRad(tp().thetaR(matched));
479 U_.boundaryField()[patchi].patchInternalField()
480 - U_.boundaryField()[patchi]
482 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
487 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
491 nWall /= (
mag(nWall) + SMALL);
497 theta += (thetaA - thetaR)*
tanh(uwall/uTheta);
511 b2[facei] =
cos(
acos(a12[facei]) - theta[facei]);
519 nHatPatch = a*AfHatPatch +
b*nHatPatch;
521 nHatPatch /= (
mag(nHatPatch) + deltaN_.value());
525Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::K
536 return -
fvc::div(tnHatfv & mesh_.Sf());
540Foam::tmp<Foam::volScalarField>
551 for (
const phase& ph : phases_)
553 tnearInt.ref() =
max(tnearInt(),
pos0(ph - 0.01)*
pos0(0.99 - ph));
560void Foam::multiphaseMixture::solveAlphas
565 static label nSolves(-1);
568 const word alphaScheme(
"div(phi,alpha)");
612 1.0/mesh_.time().deltaTValue(),
633 mesh_.newIOobject(
"sumAlpha"),
665 Info<<
"Phase-sum volume fraction, min, max = "
666 << sumAlpha.weightedAverage(mesh_.V()).value()
667 <<
' ' <<
min(sumAlpha).value()
668 <<
' ' <<
max(sumAlpha).value()
691 for (
phase& ph : phases_)
693 readOK &= ph.read(phaseData[
phasei++].
dict());
696 readEntry(
"sigmas", sigmas_);
CGAL::Exact_predicates_exact_constructions_kernel K
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & alpha1
const volScalarField & alpha2
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
dimensionedScalar deltaT() const
Return time step.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
T & first()
Access first element of the list, position [0].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
tmp< surfaceScalarField > surfaceTensionForce() const
void correct()
Correct the mixture properties.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
tmp< volScalarField > rho() const
Return the mixture density.
multiphaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
void solve()
Solve for the mixture phase-fractions.
bool read()
Read base transportProperties dictionary.
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
const Time & time() const noexcept
Return time registry.
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Lookup type of boundary radiation properties.
bool end() const
Return true if the number of sub-cycles has been reached.
Perform a subCycleTime on a field.
A class for managing temporary objects.
virtual bool read()=0
Read transportProperties dictionary.
Upwind differencing scheme class.
A class for handling words, derived from Foam::string.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
const expr V(m.psi().mesh().V())
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
Different types of constants.
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.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
tmp< volScalarField > trho
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Unit conversion functions.