Go to the documentation of this file.
32 template<
class BasePhaseSystem>
39 HeatAndMassTransferPhaseSystem<BasePhaseSystem>(
mesh),
40 volatile_(this->
lookup(
"volatile")),
41 saturationModel_(saturationModel::
New(this->subDict(
"saturationModel"))),
42 massTransfer_(this->
lookup(
"massTransfer"))
52 const phasePair& pair(phasePairIter());
68 this->mesh().time().timeName(),
83 template<
class BasePhaseSystem>
91 template<
class BasePhaseSystem>
95 return saturationModel_();
99 template<
class BasePhaseSystem>
103 typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
104 alphatPhaseChangeWallFunction;
106 autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
119 const phasePair& pair(phasePairIter());
126 const phaseModel& phase = pair.phase1();
127 const phaseModel& otherPhase = pair.phase2();
134 phase.mesh().time().timeName(),
148 "alphat." + otherPhase.name()
155 "alphat." + otherPhase.name()
165 isA<alphatPhaseChangeWallFunction>
172 refCast<const alphatPhaseChangeWallFunction>
179 label faceCelli = currPatch.faceCells()[facei];
180 mDotL[faceCelli] = patchMDotL[facei];
186 *eqns[otherPhase.name()] -= mDotL;
194 template<
class BasePhaseSystem>
199 autoPtr<phaseSystem::massTransferTable> eqnsPtr
208 const phaseModel& phase = this->phaseModels_[
phasei];
210 const PtrList<volScalarField>& Yi = phase.Y();
229 const phasePair& pair(phasePairIter());
235 const phaseModel& phase = pair.phase1();
236 const phaseModel& otherPhase = pair.phase2();
253 *eqns[otherName] += dmdt12 -
fvm::Sp(dmdt12, eqns[otherName]->
psi());
260 template<
class BasePhaseSystem>
263 typedef compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
264 alphatPhaseChangeWallFunction;
277 const phasePair& pair(phasePairIter());
284 const phaseModel&
phase1 = pair.phase1();
285 const phaseModel&
phase2 = pair.phase2();
317 this->heatTransferModels_[pair][pair.first()]->K(0)
322 this->heatTransferModels_[pair][pair.second()]->K(0)
327 scalar iDmdtRelax(this->
mesh().fieldRelaxationFactor(
"iDmdt"));
330 (1 - iDmdtRelax)*iDmdt
331 + iDmdtRelax*(H1*(Tf - T1) + H2*(Tf - T2))
335 - (
neg(iDmdt)*he1 +
pos(iDmdt)*hef1),
339 Info<<
"iDmdt." << pair.name()
340 <<
": min = " <<
min(iDmdt.internalField())
341 <<
", mean = " <<
average(iDmdt.internalField())
342 <<
", max = " <<
max(iDmdt.internalField())
351 volScalarField H1(this->heatTransferModels_[pair][pair.first()]->K());
352 volScalarField H2(this->heatTransferModels_[pair][pair.second()]->K());
355 const scalar HLimit = 1
e-4;
366 - (
neg(iDmdt)*he1 +
pos(iDmdt)*hef1)
370 Tf = (H1*T1 + H2*T2 + mDotL)/(H1 + H2);
372 Info<<
"Tf." << pair.name()
373 <<
": min = " <<
min(Tf.internalField())
374 <<
", mean = " <<
average(Tf.internalField())
375 <<
", max = " <<
max(Tf.internalField())
384 this->mesh().time().timeName(),
415 isA<alphatPhaseChangeWallFunction>
422 refCast<const alphatPhaseChangeWallFunction>
429 label faceCelli = currPatch.faceCells()[facei];
430 wDmdt[faceCelli] += patchDmdt[facei];
435 Info<<
"wDmdt." << pair.name()
436 <<
": min = " <<
min(wDmdt.internalField())
437 <<
", mean = " <<
average(wDmdt.internalField())
438 <<
", max = " <<
max(wDmdt.internalField())
443 dmdt = wDmdt + iDmdt;
445 Info<<
"dmdt." << pair.name()
446 <<
": min = " <<
min(dmdt.internalField())
447 <<
", mean = " <<
average(dmdt.internalField())
448 <<
", max = " <<
max(dmdt.internalField())
455 template<
class BasePhaseSystem>
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
virtual ~ThermalPhaseChangePhaseSystem()
Destructor.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
HashPtrTable< fvScalarMatrix, word, string::hash > heatTransferTable
const dimensionSet dimDensity
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
bool read(const char *, int32_t &)
dimensionedScalar posPart(const dimensionedScalar &ds)
virtual const volScalarField & T() const
Temperature [K].
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual bool read()
Read base phaseProperties dictionary.
HashPtrTable< fvScalarMatrix, word, string::hash > massTransferTable
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
PtrList< fvPatch > fvPatchList
container classes for fvPatch
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fvMatrix< scalar > fvScalarMatrix
ThermalPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const word & name() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
const double e
Elementary charge.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
dimensionedScalar negPart(const dimensionedScalar &ds)
const volScalarField & psi
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
const saturationModel & saturation() const
Return the saturationModel.
virtual volScalarField & p()
Pressure [Pa].
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
virtual void correctThermo()
Correct the thermodynamics.
word name(const complex &)
Return a string representation of a complex.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar pos(const dimensionedScalar &ds)