56 class integrationScheme;
58 template<
class CloudType>
61 template<
class CloudType>
62 class HeatTransferModel;
64 template<
class CloudType>
65 class CompositionModel;
79 template<
class CloudType>
83 public ThermoCloudName
112 typename parcelType::constantProperties
constProps_;
245 inline const typename parcelType::constantProperties&
249 inline typename parcelType::constantProperties&
constProps();
345 inline scalar
Tmax()
const;
348 inline scalar
Tmin()
const;
360 const label injectori
scalar hs(const scalar p, const scalar T) const
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Templated base class for dsmc cloud.
const word & cloudName() const
Return the cloud type.
const fvMesh & mesh() const
Return references to the mesh.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
Templated heat transfer model class.
autoPtr< IOobject > clone() const
Clone.
const word & name() const
Return name.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Templated base class for thermodynamic cloud.
const parcelType::constantProperties & constProps() const
Return the constant properties.
void operator=(const ThermoCloud &)=delete
Disallow default bitwise assignment.
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].
void setModels()
Set cloud sub-models.
autoPtr< CompositionModel< ThermoCloud< CloudType > > > compositionModel_
Reacting composition model.
volScalarField::Internal & radAreaPT4()
Radiation sum of parcel projected area*temperature^4 [m2K4].
void storeState()
Store the current cloud state.
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.
const volScalarField & p_
Pressure [Pa].
volScalarField::Internal & radT4()
Radiation sum of parcel temperature^4 [K4].
void relaxSources(const ThermoCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
autoPtr< volScalarField::Internal > hsCoeff_
Coefficient for carrier phase hs equation [W/K].
autoPtr< volScalarField::Internal > radT4_
Radiation sum of parcel temperature^4.
void scaleSources()
Apply scaling to (transient) cloud sources.
tmp< volScalarField::Internal > hsTrans() const
Return sensible enthalpy transfer [J/kg].
autoPtr< integrationScheme > TIntegrator_
Temperature integration.
ThermoCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const fluidThermo &carrierThermo, const bool readFields=true)
Construct given carrier fields and thermo.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Switch radiation_
Include radiation.
const integrationScheme & TIntegrator() const
Return reference to velocity integration.
volScalarField::Internal & hsTransRef()
Access sensible enthalpy transfer [J/kg].
virtual void writeFields() const
Write the field data for the cloud.
CloudType cloudType
Type of cloud this cloud was instantiated for.
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.
void cloudReset(ThermoCloud< CloudType > &c)
Reset state of cloud.
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer() const
Return reference to heat transfer model.
void evolve()
Evolve the cloud.
autoPtr< HeatTransferModel< ThermoCloud< CloudType > > > heatTransferModel_
Heat transfer model.
parcelType::constantProperties constProps_
Thermo parcel constant properties.
autoPtr< volScalarField::Internal > radAreaPT4_
Radiation sum of parcel projected areas * temperature^4.
const fluidThermo & carrierThermo_
Thermophysical properties of the carrier fluid.
const parcelThermo thermo_
SLG thermodynamics package.
ThermoCloud< CloudType > thermoCloudType
Convenience typedef for this cloud type.
const CompositionModel< ThermoCloud< CloudType > > & composition() const
Return const access to reacting composition model.
const fluidThermo & carrierThermo() const
Return const access to carrier thermo package.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void checkParcelProperties(parcelType &parcel, const label injectori)
Check parcel properties.
void preEvolve()
Pre-evolve.
void resetSourceTerms()
Reset the cloud source terms.
const volScalarField & T_
Temperature [K].
const volScalarField & p() const
Return const access to the carrier pressure field.
tmp< volScalarField > Ep() const
Return tmp equivalent particulate emission.
virtual ~ThermoCloud()
Destructor.
autoPtr< volScalarField::Internal > hsTrans_
Sensible enthalpy transfer [J/kg].
autoPtr< volScalarField::Internal > radAreaP_
Radiation sum of parcel projected areas.
tmp< volScalarField::Internal > hsCoeff() const
Return coefficient for carrier phase hs equation.
void setParcelThermoProperties(parcelType &parcel)
Set parcel thermo properties.
volScalarField::Internal & hsCoeffRef()
Access coefficient for carrier phase hs equation.
tmp< volScalarField > ap() const
Return tmp equivalent particulate absorption.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Base-class for fluid thermodynamic properties.
Mesh data needed to do the Finite Volume discretisation.
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,...
A class for managing temporary objects.
A class for handling words, derived from string.
Forward declarations of fvMatrix specialisations.
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
TemplateName(FvFaceCellWave)