44#ifndef ReactingHeterogeneousCloud_H
45#define ReactingHeterogeneousCloud_H
55template<
class CloudType>
63template<
class CloudType>
64class ReactingHeterogeneousCloud
80 typedef ReactingHeterogeneousCloud<CloudType>
95 ReactingHeterogeneousCloud(
const ReactingHeterogeneousCloud&) =
delete;
98 void operator=(
const ReactingHeterogeneousCloud&) =
delete;
112 ReactingHeterogeneousCloud<CloudType>
128 void cloudReset(ReactingHeterogeneousCloud<CloudType>& c);
136 ReactingHeterogeneousCloud
147 ReactingHeterogeneousCloud
149 ReactingHeterogeneousCloud<CloudType>& c,
const word&
name
153 ReactingHeterogeneousCloud
157 const ReactingHeterogeneousCloud<CloudType>& c
165 new ReactingHeterogeneousCloud(*
this,
name)
174 new ReactingHeterogeneousCloud(this->
mesh(), name, *
this)
188 inline const ReactingHeterogeneousCloud&
cloudCopy()
const;
191 inline label
nF()
const;
199 ReactingHeterogeneousCloud<CloudType>
205 ReactingHeterogeneousCloud<CloudType>
215 const scalar lagrangianDt
222 const scalar lagrangianDt,
223 const bool fullyDescribed
const uniformDimensionedVectorField & g
ParticleType particleType
const word & cloudName() const
const fvMesh & mesh() const
Base class for heterogeneous reacting models.
autoPtr< IOobject > clone() const
Clone.
ReactingHeterogeneousCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const SLGThermo &thermo, bool readFields=true)
Construct given carrier gas fields.
void setModels()
Set cloud sub-models.
HeterogeneousReactingModel< ReactingHeterogeneousCloud< CloudType > > & heterogeneousReaction()
void storeState()
Store the current cloud state.
ReactingHeterogeneousCloud(const fvMesh &mesh, const word &name, const ReactingHeterogeneousCloud< CloudType > &c)
Copy constructor with new name - creates bare cloud.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
virtual void readObjects(const objectRegistry &obr)
Read particle fields as objects from the obr registry.
label nF() const
Return progress variable dimension.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
autoPtr< HeterogeneousReactingModel< ReactingHeterogeneousCloud< ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicHeterogeneousReactingParcel > > > > > > > heterogeneousReactionModel_
virtual ~ReactingHeterogeneousCloud()=default
Destructor.
const ReactingHeterogeneousCloud & cloudCopy() const
Return a reference to the cloud copy.
ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicHeterogeneousReactingParcel > > > >::particleType parcelType
virtual void writeFields() const
Write the field data for the cloud.
ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicHeterogeneousReactingParcel > > > > cloudType
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.
ReactingHeterogeneousCloud(ReactingHeterogeneousCloud< CloudType > &c, const word &name)
Copy constructor with new name.
void evolve()
Evolve the cloud.
void cloudReset(ReactingHeterogeneousCloud< CloudType > &c)
Reset state of cloud.
ReactingHeterogeneousCloud< ReactingCloud< ThermoCloud< KinematicCloud< Cloud< basicHeterogeneousReactingParcel > > > > > reactingHeterogeneousCloudType
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
const HeterogeneousReactingModel< ReactingHeterogeneousCloud< CloudType > > & heterogeneousReaction() const
Return reference to model.
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.
reactingHeterogeneousCloud()=default
Null constructor.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.