33 template<
class CloudType>
36 heterogeneousReactionModel_.reset
40 this->subModelProperties(),
47 template<
class CloudType>
53 CloudType::cloudReset(
c);
54 heterogeneousReactionModel_.reset(
c.heterogeneousReactionModel_.ptr());
60 template<
class CloudType>
73 cloudCopyPtr_(nullptr),
74 heterogeneousReactionModel_(nullptr)
83 this->deleteLostParticles();
89 template<
class CloudType>
92 ReactingHeterogeneousCloud<CloudType>&
c,
97 reactingHeterogeneousCloud(),
98 cloudCopyPtr_(nullptr),
99 heterogeneousReactionModel_(
c.heterogeneousReactionModel_->clone())
103 template<
class CloudType>
113 cloudCopyPtr_(nullptr),
114 heterogeneousReactionModel_(
c.heterogeneousReactionModel_->clone())
120 template<
class CloudType>
124 const scalar lagrangianDt
127 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
135 parcel.F().setSize(heterogeneousReactionModel_->nF(), 0.0);
138 parcel.canCombust() = 1;
142 if (this->constProps_.rho0() == -1)
145 const label idLiquid = this->
composition().idLiquid();
146 const label idSolid = this->
composition().idSolid();
153 this->
composition().thermo().thermo().p()[parcel.cell()];
154 const scalar
T0 = this->constProps_.T0();
161 template<
class CloudType>
165 const scalar lagrangianDt,
166 const bool fullyDescribed
169 CloudType::checkParcelProperties(parcel, lagrangianDt,
false);
172 const label liqId = this->
composition().idLiquid();
184 <<
"The supplied composition must be : " <<
nl 185 <<
" YGasTot0 0 : " <<
nl 186 <<
" YLiquidTot0 0 : " <<
nl 187 <<
" YSolidTot0 1 : " <<
nl 188 <<
"This Cloud only works with pure solid particles." 194 <<
"The supplied composition has a liquid phase. " <<
nl 195 <<
"This Cloud only works with pure solid particles." 201 template<
class CloudType>
206 static_cast<ReactingHeterogeneousCloud<CloudType>*
> 208 clone(this->
name() +
"Copy").ptr()
214 template<
class CloudType>
217 cloudReset(cloudCopyPtr_());
218 cloudCopyPtr_.clear();
222 template<
class CloudType>
227 typename parcelType::trackingData td(*
this);
234 template<
class CloudType>
237 const mapPolyMesh& mapper
246 template<
class CloudType>
251 heterogeneousReactionModel_->info();
255 template<
class CloudType>
262 template<
class CloudType>
266 CloudType::particleType::readObjects(*
this, this->
composition(), obr);
270 template<
class CloudType>
274 CloudType::particleType::writeObjects(*
this, this->
composition(), obr);
DSMCCloud< dsmcParcel > CloudType
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
basicSpecieMixture & composition
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
constexpr char nl
The newline '\n' character (0x0a)
void restoreState()
Reset the current cloud to the previously stored state.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Templated base class for reacting heterogeneous cloud.
Virtual abstract base class for templated ReactingCloud.
void storeState()
Store the current cloud state.
void cloudReset(ReactingHeterogeneousCloud< CloudType > &c)
Reset state of cloud.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
void setModels()
Set cloud sub-models.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
errorManip< error > abort(error &err)
Base cloud calls templated on particle type.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
const uniformDimensionedVectorField & g
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
void evolve()
Evolve the cloud.
void info()
Print cloud information.
Base class for heterogeneous reacting models.
Mesh data needed to do the Finite Volume discretisation.
const dimensionedScalar c
Speed of light in a vacuum.
virtual void writeFields() const
Write the field data for the cloud.
Selector class for relaxation factors, solver type and solution.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type.
Registry of regIOobjects.
Templated base class for dsmc cloud.
const volScalarField & p0
virtual void readObjects(const objectRegistry &obr)
Read particle fields as objects from the obr registry.