Go to the documentation of this file.
35 template<
class CloudType>
38 devolatilisationModel_.reset
42 this->subModelProperties(),
47 surfaceReactionModel_.reset
51 this->subModelProperties(),
58 template<
class CloudType>
64 CloudType::cloudReset(
c);
66 devolatilisationModel_.reset(
c.devolatilisationModel_.ptr());
67 surfaceReactionModel_.reset(
c.surfaceReactionModel_.ptr());
69 dMassDevolatilisation_ =
c.dMassDevolatilisation_;
70 dMassSurfaceReaction_ =
c.dMassSurfaceReaction_;
76 template<
class CloudType>
89 cloudCopyPtr_(
nullptr),
90 constProps_(this->particleProperties()),
91 devolatilisationModel_(
nullptr),
92 surfaceReactionModel_(
nullptr),
93 dMassDevolatilisation_(0.0),
94 dMassSurfaceReaction_(0.0)
103 this->deleteLostParticles();
107 if (this->
solution().resetSourcesOnStartup())
114 template<
class CloudType>
123 cloudCopyPtr_(
nullptr),
124 constProps_(
c.constProps_),
125 devolatilisationModel_(
c.devolatilisationModel_->clone()),
126 surfaceReactionModel_(
c.surfaceReactionModel_->clone()),
127 dMassDevolatilisation_(
c.dMassDevolatilisation_),
128 dMassSurfaceReaction_(
c.dMassSurfaceReaction_)
132 template<
class CloudType>
142 cloudCopyPtr_(
nullptr),
144 devolatilisationModel_(
nullptr),
145 surfaceReactionModel_(
nullptr),
146 dMassDevolatilisation_(0.0),
147 dMassSurfaceReaction_(0.0)
153 template<
class CloudType>
160 template<
class CloudType>
164 const scalar lagrangianDt
167 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
174 parcel.YLiquid() = this->
composition().Y0(idLiquid);
179 template<
class CloudType>
183 const scalar lagrangianDt,
184 const bool fullyDescribed
187 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
195 this->checkSuppliedComposition
201 this->checkSuppliedComposition
207 this->checkSuppliedComposition
217 template<
class CloudType>
224 clone(this->
name() +
"Copy").ptr()
230 template<
class CloudType>
233 cloudReset(cloudCopyPtr_());
234 cloudCopyPtr_.clear();
238 template<
class CloudType>
241 CloudType::resetSourceTerms();
245 template<
class CloudType>
250 typename parcelType::trackingData td(*
this);
252 this->
solve(*
this, td);
257 template<
class CloudType>
269 template<
class CloudType>
274 this->devolatilisation().info(
Info);
275 this->surfaceReaction().info(
Info);
279 template<
class CloudType>
282 if (this->compositionModel_.valid())
Templated surface reaction model class.
void cloudReset(ReactingMultiphaseCloud< CloudType > &c)
Reset state of cloud.
Selector class for relaxation factors, solver type and solution.
void restoreState()
Reset the current cloud to the previously stored state.
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,...
void storeState()
Store the current cloud state.
Virtual abstract base class for templated reactingMultiphaseCloud.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
void info()
Print cloud information.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
ParticleType parcelType
Parcels are just particles.
virtual void writeFields() const
Write the field data for the cloud.
basicSpecieMixture & composition
void setModels()
Set cloud sub-models.
messageStream Info
Information stream (uses stdout - output is on the master only)
word name(const complex &c)
Return string representation of complex.
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields)
Templated devolatilisation model class.
DSMCCloud< dsmcParcel > CloudType
Mesh data needed to do the Finite Volume discretisation.
const uniformDimensionedVectorField & g
virtual ~ReactingMultiphaseCloud()
Destructor.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void evolve()
Evolve the cloud.
void resetSourceTerms()
Reset the cloud source terms.
Base cloud calls templated on particle type.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const dimensionedScalar c
Speed of light in a vacuum.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Templated base class for multiphase reacting cloud.