41 #ifndef ReactingMultiphaseCloud_H 42 #define ReactingMultiphaseCloud_H 53 template<
class CloudType>
56 template<
class CloudType>
63 template<
class CloudType>
105 typename parcelType::constantProperties
constProps_;
211 inline const typename parcelType::constantProperties&
215 inline typename parcelType::constantProperties&
constProps();
255 const scalar lagrangianDt
262 const scalar lagrangianDt,
263 const bool fullyDescribed
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
autoPtr< IOobject > clone() const
Clone.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
const word & name() const
Return name.
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.
void setModels()
Set cloud sub-models.
const word & cloudName() const
Return the cloud type.
parcelType::constantProperties constProps_
Parcel constant properties.
rhoReactionThermo & thermo
autoPtr< SurfaceReactionModel< ReactingMultiphaseCloud< CloudType > > > surfaceReactionModel_
Surface reaction model.
Templated base class for multiphase reacting cloud.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void cloudReset(ReactingMultiphaseCloud< CloudType > &c)
Reset state of cloud.
void restoreState()
Reset the current cloud to the previously stored state.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
A class for handling words, derived from string.
CloudType cloudType
Type of cloud this cloud was instantiated for.
void resetSourceTerms()
Reset the cloud source terms.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
void storeState()
Store the current cloud state.
Virtual abstract base class for templated reactingMultiphaseCloud.
virtual ~ReactingMultiphaseCloud()
Destructor.
ReactingMultiphaseCloud< CloudType > reactingMultiphaseCloudType
Convenience typedef for this cloud type.
autoPtr< DevolatilisationModel< ReactingMultiphaseCloud< CloudType > > > devolatilisationModel_
Devolatilisation model.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const fvMesh & mesh() const
Return references to the mesh.
virtual void writeFields() const
Write the field data for the cloud.
const ReactingMultiphaseCloud & cloudCopy() const
Return a reference to the cloud copy.
Mesh data needed to do the Finite Volume discretisation.
const dimensionedScalar c
Speed of light in a vacuum.
scalar dMassDevolatilisation_
Total mass transferred to continuous phase via devolatilisation.
const SurfaceReactionModel< ReactingMultiphaseCloud< CloudType > > & surfaceReaction() const
Return const access to reacting surface reaction model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Templated devolatilisation model class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
const dimensionedVector & g
void info()
Print cloud information.
void evolve()
Evolve the cloud.
scalar dMassSurfaceReaction_
Total mass transferred to continuous phase via surface.
Templated surface reaction model class.
Templated base class for dsmc cloud.
const DevolatilisationModel< ReactingMultiphaseCloud< CloudType > > & devolatilisation() const
Return const access to devolatilisation model.