33 template<
class CloudType>
37 return cloudCopyPtr_();
41 template<
class CloudType>
42 inline const typename CloudType::particleType::constantProperties&
49 template<
class CloudType>
50 inline typename CloudType::particleType::constantProperties&
57 template<
class CloudType>
61 return carrierThermo_;
65 template<
class CloudType>
72 template<
class CloudType>
79 template<
class CloudType>
86 template<
class CloudType>
90 return heatTransferModel_;
94 template<
class CloudType>
98 return compositionModel_;
102 template<
class CloudType>
110 template<
class CloudType>
117 template<
class CloudType>
124 <<
"Radiation field requested, but radiation model not active"
132 template<
class CloudType>
139 <<
"Radiation field requested, but radiation model not active"
147 template<
class CloudType>
154 <<
"Radiation field requested, but radiation model not active"
162 template<
class CloudType>
169 <<
"Radiation field requested, but radiation model not active"
177 template<
class CloudType>
184 <<
"Radiation field requested, but radiation model not active"
188 return radAreaPT4_();
192 template<
class CloudType>
199 <<
"Radiation field requested, but radiation model not active"
203 return radAreaPT4_();
207 template<
class CloudType>
215 template<
class CloudType>
223 template<
class CloudType>
231 template<
class CloudType>
239 template<
class CloudType>
245 Info<<
"hsTrans min/max = " <<
min(hsTrans()).value() <<
", "
246 <<
max(hsTrans()).value() <<
nl
247 <<
"hsCoeff min/max = " <<
min(hsCoeff()).value() <<
", "
248 <<
max(hsCoeff()).value() <<
endl;
253 if (this->
solution().semiImplicit(
"h"))
257 Vdt(this->mesh().V()*this->db().time().deltaT());
264 + (hsCoeff()/Vdt)*
hs;
271 + hsCoeff()/(
Cp*Vdt)*
hs;
279 fvm.
source() = -hsTrans()/(this->db().time().deltaT());
289 template<
class CloudType>
296 this->
name() +
":radiation:Ep",
305 const scalar dt = this->db().time().deltaTValue();
307 const scalar
epsilon = constProps_.epsilon0();
308 const scalarField& sumAreaPT4 = radAreaPT4_->primitiveField();
317 template<
class CloudType>
324 this->
name() +
":radiation:ap",
333 const scalar dt = this->db().time().deltaTValue();
335 const scalar
epsilon = constProps_.epsilon0();
336 const scalarField& sumAreaP = radAreaP_->primitiveField();
345 template<
class CloudType>
353 this->
name() +
":radiation:sigmap",
362 const scalar dt = this->db().time().deltaTValue();
364 const scalar
epsilon = constProps_.epsilon0();
365 const scalar
f = constProps_.f0();
366 const scalarField& sumAreaP = radAreaP_->primitiveField();
368 sigmap = sumAreaP*(1.0 -
f)*(1.0 -
epsilon)/V/dt;
375 template<
class CloudType>
401 template<
class CloudType>
scalar hs(const scalar p, const scalar T) const
scalar Cp(const scalar p, const scalar T) const
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Templated heat transfer model class.
Templated base class for thermodynamic cloud.
const parcelType::constantProperties & constProps() const
Return the constant properties.
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/kg/m^3/s].
volScalarField::Internal & radAreaP()
Radiation sum of parcel projected areas [m^2].
volScalarField::Internal & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
const parcelThermo & thermo() const
Return const access to thermo package.
const volScalarField & T() const
Return const access to the carrier temperature field.
bool radiation() const
Radiation flag.
scalar Tmin() const
Minimum temperature.
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
tmp< volScalarField::Internal > hsTrans() const
Return sensible enthalpy transfer [J/kg].
const integrationScheme & TIntegrator() const
Return reference to velocity integration.
volScalarField::Internal & hsTransRef()
Access sensible enthalpy transfer [J/kg].
const ThermoCloud & cloudCopy() const
Return a reference to the cloud copy.
scalar Tmax() const
Maximum temperature.
tmp< volScalarField > sigmap() const
Return tmp equivalent particulate scattering factor.
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
const CompositionModel< ThermoCloud< CloudType > > & composition() const
Return const access to reacting composition model.
const fluidThermo & carrierThermo() const
Return const access to carrier thermo package.
const volScalarField & p() const
Return const access to the carrier pressure field.
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
tmp< volScalarField::Internal > hsCoeff() const
Return coefficient for carrier phase hs equation.
volScalarField::Internal & hsCoeffRef()
Access coefficient for carrier phase hs equation.
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
const Type & value() const
Return const reference to value.
Base-class for fluid thermodynamic properties.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Base for a set of schemes which integrate simple ODEs which arise from semi-implicit rate expressions...
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Selector class for relaxation factors, solver type and solution.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the matrix for implicit and explicit sources.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionSet dimEnergy
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
fvMatrix< scalar > fvScalarMatrix
errorManip< error > abort(error &err)
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.