46#ifndef ReactingCloud_H
47#define ReactingCloud_H
58template<
class CloudType>
59class CompositionModel;
61template<
class CloudType>
62class PhaseChangeModel;
68template<
class CloudType>
259 Srho(
const label i)
const;
275 const scalar lagrangianDt
282 const scalar lagrangianDt,
283 const bool fullyDescribed
const uniformDimensionedVectorField & g
ParticleType particleType
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 reference to the mesh.
Class to hold DSMC particle constant properties.
const word & name() const noexcept
Return the object name.
autoPtr< IOobject > clone() const
Clone.
Templated phase change model class.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Templated base class for reacting cloud.
virtual ~ReactingCloud()=default
Destructor.
autoPtr< PhaseChangeModel< ReactingCloud< CloudType > > > phaseChangeModel_
Reacting phase change model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
void setModels()
Set cloud sub-models.
void storeState()
Store the current cloud state.
ReactingCloud< CloudType > reactingCloudType
Convenience typedef for this cloud type.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
void scaleSources()
Apply scaling to (transient) cloud sources.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
virtual void writeFields() const
Write the field data for the cloud.
CloudType cloudType
Type of cloud this cloud was instantiated for.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
void evolve()
Evolve the cloud.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
parcelType::constantProperties constProps_
Parcel constant properties.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
autoPtr< CompositionModel< ReactingCloud< CloudType > > > compositionModel_
Reacting composition model.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void resetSourceTerms()
Reset the cloud source terms.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Mesh data needed to do the Finite Volume discretisation.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
Virtual abstract base class for templated ReactingCloud.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
PtrList< volScalarField > & Y
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.