30#include "alphaContactAngleFvPatchScalarField.H"
39#include "surfaceInterpolate.H"
52void Foam::multiphaseMixtureThermo::calcAlphas()
57 for (
const phaseModel& phase : phases_)
59 alphas_ += level * phase;
109 sigmas_(
lookup(
"sigmas")),
110 dimSigma_(1, 0, -2, 0, 0),
117 rhoPhi_.setOriented();
128 for (phaseModel& phase : phases_)
133 auto phasei = phases_.cbegin();
159 auto phasei = phases_.cbegin();
176 if (!
phase.thermo().incompressible())
190 if (!
phase.thermo().isochoric())
260 phasei().boundaryField()[patchi]*
phasei().thermo().he(
p,
T, patchi);
269 auto phasei = phases_.cbegin();
310 auto phasei = phases_.cbegin();
328 auto phasei = phases_.cbegin();
338 phasei().boundaryField()[patchi]*
phasei().thermo().rho(patchi);
347 auto phasei = phases_.cbegin();
367 auto phasei = phases_.cbegin();
377 phasei().boundaryField()[patchi]*
phasei().thermo().Cp(
p,
T, patchi);
386 auto phasei = phases_.cbegin();
406 auto phasei = phases_.cbegin();
416 phasei().boundaryField()[patchi]*
phasei().thermo().Cv(
p,
T, patchi);
425 auto phasei = phases_.cbegin();
445 auto phasei = phases_.cbegin();
455 phasei().boundaryField()[patchi]
456 *
phasei().thermo().gamma(
p,
T, patchi);
465 auto phasei = phases_.cbegin();
485 auto phasei = phases_.cbegin();
495 phasei().boundaryField()[patchi]
496 *
phasei().thermo().Cpv(
p,
T, patchi);
505 auto phasei = phases_.cbegin();
525 auto phasei = phases_.cbegin();
535 phasei().boundaryField()[patchi]
536 *
phasei().thermo().CpByCpv(
p,
T, patchi);
545 auto phasei = phases_.cbegin();
569 return mu(patchi)/
rho(patchi);
575 auto phasei = phases_.cbegin();
593 auto phasei = phases_.cbegin();
603 phasei().boundaryField()[patchi]*
phasei().thermo().kappa(patchi);
618 talphaEff.ref() +=
phasei()*
phasei().thermo().alphahe();
634 phasei().boundaryField()[patchi]
641 phasei().boundaryField()[patchi]
642 *
phasei().thermo().alphahe(patchi);
654 auto phasei = phases_.cbegin();
660 tkappaEff.ref() +=
phasei()*
phasei().thermo().kappaEff(alphat);
673 auto phasei = phases_.cbegin();
677 phasei().boundaryField()[patchi]
684 phasei().boundaryField()[patchi]
685 *
phasei().thermo().kappaEff(alphat, patchi);
697 auto phasei = phases_.cbegin();
703 talphaEff.ref() +=
phasei()*
phasei().thermo().alphaEff(alphat);
716 auto phasei = phases_.cbegin();
720 phasei().boundaryField()[patchi]
727 phasei().boundaryField()[patchi]
728 *
phasei().thermo().alphaEff(alphat, patchi);
737 auto phasei = phases_.cbegin();
750Foam::tmp<Foam::surfaceScalarField>
759 "surfaceTensionForce",
760 mesh_.time().timeName(),
786 <<
"Cannot find interface " << interfacePair(
alpha1,
alpha2)
787 <<
" in list of sigma values"
822 !(++alphaSubCycle).
end();
838Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixtureThermo::nHatfv
860 return gradAlphaf/(
mag(gradAlphaf) + deltaN_);
864Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixtureThermo::nHatf
881void Foam::multiphaseMixtureThermo::correctContactAngle
904 mesh_.Sf().boundaryField()[patchi]
905 /mesh_.magSf().boundaryField()[patchi]
914 <<
"Cannot find interface " << interfacePair(
alpha1,
alpha2)
915 <<
"\n in table of theta properties for patch "
916 << acap.patch().name()
920 const bool matched = (tp.key().first() ==
alpha1.
name());
922 const scalar theta0 =
degToRad(tp().theta0(matched));
925 const scalar uTheta = tp().uTheta();
930 const scalar thetaA =
degToRad(tp().thetaA(matched));
931 const scalar thetaR =
degToRad(tp().thetaR(matched));
936 U_.boundaryField()[patchi].patchInternalField()
937 - U_.boundaryField()[patchi]
939 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
944 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
948 nWall /= (
mag(nWall) + SMALL);
954 theta += (thetaA - thetaR)*
tanh(uwall/uTheta);
968 b2[facei] =
cos(
acos(a12[facei]) - theta[facei]);
976 nHatPatch = a*AfHatPatch +
b*nHatPatch;
978 nHatPatch /= (
mag(nHatPatch) + deltaN_.value());
984Foam::tmp<Foam::volScalarField> Foam::multiphaseMixtureThermo::K
995 return -
fvc::div(tnHatfv & mesh_.Sf());
999Foam::tmp<Foam::volScalarField>
1020void Foam::multiphaseMixtureThermo::solveAlphas
1025 static label nSolves(-1);
1028 const word alphaScheme(
"div(phi,alpha)");
1072 1.0/mesh_.time().deltaTValue(),
1096 mesh_.time().timeName(),
1119 mesh_.time().timeName(),
1131 mesh_.time().timeName(),
1144 if (dgdt[celli] < 0.0 &&
alpha[celli] > 0.0)
1146 Sp[celli] += dgdt[celli]*
alpha[celli];
1147 Su[celli] -= dgdt[celli]*
alpha[celli];
1149 else if (dgdt[celli] > 0.0 &&
alpha[celli] < 1.0)
1151 Sp[celli] -= dgdt[celli]*(1.0 -
alpha[celli]);
1164 if (dgdt2[celli] > 0.0 &&
alpha2[celli] < 1.0)
1166 Sp[celli] -= dgdt2[celli]*(1.0 -
alpha2[celli]);
1167 Su[celli] += dgdt2[celli]*
alpha[celli];
1169 else if (dgdt2[celli] < 0.0 &&
alpha2[celli] > 0.0)
1171 Sp[celli] += dgdt2[celli]*
alpha2[celli];
1198 Info<<
"Phase-sum volume fraction, min, max = "
1199 << sumAlpha.weightedAverage(mesh_.V()).value()
1200 <<
' ' <<
min(sumAlpha).value()
1201 <<
' ' <<
max(sumAlpha).value()
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
const dimensionSet & dimensions() const noexcept
Return dimensions.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
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.
@ 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.
An STL-conforming const_iterator.
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].
const_iterator cend() const noexcept
Return const_iterator to end traversing the constant UList.
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
volScalarField alpha_
Laminar thermal diffusivity [kg/m/s].
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.
virtual word thermoName() const
Return the name of the thermo physics.
tmp< surfaceScalarField > surfaceTensionForce() const
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
virtual bool incompressible() const
Return true if the equation of state is incompressible.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
virtual void correct()
Update properties.
multiphaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
void solve()
Solve for the mixture phase-fractions.
virtual bool isochoric() const
Return true if the equation of state is isochoric.
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...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const dimensionedScalar & rho() const
Return const-access to phase1 density.
Basic thermodynamic properties based on compressibility.
volScalarField mu_
Dynamic viscosity [kg/m/s].
volScalarField psi_
Compressibility [s^2/m^2].
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.
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 ...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const tmp< volScalarField > & tCv
const volScalarField & Cv
const tmp< volScalarField > & tCp
const volScalarField & Cp
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
#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 mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the snGrad of the given volField.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
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)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
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)
List< label > labelList
A List of labels.
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).
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
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)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
tmp< volScalarField > trho
psiReactionThermo & thermo
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.