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>
162 heatTransferModelIter
167 this->phasePairs_[heatTransferModelIter.key()]
174 heatTransferModelIter().first()->
K(),
175 heatTransferModelIter().second()->
K()
224 if (this->heatTransferModels_.found(phasePairIter.key()))
246template<
class BasePhaseSystem>
250 BasePhaseSystem::correctEnergyTransport();
252 correctInterfaceThermo();
256template<
class BasePhaseSystem>
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);
314template<
class BasePhaseSystem>
317 if (BasePhaseSystem::read())
CGAL::Exact_predicates_exact_constructions_kernel K
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
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, overwriting existing entries.
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.
Class which models interfacial heat transfer between a number of phases. Two heat transfer models are...
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.
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 const volScalarField & T() const
Temperature [K].
virtual volScalarField & p()
Pressure [Pa].
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
A special matrix type and solver, designed for finite volume solutions of scalar equations....
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
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
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"))