31#include "BlendedInterfacialModel.H"
33#include "virtualMassModel.H"
35#include "wallLubricationModel.H"
36#include "turbulentDispersionModel.H"
53template<
class BasePhaseSystem>
55Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::Kd
57 const phasePairKey& key
60 if (dragModels_.found(key))
62 return dragModels_[
key]->K();
65 return volScalarField::New
67 IOobject::scopedName(dragModel::typeName,
"K"),
68 IOobject::NO_REGISTER,
75template<
class BasePhaseSystem>
77Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::Kdf
79 const phasePairKey& key
82 if (dragModels_.found(key))
84 return dragModels_[
key]->Kf();
87 return surfaceScalarField::New
89 IOobject::scopedName(dragModel::typeName,
"K"),
90 IOobject::NO_REGISTER,
97template<
class BasePhaseSystem>
99Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::Vm
101 const phasePairKey& key
104 if (virtualMassModels_.found(key))
106 return virtualMassModels_[
key]->K();
109 return volScalarField::New
111 IOobject::scopedName(virtualMassModel::typeName,
"K"),
112 IOobject::NO_REGISTER,
119template<
class BasePhaseSystem>
125 phaseSystem::phasePairTable,
130 const phasePair& pair(phasePairIter());
143 if (!pair.phase1().stationary())
148 eqn += dmdt21*pair.phase2().U() - fvm::Sp(dmdt21, eqn.psi());
151 if (!pair.phase2().stationary())
156 eqn -= dmdt12*pair.phase1().U() - fvm::Sp(dmdt12, eqn.psi());
164template<
class BasePhaseSystem>
171 BasePhaseSystem(
mesh)
173 this->generatePairsAndSubModels
179 this->generatePairsAndSubModels
185 this->generatePairsAndSubModels
191 this->generatePairsAndSubModels
194 wallLubricationModels_
197 this->generatePairsAndSubModels
199 "turbulentDispersion",
200 turbulentDispersionModels_
210 const phasePair& pair(this->phasePairs_[dragModelIter.key()]);
228 dragModelIter()->Kf()
240 const phasePair& pair(this->phasePairs_[virtualMassModelIter.key()]);
248 virtualMassModelIter()->
K()
258 virtualMassModelIter()->Kf()
267template<
class BasePhaseSystem>
275template<
class BasePhaseSystem>
287 forAll(this->movingPhases(), movingPhasei)
306 *Kds_[dragModelIter.key()] = dragModelIter()->K();
307 *Kdfs_[dragModelIter.key()] = dragModelIter()->Kf();
314 const phasePair& pair(this->phasePairs_[KdIter.key()]);
318 if (!iter().stationary())
330 virtualMassModelTable,
335 *Vms_[virtualMassModelIter.key()] = virtualMassModelIter()->K();
336 *Vmfs_[virtualMassModelIter.key()] = virtualMassModelIter()->Kf();
343 const phasePair& pair(this->phasePairs_[VmIter.key()]);
350 if (!
phase.stationary())
366 + this->MRF_.DDt(Vm,
U - otherPhase.U());
372 addMassTransferMomentumTransfer(eqns);
378template<
class BasePhaseSystem>
390 forAll(this->movingPhases(), movingPhasei)
402 PtrList<fvVectorMatrix> UgradUs(this->phaseModels_.size());
405 const phaseModel& phase = this->phaseModels_[
phasei];
407 if (!phase.stationary())
409 const tmp<volVectorField> tU(phase.U());
429 const phasePair& pair(this->phasePairs_[VmIter.key()]);
436 if (!
phase.stationary())
441 UgradUs[
phase.index()]
442 - (UgradUs[otherPhase.index()] & otherPhase.U())
449 addMassTransferMomentumTransfer(eqns);
455template<
class BasePhaseSystem>
465 const phasePair& pair(this->phasePairs_[KdfIter.key()]);
469 this->addField(iter(),
"AFf", Kf,
AFfs);
477 const phasePair& pair(this->phasePairs_[VmfIter.key()]);
485 if (this->fillFields_)
494template<
class BasePhaseSystem>
512 const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
533 wallLubricationModelTable,
534 wallLubricationModels_,
535 wallLubricationModelIter
540 pair(this->phasePairs_[wallLubricationModelIter.key()]);
576 const bool implicitPhasePressure =
577 this->mesh_.solverDict(
phase.volScalarField::name()).
578 template getOrDefault<Switch>
580 "implicitPhasePressure",
584 if (implicitPhasePressure)
586 this->addField(
phase,
"DByAf", pPrimeByAf, DByAfs_);
593 turbulentDispersionModelTable,
594 turbulentDispersionModels_,
595 turbulentDispersionModelIter
599 pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
620 if (DByAfs_.found(pair.phase1().name()))
622 this->addField(pair.phase1(),
"DByAf", DByA1f, DByAfs_);
626 if (this->fillFields_)
635template<
class BasePhaseSystem>
648 const phasePair& pair(this->phasePairs_[VmfIter.key()]);
659 + iter.otherPhase().DUDtf()
675 const phasePair& pair(this->phasePairs_[liftModelIter.key()]);
681 rAUfs[pair.phase1().index()]*
Ff,
688 -
rAUfs[pair.phase2().index()]*
Ff,
696 wallLubricationModelTable,
697 wallLubricationModels_,
698 wallLubricationModelIter
703 pair(this->phasePairs_[wallLubricationModelIter.key()]);
709 rAUfs[pair.phase1().index()]*
Ff,
716 -
rAUfs[pair.phase2().index()]*
Ff,
739 const bool implicitPhasePressure =
740 this->mesh_.solverDict(
phase.volScalarField::name()).
741 template getOrDefault<Switch>
743 "implicitPhasePressure",
747 if (implicitPhasePressure)
749 this->addField(
phase,
"DByAf", pPrimeByAf, DByAfs_);
756 turbulentDispersionModelTable,
757 turbulentDispersionModels_,
758 turbulentDispersionModelIter
762 pair(this->phasePairs_[turbulentDispersionModelIter.key()]);
783 if (DByAfs_.found(pair.phase1().name()))
785 this->addField(pair.phase1(),
"DByAf", DByAf1, DByAfs_);
789 if (this->fillFields_)
798template<
class BasePhaseSystem>
811 const phasePair& pair(this->phasePairs_[KdIter.key()]);
820 *this->
MRF().absolute(iter.otherPhase().phi()),
826 if (this->fillFields_)
840template<
class BasePhaseSystem>
853 const phasePair& pair(this->phasePairs_[KdfIter.key()]);
861 -
rAUfs[iter().index()]*Kf
862 *this->
MRF().absolute(iter.otherPhase().phi()),
868 if (this->fillFields_)
882template<
class BasePhaseSystem>
895 const phasePair& pair(this->phasePairs_[KdIter.key()]);
903 -
rAUs[iter().index()]*
K*iter.otherPhase().U(),
909 if (this->fillFields_)
918template<
class BasePhaseSystem>
923 const bool includeVirtualMass
945 const phaseModel& phase = this->phaseModels_[
phasei];
954 tmp<surfaceScalarField> phiCorrCoeff =
pos0(alphafBar - 0.99);
957 phiCorrCoeff.ref().boundaryFieldRef();
959 forAll(this->mesh_.boundary(), patchi)
964 !this->mesh_.boundary()[patchi].coupled()
968 phiCorrCoeffBf[patchi] = 0;
985 if (includeVirtualMass)
990 const phasePair& pair(this->phasePairs_[VmIter.key()]);
1003 phiCorrs[
phase.index()]
1004 + this->MRF().absolute(otherPhase.phi())
1006 - phiCorrs[otherPhase.index()]
1018template<
class BasePhaseSystem>
1021 const PtrList<volScalarField>&
rAUs
1024 Info<<
"Inverting drag systems: ";
1029 PtrList<PtrList<volScalarField>> KdByAs(
phases.size());
1030 PtrList<PtrList<surfaceScalarField>> phiKds(
phases.size());
1037 new PtrList<volScalarField>(
phases.size())
1043 new PtrList<surfaceScalarField>(
phases.size())
1050 const phasePair& pair(this->phasePairs_[KdIter.key()]);
1052 const label phase1i = pair.phase1().index();
1053 const label phase2i = pair.phase2().index();
1096 for (label i = 0; i <
phases.size(); ++ i)
1098 for (label j = i + 1; j <
phases.size(); ++ j)
1100 KdByAs[i][j] /= KdByAs[i][i];
1101 phiKds[i][j] /= phiKds[i][i];
1102 for (label
k = i + 1;
k <
phases.size(); ++
k)
1104 KdByAs[j][
k] -= KdByAs[j][i]*KdByAs[i][
k];
1105 phiKds[j][
k] -= phiKds[j][i]*phiKds[i][
k];
1112 for (label i = 1; i <
phases.size(); ++ i)
1114 detKdByAs *= KdByAs[i][i];
1115 detPhiKdfs *= phiKds[i][i];
1117 Info<<
"Min cell/face det = " <<
gMin(detKdByAs.primitiveField())
1118 <<
"/" <<
gMin(detPhiKdfs.primitiveField()) <<
endl;
1122 for (label i = 1; i <
phases.size(); ++ i)
1124 if (!
phases[i].stationary())
1126 for (label j = 0; j < i; j ++)
1133 for (label i =
phases.size() - 1; i >= 0; i --)
1135 if (!
phases[i].stationary())
1137 for (label j =
phases.size() - 1; j > i; j --)
1142 phases[i].URef() /= KdByAs[i][i];
1149template<
class BasePhaseSystem>
1155 Info<<
"Inverting drag system: ";
1174 const phasePair& pair(this->phasePairs_[KdfIter.key()]);
1176 const label phase1i = pair.phase1().index();
1177 const label phase2i = pair.phase2().index();
1203 for (label i = 0; i <
phases.size(); ++ i)
1205 for (label j = i + 1; j <
phases.size(); ++ j)
1207 phiKdfs[i][j] /= phiKdfs[i][i];
1208 for (label
k = i + 1;
k <
phases.size(); ++
k)
1210 phiKdfs[j][
k] -= phiKdfs[j][i]*phiKdfs[i][
k];
1216 for (label i = 1; i <
phases.size(); ++ i)
1218 detPhiKdfs *= phiKdfs[i][i];
1220 Info<<
"Min face det = " <<
gMin(detPhiKdfs.primitiveField()) <<
endl;
1224 for (label i = 1; i <
phases.size(); ++ i)
1226 if (!
phases[i].stationary())
1228 for (label j = 0; j < i; j ++)
1234 for (label i =
phases.size() - 1; i >= 0; i --)
1236 if (!
phases[i].stationary())
1238 for (label j =
phases.size() - 1; j > i; j --)
1242 phases[i].phiRef() /= phiKdfs[i][i];
1248template<
class BasePhaseSystem>
1249const Foam::HashPtrTable<Foam::surfaceScalarField>&
1256template<
class BasePhaseSystem>
1259 if (BasePhaseSystem::read())
CGAL::Exact_predicates_exact_constructions_kernel K
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Class which models interfacial momenum transfer between a number of phases. Drag, virtual mass,...
virtual PtrList< surfaceScalarField > AFfs() const
Return implicit force coefficients on the faces, for the face-based.
HashTable< autoPtr< BlendedInterfacialModel< dragModel > >, phasePairKey, phasePairKey::hash > dragModelTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > KdTable
virtual PtrList< surfaceScalarField > phiFs(const PtrList< volScalarField > &rAUs)
Return the explicit force fluxes for the cell-based algorithm, that.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual PtrList< surfaceScalarField > ddtCorrByAs(const PtrList< volScalarField > &rAUs, const bool includeVirtualMass=false) const
Return the flux corrections for the cell-based algorithm. These.
virtual void partialElimination(const PtrList< volScalarField > &rAUs)
Solve the drag system for the velocities and fluxes.
virtual PtrList< surfaceScalarField > phiKdPhifs(const PtrList< surfaceScalarField > &rAUfs) const
As phiKdPhis, but for the face-based algorithm.
HashTable< autoPtr< BlendedInterfacialModel< turbulentDispersionModel > >, phasePairKey, phasePairKey::hash > turbulentDispersionModelTable
HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > VmfTable
HashPtrTable< surfaceScalarField, phasePairKey, phasePairKey::hash > KdfTable
virtual ~MomentumTransferPhaseSystem()
Destructor.
HashTable< autoPtr< BlendedInterfacialModel< virtualMassModel > >, phasePairKey, phasePairKey::hash > virtualMassModelTable
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > VmTable
virtual PtrList< volVectorField > KdUByAs(const PtrList< volScalarField > &rAUs) const
Return the explicit part of the drag force for the cell-based.
virtual const HashPtrTable< surfaceScalarField > & DByAfs() const
Return the phase diffusivities divided by the momentum coefficients.
virtual void partialEliminationf(const PtrList< surfaceScalarField > &rAUfs)
As partialElimination, but for the face-based algorithm. Only solves.
virtual PtrList< surfaceScalarField > phiFfs(const PtrList< surfaceScalarField > &rAUfs)
As phiFs, but for the face-based algorithm.
MomentumTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
As momentumTransfer, but for the face-based algorithm.
virtual PtrList< surfaceScalarField > phiKdPhis(const PtrList< volScalarField > &rAUs) const
Return the explicit drag force fluxes for the cell-based algorithm.
HashTable< autoPtr< BlendedInterfacialModel< wallLubricationModel > >, phasePairKey, phasePairKey::hash > wallLubricationModelTable
HashTable< autoPtr< BlendedInterfacialModel< liftModel > >, phasePairKey, phasePairKey::hash > liftModelTable
virtual bool read()
Read base phaseProperties dictionary.
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,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Mesh data needed to do the Finite Volume discretisation.
const word & name() const
The name of this phase.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const volVectorField & U() const
const surfaceScalarField & phi() const
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
virtual tmp< volVectorField > DUDt() const =0
Return the substantive acceleration.
label index() const
Return the index of the phase.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
virtual word name() const
Pair name.
const multiphaseInter::phaseModel & phase1() const
const multiphaseInter::phaseModel & phase2() const
PtrListDictionary< phaseModel > phaseModelList
HashPtrTable< fvVectorMatrix > momentumTransferTable
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const word & name() const
const dimensionedScalar & rho() const
Return const-access to phase1 density.
A class for managing temporary objects.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime(); *Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Area-weighted average a surfaceField creating a volField.
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 matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the finiteVolume matrix for implicit and explicit sources.
PtrList< surfaceScalarField > phiFs(fluid.phiFs(rAUs))
fluid fillFields("rAU", dimTime/dimDensity, rAUs)
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
PtrList< volScalarField > rAUs
PtrList< surfaceScalarField > rAUfs
PtrList< surfaceScalarField > AFfs(fluid.AFfs())
volVectorField F(fluid.F())
surfaceScalarField snGradAlpha1(fvc::snGrad(alpha1) *mesh.magSf())
surfaceScalarField Vmf("Vmf", fluid.Vmf())
surfaceScalarField Ff(fluid.Ff())
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, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
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< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< volScalarField > byDt(const volScalarField &vf)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar negPart(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimForce
Type gMin(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensionedScalar posPart(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
const dimensionedScalar & D
multiphaseSystem::phaseModelList & phases
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.