30#include "BlendedInterfacialModel.H"
31#include "heatTransferModel.H"
42template<
class BasePhaseSystem>
51 this->generatePairsAndSubModels
66 const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
71 <<
"A heat transfer model for the " << pair.
phase1().
name()
72 <<
" side of the " << pair <<
" pair is not specified"
78 <<
"A heat transfer model for the " << pair.
phase2().
name()
79 <<
" side of the " << pair <<
" pair is not specified"
92 const phasePair& pair = this->phasePairs_[heatTransferModelIter.key()];
112 this->mesh().time().timeName(),
117 (H1*T1 + H2*T2)/
max(H1 + H2, HSmall)
126template<
class BasePhaseSystem>
134template<
class BasePhaseSystem>
160 heatTransferModelTable,
162 heatTransferModelIter
165 const phasePair& pair
167 this->phasePairs_[heatTransferModelIter.key()]
172 const Pair<tmp<volScalarField>> Ks
174 heatTransferModelIter().first()->
K(),
175 heatTransferModelIter().second()->
K()
180 const phaseModel& phase = iter();
186 *eqns[phase.name()] +=
187 K*(Tf - phase.thermo().
T())
224 if (this->heatTransferModels_.found(phasePairIter.key()))
246template<
class BasePhaseSystem>
250 BasePhaseSystem::correctEnergyTransport();
256template<
class BasePhaseSystem>
262 heatTransferModelTable,
264 heatTransferModelIter
269 this->phasePairs_[heatTransferModelIter.key()]
291 this->heatTransferModels_[pair].first()->
K()
296 this->heatTransferModels_[pair].second()->
K()
303 Tf = (H1*T1 + H2*T2 + dmdt*
L)/(H1 + H2);
305 Info<<
"Tf." << pair.name()
306 <<
": min = " <<
min(Tf.primitiveField())
307 <<
", mean = " <<
average(Tf.primitiveField())
308 <<
", max = " <<
max(Tf.primitiveField())
314template<
class BasePhaseSystem>
317 if (BasePhaseSystem::read())
CGAL::Exact_predicates_exact_constructions_kernel K
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place).
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An ordered pair of two objects of type <T> with first() and second() elements.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat and Tf.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
HashTable< Pair< autoPtr< BlendedInterfacialModel< heatTransferModel > > >, phasePairKey, phasePairKey::hash > heatTransferModelTable
TwoResistanceHeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
heatTransferModelTable heatTransferModels_
Heat transfer models.
virtual ~TwoResistanceHeatTransferPhaseSystem()
Destructor.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > Tf_
Interface temperatures.
virtual bool read()
Read base phaseProperties dictionary.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
virtual volScalarField & p()
Pressure [Pa].
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
Mesh data needed to do the Finite Volume discretisation.
static const dimensionSet dimK
Coefficient dimensions.
const word & name() const
The name of this phase.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
virtual tmp< volScalarField > K() const =0
Return the phase kinetic energy.
const word & name() const
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
bool ordered() const noexcept
Return the ordered flag.
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
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
HashPtrTable< fvScalarMatrix > heatTransferTable
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const word & name() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the divergence of the given field.
Calculate the finiteVolume matrix for implicit and explicit sources.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimEnergy
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).
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar negPart(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
dimensionedScalar posPart(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
const vector L(dict.get< vector >("L"))