Go to the documentation of this file.
33 template<
class CloudType>
36 compositionModel_.reset
40 this->subModelProperties(),
45 phaseChangeModel_.reset
49 this->subModelProperties(),
56 template<
class CloudType>
64 if (YSupplied.size() !=
Y.size())
67 << YName <<
" supplied, but size is not compatible with "
68 <<
"parcel composition: " <<
nl <<
" "
69 << YName <<
"(" << YSupplied.size() <<
") vs required composition "
70 << YName <<
"(" <<
Y.size() <<
")" <<
nl
76 template<
class CloudType>
79 CloudType::cloudReset(
c);
81 compositionModel_.reset(
c.compositionModel_.ptr());
82 phaseChangeModel_.reset(
c.phaseChangeModel_.ptr());
88 template<
class CloudType>
102 constProps_(this->particleProperties()),
103 compositionModel_(NULL),
104 phaseChangeModel_(NULL),
105 rhoTrans_(
thermo.carrier().species().size())
120 const word& specieName =
thermo.carrier().species()[i];
128 this->
name() +
":rhoTrans_" + specieName,
131 IOobject::READ_IF_PRESENT,
140 if (this->
solution().resetSourcesOnStartup())
147 template<
class CloudType>
157 constProps_(
c.constProps_),
158 compositionModel_(
c.compositionModel_->clone()),
159 phaseChangeModel_(
c.phaseChangeModel_->clone()),
160 rhoTrans_(
c.rhoTrans_.size())
164 const word& specieName = this->
thermo().carrier().species()[i];
172 this->
name() +
":rhoTrans_" + specieName,
186 template<
class CloudType>
198 compositionModel_(
c.compositionModel_->clone()),
200 phaseChangeModel_(NULL),
207 template<
class CloudType>
214 template<
class CloudType>
218 const scalar lagrangianDt
221 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
223 parcel.pc() = this->
thermo().thermo().p()[parcel.cell()];
228 template<
class CloudType>
232 const scalar lagrangianDt,
233 const bool fullyDescribed
236 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
240 checkSuppliedComposition
249 parcel.mass0() = parcel.mass();
253 template<
class CloudType>
260 clone(this->
name() +
"Copy").ptr()
266 template<
class CloudType>
269 cloudReset(cloudCopyPtr_());
270 cloudCopyPtr_.clear();
274 template<
class CloudType>
277 CloudType::resetSourceTerms();
280 rhoTrans_[i].field() = 0.0;
285 template<
class CloudType>
291 CloudType::relaxSources(cloudOldTime);
297 dsfType& rhoT = rhoTrans_[fieldI];
298 const dsfType& rhoT0 = cloudOldTime.
rhoTrans()[fieldI];
299 this->
relax(rhoT, rhoT0,
"rho");
304 template<
class CloudType>
307 CloudType::scaleSources();
313 dsfType& rhoT = rhoTrans_[fieldI];
314 this->scale(rhoT,
"rho");
319 template<
class CloudType>
324 typename parcelType::template
325 TrackingData<ReactingCloud<CloudType> > td(*
this);
332 template<
class CloudType>
345 template<
class CloudType>
350 this->phaseChange().info(
Info);
354 template<
class CloudType>
359 CloudType::particleType::writeFields(*
this, this->
composition());
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Templated phase change model class.
void resetSourceTerms()
Reset the cloud source terms.
virtual void writeFields() const
Write the field data for the cloud.
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
#define forAll(list, i)
Loop across all elements in list.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
void evolve()
Evolve the cloud.
virtual ~ReactingCloud()
Destructor.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
const dimensionedVector & g
This function object reads fields from the time directories and adds them to the mesh database for fu...
ReactingCloud(const ReactingCloud &)
Disallow default bitwise copy construct.
void scaleSources()
Apply scaling to (transient) cloud sources.
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Pre-declare SubField and related Field type.
DimensionedField< scalar, volMesh > & rhoTrans(const label i)
Mass.
void setModels()
Set cloud sub-models.
Virtual abstract base class for templated ReactingCloud.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
DSMCCloud< dsmcParcel > CloudType
void storeState()
Store the current cloud state.
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
errorManip< error > abort(error &err)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Base cloud calls templated on particle type.
basicMultiComponentMixture & composition
const word cloudName(propsDict.lookup("cloudName"))
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const dimensionedScalar c
Speed of light in a vacuum.
void restoreState()
Reset the current cloud to the previously stored state.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
void info()
Print cloud information.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
Generic GeometricField class.
word name(const complex &)
Return a string representation of a complex.
PtrList< volScalarField > & Y
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Templated base class for reacting cloud.
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.