ThermoCloud< CloudType > Class Template Reference

Templated base class for thermodynamic cloud. More...

Inheritance diagram for ThermoCloud< CloudType >:
[legend]
Collaboration diagram for ThermoCloud< CloudType >:
[legend]

Public Types

typedef CloudType cloudType
 Type of cloud this cloud was instantiated for. More...
 
typedef CloudType::particleType parcelType
 Type of parcel the cloud was instantiated for. More...
 
typedef ThermoCloud< CloudTypethermoCloudType
 Convenience typedef for this cloud type. More...
 
- Public Types inherited from DSMCCloud< ParcelType >
typedef ParcelType parcelType
 Type of parcel the cloud was instantiated for. More...
 
- Public Types inherited from Cloud< ParticleType >
typedef ParticleType particleType
 
typedef ParticleType parcelType
 Parcels are just particles. More...
 

Public Member Functions

 ThermoCloud (const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const SLGThermo &thermo, bool readFields=true)
 Construct given carrier gas fields. More...
 
 ThermoCloud (ThermoCloud< CloudType > &c, const word &name)
 Copy constructor with new name. More...
 
 ThermoCloud (const fvMesh &mesh, const word &name, const ThermoCloud< CloudType > &c)
 Copy constructor with new name - creates bare cloud. More...
 
virtual autoPtr< Cloud< parcelType > > clone (const word &name)
 Construct and return clone based on (this) with new name. More...
 
virtual autoPtr< Cloud< parcelType > > cloneBare (const word &name) const
 Construct and return bare clone based on (this) with new name. More...
 
virtual ~ThermoCloud ()=default
 Destructor. More...
 
const ThermoCloudcloudCopy () const
 Return a reference to the cloud copy. More...
 
const parcelType::constantProperties & constProps () const
 Return the constant properties. More...
 
parcelType::constantProperties & constProps ()
 Return access to the constant properties. More...
 
const SLGThermothermo () const
 Return const access to thermo package. More...
 
const volScalarFieldT () const
 Return const access to the carrier temperature field. More...
 
const volScalarFieldp () const
 Return const access to the carrier pressure field. More...
 
const HeatTransferModel< ThermoCloud< CloudType > > & heatTransfer () const
 Return reference to heat transfer model. More...
 
const integrationSchemeTIntegrator () const
 Return reference to velocity integration. More...
 
bool radiation () const
 Radiation flag. More...
 
volScalarField::InternalradAreaP ()
 Radiation sum of parcel projected areas [m2]. More...
 
const volScalarField::InternalradAreaP () const
 Radiation sum of parcel projected areas [m2]. More...
 
volScalarField::InternalradT4 ()
 Radiation sum of parcel temperature^4 [K4]. More...
 
const volScalarField::InternalradT4 () const
 Radiation sum of parcel temperature^4 [K4]. More...
 
volScalarField::InternalradAreaPT4 ()
 Radiation sum of parcel projected area*temperature^4 [m2K4]. More...
 
const volScalarField::InternalradAreaPT4 () const
 Radiation sum of parcel temperature^4 [m2K4]. More...
 
volScalarField::InternalhsTrans ()
 Sensible enthalpy transfer [J/kg]. More...
 
const volScalarField::InternalhsTrans () const
 Sensible enthalpy transfer [J/kg]. More...
 
volScalarField::InternalhsCoeff ()
 Return coefficient for carrier phase hs equation. More...
 
const volScalarField::InternalhsCoeff () const
 Return const coefficient for carrier phase hs equation. More...
 
tmp< fvScalarMatrixSh (volScalarField &hs) const
 Return sensible enthalpy source term [J/kg/m3/s]. More...
 
tmp< volScalarFieldEp () const
 Return tmp equivalent particulate emission. More...
 
tmp< volScalarFieldap () const
 Return tmp equivalent particulate absorption. More...
 
tmp< volScalarFieldsigmap () const
 Return tmp equivalent particulate scattering factor. More...
 
scalar Tmax () const
 Maximum temperature. More...
 
scalar Tmin () const
 Minimum temperature. More...
 
void setParcelThermoProperties (parcelType &parcel, const scalar lagrangianDt)
 Set parcel thermo properties. More...
 
void checkParcelProperties (parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
 Check parcel properties. More...
 
void storeState ()
 Store the current cloud state. More...
 
void restoreState ()
 Reset the current cloud to the previously stored state. More...
 
void resetSourceTerms ()
 Reset the cloud source terms. More...
 
void relaxSources (const ThermoCloud< CloudType > &cloudOldTime)
 Apply relaxation to (steady state) cloud sources. More...
 
void scaleSources ()
 Apply scaling to (transient) cloud sources. More...
 
void preEvolve (const typename parcelType::trackingData &td)
 Pre-evolve. More...
 
void evolve ()
 Evolve the cloud. More...
 
virtual void autoMap (const mapPolyMesh &)
 Remap the cells of particles corresponding to the. More...
 
void info ()
 Print cloud information. More...
 
- Public Member Functions inherited from DSMCCloud< ParcelType >
 DSMCCloud (const word &cloudName, const fvMesh &mesh, bool readFields=true)
 Construct given name and mesh, will read Parcels and fields from. More...
 
 DSMCCloud (const word &cloudName, const fvMesh &mesh, const IOdictionary &dsmcInitialiseDict)
 Construct given name, mesh and initialisation dictionary. More...
 
virtual ~DSMCCloud ()
 Destructor. More...
 
const wordcloudName () const
 Return the cloud type. More...
 
const fvMeshmesh () const
 Return reference to the mesh. More...
 
const IOdictionaryparticleProperties () const
 Return particle properties dictionary. More...
 
const List< word > & typeIdList () const
 Return the idList. More...
 
scalar nParticle () const
 Return the number of real particles represented by one. More...
 
const List< DynamicList< ParcelType * > > & cellOccupancy () const
 Return the cell occupancy addressing. More...
 
volScalarFieldsigmaTcRMax ()
 Return the sigmaTcRMax field. non-const access to allow. More...
 
scalarFieldcollisionSelectionRemainder ()
 Return the collision selection remainder field. non-const. More...
 
const List< typename ParcelType::constantProperties > & constProps () const
 Return all of the constant properties. More...
 
const ParcelType::constantProperties & constProps (label typeId) const
 Return the constant properties of the given typeId. More...
 
RandomrndGen ()
 Return reference to the random object. More...
 
volScalarField::Boundary & qBF ()
 Return non-const heat flux boundary field reference. More...
 
volVectorField::Boundary & fDBF ()
 Return non-const force density at boundary field reference. More...
 
volScalarField::Boundary & rhoNBF ()
 Return non-const number density boundary field reference. More...
 
volScalarField::Boundary & rhoMBF ()
 Return non-const mass density boundary field reference. More...
 
volScalarField::Boundary & linearKEBF ()
 Return non-const linear kinetic energy density boundary. More...
 
volScalarField::Boundary & internalEBF ()
 Return non-const internal energy density boundary field. More...
 
volScalarField::Boundary & iDofBF ()
 Return non-const internal degree of freedom density boundary. More...
 
volVectorField::Boundary & momentumBF ()
 Return non-const momentum density boundary field reference. More...
 
const volScalarFieldboundaryT () const
 Return macroscopic temperature. More...
 
const volVectorFieldboundaryU () const
 Return macroscopic velocity. More...
 
const volScalarFieldq () const
 Return heat flux at surface field. More...
 
const volVectorFieldfD () const
 Return force density at surface field. More...
 
const volScalarFieldrhoN () const
 Return the real particle number density field. More...
 
const volScalarFieldrhoM () const
 Return the particle mass density field. More...
 
const volScalarFielddsmcRhoN () const
 Return the field of number of DSMC particles. More...
 
const volScalarFieldlinearKE () const
 Return the total linear kinetic energy (translational and. More...
 
const volScalarFieldinternalE () const
 Return the internal energy density field. More...
 
const volScalarFieldiDof () const
 Return the average internal degrees of freedom field. More...
 
const volVectorFieldmomentum () const
 Return the momentum density field. More...
 
vector equipartitionLinearVelocity (scalar temperature, scalar mass)
 Generate a random velocity sampled from the Maxwellian speed. More...
 
scalar equipartitionInternalEnergy (scalar temperature, direction internalDegreesOfFreedom)
 Generate a random internal energy, sampled from the. More...
 
scalar maxwellianAverageSpeed (scalar temperature, scalar mass) const
 Average particle speed. More...
 
scalarField maxwellianAverageSpeed (scalarField temperature, scalar mass) const
 
scalar maxwellianRMSSpeed (scalar temperature, scalar mass) const
 RMS particle speed. More...
 
scalarField maxwellianRMSSpeed (scalarField temperature, scalar mass) const
 
scalar maxwellianMostProbableSpeed (scalar temperature, scalar mass) const
 Most probable speed. More...
 
scalarField maxwellianMostProbableSpeed (scalarField temperature, scalar mass) const
 
const BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision () const
 Return reference to binary elastic collision model. More...
 
BinaryCollisionModel< DSMCCloud< ParcelType > > & binaryCollision ()
 Return non-const reference to binary elastic collision model. More...
 
const WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction () const
 Return reference to wall interaction model. More...
 
WallInteractionModel< DSMCCloud< ParcelType > > & wallInteraction ()
 Return non-const reference to wall interaction model. More...
 
const InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary () const
 Return reference to wall interaction model. More...
 
InflowBoundaryModel< DSMCCloud< ParcelType > > & inflowBoundary ()
 Return non-const reference to wall interaction model. More...
 
scalar massInSystem () const
 Total mass in system. More...
 
vector linearMomentumOfSystem () const
 Total linear momentum of the system. More...
 
scalar linearKineticEnergyOfSystem () const
 Total linear kinetic energy in the system. More...
 
scalar internalEnergyOfSystem () const
 Total internal energy in the system. More...
 
void info () const
 Print cloud information. More...
 
void dumpParticlePositions () const
 Dump particle positions to .obj file. More...
 
void addNewParcel (const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
 Add new parcel. More...
 
void evolve ()
 Evolve the cloud (move, collide) More...
 
void clear ()
 Clear the Cloud. More...
 
- Public Member Functions inherited from Cloud< ParticleType >
 TypeName ("Cloud")
 Runtime type information. More...
 
 Cloud (const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
 Construct from mesh and a list of particles. More...
 
 Cloud (const polyMesh &pMesh, const word &cloudName, const bool checkClass=true)
 Construct from mesh by reading from file with given cloud instance. More...
 
const polyMeshpMesh () const
 Return the polyMesh reference. More...
 
virtual label nParcels () const
 Return the number of particles in the cloud. More...
 
DynamicList< label > & labels () const
 Return temporary addressing. More...
 
void addParticle (ParticleType *pPtr)
 Transfer particle to cloud. More...
 
void deleteParticle (ParticleType &p)
 Remove particle from cloud and delete. More...
 
void deleteLostParticles ()
 Remove lost particles from cloud and delete. More...
 
void cloudReset (const Cloud< ParticleType > &c)
 Reset the particles. More...
 
template<class TrackCloudType >
void move (TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
 Move the particles. More...
 
void autoMap (const mapPolyMesh &)
 Remap the cells of particles corresponding to the. More...
 
IOobject fieldIOobject (const word &fieldName, const IOobject::readOption r) const
 Helper to construct IOobject for field and current time. More...
 
template<class DataType >
void checkFieldIOobject (const Cloud< ParticleType > &c, const IOField< DataType > &data) const
 Check lagrangian data field. More...
 
template<class DataType >
void checkFieldFieldIOobject (const Cloud< ParticleType > &c, const CompactIOField< Field< DataType >, DataType > &data) const
 Check lagrangian data fieldfield. More...
 
template<class Type >
bool readStoreFile (const IOobject &io, const IOobject &ioNew) const
 Helper function to store a cloud field on its registry. More...
 
void readFromFiles (objectRegistry &obr, const wordRes &selectFields) const
 Read from files into objectRegistry. More...
 
virtual void writeFields () const
 Write the field data for the cloud of particles Dummy at. More...
 
virtual bool writeObject (IOstreamOption streamOpt, const bool valid) const
 Write using stream options. More...
 
void writePositions () const
 Write positions to <cloudName>_positions.obj file. More...
 
void storeGlobalPositions () const
 Call this before a topology change. More...
 
- Public Member Functions inherited from DSMCBaseCloud
 TypeName ("DSMCBaseCloud")
 Runtime type information. More...
 
 DSMCBaseCloud ()=default
 Null constructor. More...
 
virtual ~DSMCBaseCloud ()=default
 Destructor. More...
 
- Public Member Functions inherited from thermoCloud
 TypeName ("thermoCloud")
 Runtime type information. More...
 
 thermoCloud ()=default
 Null constructor. More...
 
virtual ~thermoCloud ()=default
 Destructor. More...
 

Protected Member Functions

void setModels ()
 Set cloud sub-models. More...
 
void cloudReset (ThermoCloud< CloudType > &c)
 Reset state of cloud. More...
 

Protected Attributes

parcelType::constantProperties constProps_
 Thermo parcel constant properties. More...
 
const SLGThermothermo_
 SLG thermodynamics package. More...
 
const volScalarFieldT_
 Temperature [K]. More...
 
const volScalarFieldp_
 Pressure [Pa]. More...
 
autoPtr< HeatTransferModel< ThermoCloud< CloudType > > > heatTransferModel_
 Heat transfer model. More...
 
autoPtr< integrationSchemeTIntegrator_
 Temperature integration. More...
 
Switch radiation_
 Include radiation. More...
 
autoPtr< volScalarField::InternalradAreaP_
 Radiation sum of parcel projected areas. More...
 
autoPtr< volScalarField::InternalradT4_
 Radiation sum of parcel temperature^4. More...
 
autoPtr< volScalarField::InternalradAreaPT4_
 Radiation sum of parcel projected areas * temperature^4. More...
 
autoPtr< volScalarField::InternalhsTrans_
 Sensible enthalpy transfer [J/kg]. More...
 
autoPtr< volScalarField::InternalhsCoeff_
 Coefficient for carrier phase hs equation [W/K]. More...
 
- Protected Attributes inherited from Cloud< ParticleType >
cloud::geometryType geometryType_
 Geometry type. More...
 

Additional Inherited Members

- Static Public Attributes inherited from Cloud< ParticleType >
static word cloudPropertiesName
 Name of cloud properties dictionary. More...
 

Detailed Description

template<class CloudType>
class Foam::ThermoCloud< CloudType >

Templated base class for thermodynamic cloud.

- Adds to kinematic cloud

  • Heat transfer
Source files

Definition at line 66 of file ThermoCloud.H.

Member Typedef Documentation

◆ cloudType

Type of cloud this cloud was instantiated for.

Definition at line 76 of file ThermoCloud.H.

◆ parcelType

Type of parcel the cloud was instantiated for.

Definition at line 79 of file ThermoCloud.H.

◆ thermoCloudType

Convenience typedef for this cloud type.

Definition at line 82 of file ThermoCloud.H.

Constructor & Destructor Documentation

◆ ThermoCloud() [1/3]

ThermoCloud ( const word cloudName,
const volScalarField rho,
const volVectorField U,
const dimensionedVector g,
const SLGThermo thermo,
bool  readFields = true 
)

Construct given carrier gas fields.

Definition at line 131 of file ThermoCloud.C.

◆ ThermoCloud() [2/3]

ThermoCloud ( ThermoCloud< CloudType > &  c,
const word name 
)

Copy constructor with new name.

Definition at line 214 of file ThermoCloud.C.

◆ ThermoCloud() [3/3]

ThermoCloud ( const fvMesh mesh,
const word name,
const ThermoCloud< CloudType > &  c 
)

Copy constructor with new name - creates bare cloud.

Definition at line 323 of file ThermoCloud.C.

◆ ~ThermoCloud()

virtual ~ThermoCloud ( )
virtualdefault

Destructor.

Member Function Documentation

◆ setModels()

void setModels ( )
protected

Set cloud sub-models.

Definition at line 37 of file ThermoCloud.C.

◆ cloudReset()

void cloudReset ( ThermoCloud< CloudType > &  c)
protected

Reset state of cloud.

Definition at line 116 of file ThermoCloud.C.

◆ clone()

virtual autoPtr<Cloud<parcelType> > clone ( const word name)
inlinevirtual

Construct and return clone based on (this) with new name.

Definition at line 200 of file ThermoCloud.H.

◆ cloneBare()

virtual autoPtr<Cloud<parcelType> > cloneBare ( const word name) const
inlinevirtual

Construct and return bare clone based on (this) with new name.

Definition at line 209 of file ThermoCloud.H.

◆ cloudCopy()

const Foam::ThermoCloud< CloudType > & cloudCopy ( ) const
inline

Return a reference to the cloud copy.

Definition at line 37 of file ThermoCloudI.H.

◆ constProps() [1/2]

const CloudType::particleType::constantProperties & constProps ( ) const
inline

Return the constant properties.

Definition at line 45 of file ThermoCloudI.H.

◆ constProps() [2/2]

CloudType::particleType::constantProperties & constProps ( )
inline

Return access to the constant properties.

Definition at line 53 of file ThermoCloudI.H.

◆ thermo()

const Foam::SLGThermo & thermo ( ) const
inline

Return const access to thermo package.

Definition at line 60 of file ThermoCloudI.H.

◆ T()

const Foam::volScalarField & T ( ) const
inline

Return const access to the carrier temperature field.

Definition at line 67 of file ThermoCloudI.H.

◆ p()

const Foam::volScalarField & p ( ) const
inline

Return const access to the carrier pressure field.

Definition at line 74 of file ThermoCloudI.H.

◆ heatTransfer()

const Foam::HeatTransferModel< Foam::ThermoCloud< CloudType > > & heatTransfer ( ) const
inline

Return reference to heat transfer model.

Definition at line 82 of file ThermoCloudI.H.

◆ TIntegrator()

const Foam::integrationScheme & TIntegrator ( ) const
inline

Return reference to velocity integration.

Definition at line 90 of file ThermoCloudI.H.

◆ radiation()

bool radiation ( ) const
inline

Radiation flag.

Definition at line 97 of file ThermoCloudI.H.

◆ radAreaP() [1/2]

Foam::DimensionedField< Foam::scalar, Foam::volMesh > & radAreaP ( )
inline

Radiation sum of parcel projected areas [m2].

Definition at line 105 of file ThermoCloudI.H.

Referenced by ThermoCloud< Foam::DSMCCloud >::relaxSources().

Here is the caller graph for this function:

◆ radAreaP() [2/2]

const Foam::DimensionedField< Foam::scalar, Foam::volMesh > & radAreaP ( ) const
inline

Radiation sum of parcel projected areas [m2].

Definition at line 120 of file ThermoCloudI.H.

◆ radT4() [1/2]

Foam::DimensionedField< Foam::scalar, Foam::volMesh > & radT4 ( )
inline

Radiation sum of parcel temperature^4 [K4].

Definition at line 135 of file ThermoCloudI.H.

Referenced by ThermoCloud< Foam::DSMCCloud >::relaxSources().

Here is the caller graph for this function:

◆ radT4() [2/2]

const Foam::DimensionedField< Foam::scalar, Foam::volMesh > & radT4 ( ) const
inline

Radiation sum of parcel temperature^4 [K4].

Definition at line 150 of file ThermoCloudI.H.

◆ radAreaPT4() [1/2]

Foam::DimensionedField< Foam::scalar, Foam::volMesh > & radAreaPT4 ( )
inline

Radiation sum of parcel projected area*temperature^4 [m2K4].

Definition at line 165 of file ThermoCloudI.H.

Referenced by ThermoCloud< Foam::DSMCCloud >::relaxSources().

Here is the caller graph for this function:

◆ radAreaPT4() [2/2]

const Foam::DimensionedField< Foam::scalar, Foam::volMesh > & radAreaPT4 ( ) const
inline

Radiation sum of parcel temperature^4 [m2K4].

Definition at line 180 of file ThermoCloudI.H.

◆ hsTrans() [1/2]

Foam::DimensionedField< Foam::scalar, Foam::volMesh > & hsTrans ( )
inline

Sensible enthalpy transfer [J/kg].

Definition at line 195 of file ThermoCloudI.H.

Referenced by ThermoCloud< Foam::DSMCCloud >::relaxSources().

Here is the caller graph for this function:

◆ hsTrans() [2/2]

const Foam::DimensionedField< Foam::scalar, Foam::volMesh > & hsTrans ( ) const
inline

Sensible enthalpy transfer [J/kg].

Definition at line 203 of file ThermoCloudI.H.

◆ hsCoeff() [1/2]

Foam::DimensionedField< Foam::scalar, Foam::volMesh > & hsCoeff ( )
inline

Return coefficient for carrier phase hs equation.

Definition at line 211 of file ThermoCloudI.H.

Referenced by ThermoCloud< Foam::DSMCCloud >::relaxSources().

Here is the caller graph for this function:

◆ hsCoeff() [2/2]

const Foam::DimensionedField< Foam::scalar, Foam::volMesh > & hsCoeff ( ) const
inline

Return const coefficient for carrier phase hs equation.

Definition at line 219 of file ThermoCloudI.H.

◆ Sh()

Foam::tmp< Foam::fvScalarMatrix > Sh ( volScalarField hs) const
inline

Return sensible enthalpy source term [J/kg/m3/s].

Definition at line 227 of file ThermoCloudI.H.

◆ Ep()

Foam::tmp< Foam::volScalarField > Ep ( ) const
inlinevirtual

Return tmp equivalent particulate emission.

Implements thermoCloud.

Definition at line 266 of file ThermoCloudI.H.

◆ ap()

Foam::tmp< Foam::volScalarField > ap ( ) const
inlinevirtual

Return tmp equivalent particulate absorption.

Implements thermoCloud.

Definition at line 302 of file ThermoCloudI.H.

◆ sigmap()

Foam::tmp< Foam::volScalarField > sigmap ( ) const
inlinevirtual

Return tmp equivalent particulate scattering factor.

Implements thermoCloud.

Definition at line 339 of file ThermoCloudI.H.

◆ Tmax()

Foam::scalar Tmax ( ) const
inline

Maximum temperature.

Definition at line 376 of file ThermoCloudI.H.

◆ Tmin()

Foam::scalar Tmin ( ) const
inline

Minimum temperature.

Definition at line 397 of file ThermoCloudI.H.

◆ setParcelThermoProperties()

void setParcelThermoProperties ( parcelType parcel,
const scalar  lagrangianDt 
)

Set parcel thermo properties.

Definition at line 351 of file ThermoCloud.C.

◆ checkParcelProperties()

void checkParcelProperties ( parcelType parcel,
const scalar  lagrangianDt,
const bool  fullyDescribed 
)

Check parcel properties.

Definition at line 365 of file ThermoCloud.C.

◆ storeState()

void storeState ( )

Store the current cloud state.

Definition at line 376 of file ThermoCloud.C.

◆ restoreState()

void restoreState ( )

Reset the current cloud to the previously stored state.

Definition at line 389 of file ThermoCloud.C.

◆ resetSourceTerms()

void resetSourceTerms ( )

Reset the cloud source terms.

Definition at line 397 of file ThermoCloud.C.

◆ relaxSources()

void relaxSources ( const ThermoCloud< CloudType > &  cloudOldTime)

Apply relaxation to (steady state) cloud sources.

Definition at line 414 of file ThermoCloud.C.

◆ scaleSources()

void scaleSources ( )

Apply scaling to (transient) cloud sources.

Definition at line 433 of file ThermoCloud.C.

◆ preEvolve()

void preEvolve ( const typename parcelType::trackingData &  td)

Pre-evolve.

Definition at line 451 of file ThermoCloud.C.

◆ evolve()

void evolve ( )

Evolve the cloud.

Definition at line 462 of file ThermoCloud.C.

◆ autoMap()

void autoMap ( const mapPolyMesh mapper)
virtual

Remap the cells of particles corresponding to the.

mesh topology change with a default tracking data object

Reimplemented from DSMCCloud< ParcelType >.

Definition at line 474 of file ThermoCloud.C.

◆ info()

void info ( )

Print cloud information.

Definition at line 483 of file ThermoCloud.C.

Member Data Documentation

◆ constProps_

parcelType::constantProperties constProps_
protected

Thermo parcel constant properties.

Definition at line 107 of file ThermoCloud.H.

◆ thermo_

const SLGThermo& thermo_
protected

SLG thermodynamics package.

Definition at line 113 of file ThermoCloud.H.

◆ T_

const volScalarField& T_
protected

Temperature [K].

Definition at line 116 of file ThermoCloud.H.

◆ p_

const volScalarField& p_
protected

Pressure [Pa].

Definition at line 119 of file ThermoCloud.H.

◆ heatTransferModel_

autoPtr<HeatTransferModel<ThermoCloud<CloudType> > > heatTransferModel_
protected

Heat transfer model.

Definition at line 126 of file ThermoCloud.H.

◆ TIntegrator_

autoPtr<integrationScheme> TIntegrator_
protected

Temperature integration.

Definition at line 132 of file ThermoCloud.H.

◆ radiation_

Switch radiation_
protected

Include radiation.

Definition at line 138 of file ThermoCloud.H.

◆ radAreaP_

autoPtr<volScalarField::Internal> radAreaP_
protected

Radiation sum of parcel projected areas.

Definition at line 141 of file ThermoCloud.H.

◆ radT4_

autoPtr<volScalarField::Internal> radT4_
protected

Radiation sum of parcel temperature^4.

Definition at line 144 of file ThermoCloud.H.

◆ radAreaPT4_

autoPtr<volScalarField::Internal> radAreaPT4_
protected

Radiation sum of parcel projected areas * temperature^4.

Definition at line 147 of file ThermoCloud.H.

◆ hsTrans_

autoPtr<volScalarField::Internal> hsTrans_
protected

Sensible enthalpy transfer [J/kg].

Definition at line 153 of file ThermoCloud.H.

◆ hsCoeff_

autoPtr<volScalarField::Internal> hsCoeff_
protected

Coefficient for carrier phase hs equation [W/K].

Definition at line 156 of file ThermoCloud.H.


The documentation for this class was generated from the following files: