Go to the documentation of this file.
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)
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>
250 heterogeneousReactionModel_->info(
Info);
254 template<
class CloudType>
261 template<
class CloudType>
265 CloudType::particleType::readObjects(*
this, this->
composition(), obr);
269 template<
class CloudType>
273 CloudType::particleType::writeObjects(*
this, this->
composition(), obr);
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
PtrList< volScalarField > & Y()
virtual void writeFields() const
Selector class for relaxation factors, solver type and solution.
Base class for heterogeneous reacting models.
virtual void readObjects(const objectRegistry &obr)
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
const word cloudName(propsDict.get< word >("cloud"))
A class for handling words, derived from Foam::string.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Templated base class for reacting heterogeneous cloud.
Virtual abstract base class for templated ReactingCloud.
basicSpecieMixture & composition
void cloudReset(ReactingHeterogeneousCloud< CloudType > &c)
Registry of regIOobjects.
Generic templated field type.
virtual void writeObjects(objectRegistry &obr) const
void autoMap(const mapPolyMesh &)
virtual void autoMap(const mapPolyMesh &)
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
DSMCCloud< dsmcParcel > CloudType
void deleteLostParticles()
virtual void writeFields() const
Templated base class for dsmc cloud.
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const uniformDimensionedVectorField & g
errorManip< error > abort(error &err)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
#define FatalErrorInFunction
Base cloud calls templated on particle type.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const dimensionedScalar c
word name(const expressions::valueTypeCode typeCode)
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
const volScalarField & p0
Generic GeometricField class.
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)