Go to the documentation of this file.
29 #include "interfaceHeatResistance.H"
43 namespace temperaturePhaseChangeTwoPhaseMixtures
48 temperaturePhaseChangeTwoPhaseMixture,
60 const thermoIncompressibleTwoPhaseMixture&
mixture,
64 temperaturePhaseChangeTwoPhaseMixture(
mixture,
mesh),
143 optionalSubDict(
type() +
"Coeffs").
get<scalar>(
"spread")
158 return Pair<tmp<volScalarField>>
160 (alphalCoeff*mDotc_)/(mixture_.alpha2() + SMALL),
161 -(alphalCoeff*mDote_)/(mixture_.alpha1() + SMALL)
172 min(
max(mixture_.alpha1(), scalar(0)), scalar(1))
177 min(
max(mixture_.alpha2(), scalar(0)), scalar(1))
180 return Pair<tmp<volScalarField>>
182 (mDotc_/(limitedAlpha2 + SMALL)),
183 -(mDote_/(limitedAlpha1 + SMALL))
192 return Pair<tmp<volScalarField>>
194 tmp<volScalarField>(mDotc_),
195 tmp<volScalarField>(mDote_)
204 const twoPhaseMixtureEThermo&
thermo =
205 refCast<const twoPhaseMixtureEThermo>
214 Pair<tmp<volScalarField>> mDotce(mDot());
216 return Pair<tmp<volScalarField>>
218 mDotc_*
pos(TSat -
T.oldTime())/(TSat -
T.oldTime()),
219 -mDote_*
pos(
T.oldTime() - TSat)/(
T.oldTime() - TSat)
231 auto& TSource = tTSource.ref();
233 const twoPhaseMixtureEThermo&
thermo =
234 refCast<const twoPhaseMixtureEThermo>
244 TSource =
fvm::Sp(IHRcoeff,
T) - IHRcoeff*TSat;
259 const twoPhaseMixtureEThermo&
thermo =
260 refCast<const twoPhaseMixtureEThermo>
271 mDotc_ = interfaceArea_*R_*
max(TSat -
T,
T0)/
L;
272 mDote_ = interfaceArea_*R_*
max(
T - TSat,
T0)/
L;
277 scalar rhobyDt = mixture_.rho1().value()/mesh_.time().deltaTValue();
278 scalar maxEvap = mixture_.alpha1()[celli]*rhobyDt;
279 scalar maxCond = -mixture_.alpha2()[celli]*rhobyDt;
280 mDotc_[celli] =
min(
max(mDotc_[celli], maxCond), maxEvap);
298 if (
max(mDotc_) > MDotMin)
311 if (
max(mDote_) > MDotMin)
326 void Foam::temperaturePhaseChangeTwoPhaseMixtures::interfaceHeatResistance::
337 cutCellIso cutCell(mesh_, ap);
339 forAll(interfaceArea_, celli)
341 label status = cutCell.calcSubCell(celli, 0.5);
342 interfaceArea_[celli] = 0;
345 interfaceArea_[celli] =
346 mag(cutCell.faceArea())/mesh_.V()[celli];
359 return Pair<tmp<volScalarField>>
372 optionalSubDict(
type() +
"Coeffs").readEntry(
"R", R_);
373 optionalSubDict(
type() +
"Coeffs").readEntry(
"spread", spread_);
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
virtual Pair< tmp< volScalarField > > vDot() const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const vector L(dict.get< vector >("L"))
interfaceHeatResistance(const thermoIncompressibleTwoPhaseMixture &mixture, const fvMesh &mesh)
virtual Pair< tmp< volScalarField > > mDot() const
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
const dimensionSet dimEnergy
const dimensionSet dimDensity
Type gAverage(const FieldField< Field, Type > &f)
const volScalarField & alpha2
static const word dictName
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
virtual Pair< tmp< volScalarField > > vDotAlphal() const
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const volScalarField & alpha1
virtual tmp< fvScalarMatrix > TSource() const
void spreadSource(volScalarField &mDotOut, const volScalarField &mDotIn, const volScalarField &alpha1, const volScalarField &alpha2, const dimensionedScalar &D, const scalar cutoff)
label min(const labelHashSet &set, label minValue=labelMax)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Interface Heat Resistance type of condensation/saturation model using spread source distribution foll...
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
const dimensionSet dimArea(sqr(dimLength))
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
const dimensionSet dimPower
virtual Pair< tmp< volScalarField > > mDotAlphal() const
Macros for easy insertion into run-time selection tables.
Calculate the matrix for implicit and explicit sources.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An ordered pair of two objects of type <T> with first() and second() elements.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
const dimensionedScalar e
virtual Pair< tmp< volScalarField > > mDotDeltaT() const
static tmp< T > New(Args &&... args)
const dimensionedScalar & D
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
defineTypeNameAndDebug(combustionModel, 0)
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
const dimensionSet dimless
dimensionedScalar pos(const dimensionedScalar &ds)