30#include "alphaContactAngleFvPatchScalarField.H"
35#include "surfaceInterpolate.H"
46void Foam::multiphaseSystem::calcAlphas()
51 for (
const phaseModel& phase : phases_)
53 alphas_ += level * phase;
59void Foam::multiphaseSystem::solveAlphas()
93 const auto cAlpha = cAlphas_.cfind(interfacePair(
phase,
phase2));
105 const word phirScheme
118 phase.correctInflowOutflow(alphaPhiCorr);
122 1.0/mesh_.time().deltaTValue(),
144 mesh_.time().timeName(),
179 Info<<
"Phase-sum volume fraction, min, max = "
180 << sumAlpha.weightedAverage(mesh_.V()).value()
181 <<
' ' <<
min(sumAlpha).value()
182 <<
' ' <<
max(sumAlpha).value()
197Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseSystem::nHatfv
219 return gradAlphaf/(
mag(gradAlphaf) + deltaN_);
223Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseSystem::nHatf
240void Foam::multiphaseSystem::correctContactAngle
275 mesh_.Sf().boundaryField()[patchi]
276 /mesh_.magSf().boundaryField()[patchi]
285 <<
"Cannot find interface " << interfacePair(
phase1,
phase2)
286 <<
"\n in table of theta properties for patch "
287 << acap.patch().name()
291 bool matched = (tp.key().first() ==
phase1.
name());
293 const scalar theta0 =
degToRad(tp().theta0(matched));
296 scalar uTheta = tp().uTheta();
301 const scalar thetaA =
degToRad(tp().thetaA(matched));
302 const scalar thetaR =
degToRad(tp().thetaR(matched));
310 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
315 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch
319 nWall /= (
mag(nWall) + SMALL);
325 theta += (thetaA - thetaR)*
tanh(uwall/uTheta);
339 b2[facei] =
cos(
acos(a12[facei]) - theta[facei]);
347 nHatPatch = a*AfHatPatch +
b*nHatPatch;
349 nHatPatch /= (
mag(nHatPatch) + deltaN_.value());
355Foam::tmp<Foam::volScalarField> Foam::multiphaseSystem::K
366 return -
fvc::div(tnHatfv & mesh_.Sf());
382 "transportProperties",
409 sigmas_(
lookup(
"sigmas")),
410 dimSigma_(1, 0, -2, 0, 0),
411 cAlphas_(
lookup(
"interfaceCompression")),
412 Cvms_(
lookup(
"virtualMass")),
422 interfaceDictTable dragModelsDict(
lookup(
"drag"));
432 *phases_.lookup(iter.key().first()),
433 *phases_.lookup(iter.key().second())
438 for (
const phaseModel&
phase1 : phases_)
440 for (
const phaseModel&
phase2 : phases_)
449 if (sigmas_.found(key) && !cAlphas_.found(key))
452 <<
"Compression coefficient not specified for phase pair ("
454 <<
") for which a surface tension coefficient is specified"
466 auto iter = phases_.cbegin();
471 for (++iter; iter != phases_.cend(); ++iter)
473 rho += iter()*iter().rho();
483 auto iter = phases_.cbegin();
485 tmp<scalarField>
trho = iter().boundaryField()[patchi]*iter().rho().value();
488 for (++iter; iter != phases_.cend(); ++iter)
490 rho += iter().boundaryField()[patchi]*iter().rho().value();
499 auto iter = phases_.cbegin();
501 tmp<volScalarField> tmu = iter()*(iter().rho()*iter().nu());
504 for (++iter; iter != phases_.cend(); ++iter)
506 mu += iter()*(iter().rho()*iter().nu());
516 auto iter = phases_.cbegin();
518 tmp<scalarField> tmu =
519 iter().boundaryField()[patchi]
520 *(iter().rho().value()*iter().nu().value());
523 for (++iter; iter != phases_.cend(); ++iter)
526 iter().boundaryField()[patchi]
527 *(iter().rho().value()*iter().nu().value());
530 return tmu/
rho(patchi);
536 const phaseModel& phase
547 for (
const phaseModel&
phase2 : phases_)
554 auto iterCvm = Cvms_.cfind(interfacePair(
phase,
phase2));
562 iterCvm = Cvms_.cfind(interfacePair(
phase2,
phase));
599 auto Cvm = Cvms_.cfind(interfacePair(
phase,
phase2));
617 tSvm.ref().boundaryFieldRef();
630 SvmBf[patchi] =
Zero;
645 const multiphaseEuler::dragModel& dm = *iter();
653 dm.phase1()*dm.phase2(),
654 dm.residualPhaseFraction()
660 mag(dm.phase1().U() - dm.phase2().U()),
669 forAll(dm.phase1().phi().boundaryField(), patchi)
675 dm.phase1().phi().boundaryField()[patchi]
683 dragCoeffsPtr().
set(iter.key(), Kptr);
686 return dragCoeffsPtr;
692 const phaseModel& phase,
704 dragModelTable::const_iterator dmIter = dragModels_.begin();
705 dragCoeffFields::const_iterator dcIter =
dragCoeffs.begin();
709 dmIter.good() && dcIter.good();
715 &phase == &dmIter()->
phase1()
716 || &phase == &dmIter()->
phase2()
719 tdragCoeff.ref() += *dcIter();
739 tSurfaceTension.ref().setOriented();
752 tSurfaceTension.ref() +=
762 return tSurfaceTension;
777 for (
const phaseModel& phase : phases_)
789 for (phaseModel& phase : phases_)
843 !(++alphaSubCycle).
end();
888 PtrList<entry> phaseData(lookup(
"phases"));
891 for (phaseModel& phase : phases_)
893 readOK &= phase.read(phaseData[
phasei++].
dict());
896 lookup(
"sigmas") >> sigmas_;
897 lookup(
"interfaceCompression") >> cAlphas_;
898 lookup(
"virtualMass") >> Cvms_;
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, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
label timeIndex() const noexcept
The time index of the field.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
typename parent_type::const_iterator const_iterator
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
const word & name() const noexcept
Return the object name.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
label timeIndex() const noexcept
Return the current time index.
dimensionedScalar deltaT() const
Return time step.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
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.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
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.
const phaseModel & phase2() const
const dimensionedScalar & residualSlip() const
static autoPtr< dragModel > New(const dictionary &dict, const phaseModel &phase1, const phaseModel &phase2)
virtual tmp< volScalarField > K(const volScalarField &Ur) const =0
The drag function K used in the momentum eq.
const phaseModel & phase1() const
const dimensionedScalar & residualPhaseFraction() const
Name pair for the interface.
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
HashPtrTable< volScalarField, interfacePair, interfacePair::symmHash > dragCoeffFields
autoPtr< dragCoeffFields > dragCoeffs() const
Return the drag coefficients for all of the interfaces.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient for a pair.
multiphaseSystem(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > nu() const
Return the mixture laminar viscosity.
tmp< surfaceScalarField > surfaceTension(const phaseModel &phase) const
tmp< volScalarField > rho() const
Return the mixture density.
tmp< volScalarField > dragCoeff(const phaseModel &phase, const dragCoeffFields &dragCoeffs) const
Return the sum of the drag coefficients for the given phase.
tmp< volScalarField > Cvm(const phaseModel &phase) const
Return the virtual-mass coefficient for the given phase.
tmp< volVectorField > Svm(const phaseModel &phase) const
Return the virtual-mass source for the given phase.
void solve()
Solve for the mixture phase-fractions.
bool read()
Read base transportProperties dictionary.
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...
const volVectorField & U() const
const surfaceScalarField & phi() const
const dimensionedScalar & rho() const
const volVectorField & DDtU() const
const word & name() const
const surfaceScalarField & phi() const
Return the mixture flux.
tmp< volVectorField > U() const
Return the mixture velocity.
const fvMesh & mesh() const
Return the mesh.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
void correct()
Correct the phase properties.
const word & name() const
const dimensionedScalar & rho() const
Return const-access to phase1 density.
bool read(const dictionary &phaseDict)
Read base transportProperties dictionary.
Lookup type of boundary radiation properties.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
virtual bool read()
Read object.
A class for managing sub-cycling times.
bool end() const
Return true if the number of sub-cycles has been reached.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
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()))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Area-weighted average a surfaceField creating a volField.
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.
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
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.
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
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)
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.
const dimensionSet dimDensity
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)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
dimensionedScalar acos(const dimensionedScalar &ds)
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.