31#include "virtualMassModel.H"
64 phase1_(phaseModels_[0]),
65 phase2_(phaseModels_[1])
134 word alphaScheme(
"div(phi," +
alpha1.name() +
')');
145 if (phase1_.divU().valid() && phase2_.divU().valid())
155 else if (phase1_.divU().valid())
163 else if (phase2_.divU().valid())
177 if (DByAfs().
found(phase1_.name()) && DByAfs().
found(phase2_.name()))
181 *DByAfs()[phase1_.name()] + *DByAfs()[phase2_.name()]
196 mesh_.newIOobject(
"Sp"),
203 mesh_.newIOobject(
"Su"),
215 if (dgdt[celli] > 0.0)
220 else if (dgdt[celli] < 0.0)
258 !(++alphaSubCycle).
end();
269 (alphaSubCycle.index()*
Sp)(),
270 (
Su - (alphaSubCycle.index() - 1)*
Sp*
alpha1)(),
275 if (alphaSubCycle.index() == 1)
304 if (alphaDbyA.valid())
315 phase1_.alphaPhiRef() += alpha1Eqn.flux();
318 phase1_.alphaRhoPhiRef() =
321 phase2_.alphaPhiRef() = phi_ - phase1_.alphaPhi();
322 phase2_.correctInflowOutflow(phase2_.alphaPhiRef());
323 phase2_.alphaRhoPhiRef() =
328 <<
" Min(alpha1) = " <<
min(
alpha1).value()
329 <<
" Max(alpha1) = " <<
max(
alpha1).value()
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & alpha1
surfaceScalarField & phi2
const volScalarField & alpha2
surfaceScalarField & phi1
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weights, const label comm=UPstream::worldComm) const
Return the global weighted average.
DimensionedField< scalar, volMesh > Internal
void clamp_range(const dimensioned< MinMax< Type > > &range)
Clamp field values (in-place) to the specified range.
void correctBoundaryConditions()
Correct boundary field.
const word & name() const noexcept
Return the object name.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
An ordered or unorder pair of phase names. Typically specified as follows.
Class to represent a system of phases and model interfacial transfers between them.
phaseModelList phaseModels_
Phase models.
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const =0
Return the phase diffusivities divided by the momentum.
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
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.
bool valid() const noexcept
Identical to good(), or bool operator.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Class which solves the volume fraction equations for two phases.
tmp< volScalarField > Kd() const
Return the drag coefficient.
virtual ~twoPhaseSystem()
Destructor.
phaseModel & phase1_
Phase model 1.
tmp< surfaceScalarField > Kdf() const
Return the face drag coefficient.
const volScalarField & dgdt() const
Return the dilatation term.
phaseModel & phase2_
Phase model 2.
const fvMesh & mesh() const
Return the mesh.
tmp< volScalarField > sigma() const
Return the surface tension coefficient.
twoPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual void solve()
Solve for the phase fractions.
tmp< volScalarField > Vm() const
Return the virtual mass coefficient.
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.
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
word alpharScheme("div(phirb,alpha)")
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the snGrad of the given volField.
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
const volScalarField Kd(fluid.Kd())
const surfaceScalarField Kdf("Kdf", fluid.Kdf())
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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< 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< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
fvMatrix< scalar > fvScalarMatrix
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0).
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.