35template<
class CloudType>
38 compositionModel_.reset
42 this->subModelProperties(),
47 phaseChangeModel_.reset
51 this->subModelProperties(),
58template<
class CloudType>
66 if (YSupplied.
size() !=
Y.size())
69 << YName <<
" supplied, but size is not compatible with "
70 <<
"parcel composition: " <<
nl <<
" "
71 << YName <<
"(" << YSupplied.
size() <<
") vs required composition "
72 << YName <<
"(" <<
Y.size() <<
")" <<
nl
78template<
class CloudType>
81 CloudType::cloudReset(c);
83 compositionModel_.reset(c.compositionModel_.ptr());
84 phaseChangeModel_.reset(c.phaseChangeModel_.ptr());
90template<
class CloudType>
103 cloudCopyPtr_(nullptr),
104 constProps_(this->particleProperties()),
105 compositionModel_(nullptr),
106 phaseChangeModel_(nullptr),
107 rhoTrans_(
thermo.carrier().species().size())
123 const word& specieName =
thermo.carrier().species()[i];
131 this->
name() +
":rhoTrans_" + specieName,
143 if (this->
solution().resetSourcesOnStartup())
150template<
class CloudType>
159 cloudCopyPtr_(nullptr),
160 constProps_(c.constProps_),
161 compositionModel_(c.compositionModel_->clone()),
162 phaseChangeModel_(c.phaseChangeModel_->clone()),
163 rhoTrans_(c.rhoTrans_.size())
175 this->
name() +
":rhoTrans_" + specieName,
189template<
class CloudType>
199 cloudCopyPtr_(nullptr),
201 compositionModel_(c.compositionModel_->clone()),
202 phaseChangeModel_(nullptr),
209template<
class CloudType>
213 const scalar lagrangianDt
222template<
class CloudType>
226 const scalar lagrangianDt,
227 const bool fullyDescribed
234 checkSuppliedComposition
243 parcel.mass0() = parcel.mass();
247template<
class CloudType>
254 clone(this->
name() +
"Copy").ptr()
260template<
class CloudType>
263 cloudReset(cloudCopyPtr_());
264 cloudCopyPtr_.clear();
268template<
class CloudType>
274 rhoTrans_[i].field() = 0.0;
279template<
class CloudType>
291 dsfType& rhoT = rhoTrans_[fieldi];
292 const dsfType& rhoT0 = cloudOldTime.
rhoTrans()[fieldi];
293 this->
relax(rhoT, rhoT0,
"rho");
298template<
class CloudType>
307 dsfType& rhoT = rhoTrans_[fieldi];
308 this->scale(rhoT,
"rho");
313template<
class CloudType>
325template<
class CloudType>
334template<
class CloudType>
339 this->phaseChange().info(
Info);
343template<
class CloudType>
346 if (compositionModel_)
353template<
class CloudType>
const uniformDimensionedVectorField & g
void deleteLostParticles()
Remove lost particles from cloud and delete.
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Templated base class for dsmc cloud.
const fvMesh & mesh() const
Return reference to the mesh.
void info() const
Print cloud information.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
void scaleSources()
Apply scaling to (transient) cloud sources.
void resetSourceTerms()
Reset the cloud source terms.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Templated phase change model class.
Templated base class for reacting cloud.
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
void setModels()
Set cloud sub-models.
void storeState()
Store the current cloud state.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
void cloudReset(ReactingCloud< CloudType > &c)
Reset state of cloud.
void scaleSources()
Apply scaling to (transient) cloud sources.
void checkSuppliedComposition(const scalarField &YSupplied, const scalarField &Y, const word &YName)
Check that size of a composition field is valid.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
volScalarField::Internal & rhoTrans(const label i)
Mass.
virtual void writeFields() const
Write the field data for the cloud.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
void evolve()
Evolve the cloud.
void relaxSources(const ReactingCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
void resetSourceTerms()
Reset the cloud source terms.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
const basicSpecieMixture & carrier() const
Return reference to the gaseous components.
void size(const label n)
Older name for setAddressableSize.
const speciesTable & species() const
Return the table of species.
Class used to pass data into container.
void writeFields() const
Write fields.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Allows specification of different writing frequency of objects registered to the database.
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.
const Time & time() const noexcept
Return time registry.
Virtual abstract base class for templated ReactingCloud.
Selector class for relaxation factors, solver type and solution.
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.
basicSpecieMixture & composition
PtrList< volScalarField > & Y
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.
const word cloudName(propsDict.get< word >("cloud"))