32 template<
class CloudType>
35 phaseChangeModel_.reset
39 this->subModelProperties(),
46 template<
class CloudType>
54 if (YSupplied.
size() !=
Y.size())
57 << YName <<
" supplied, but size is not compatible with "
58 <<
"parcel composition: " <<
nl <<
" "
59 << YName <<
"(" << YSupplied.
size() <<
") vs required composition "
60 << YName <<
"(" <<
Y.size() <<
")" <<
nl
66 template<
class CloudType>
69 CloudType::cloudReset(
c);
71 phaseChangeModel_.reset(
c.phaseChangeModel_.ptr());
77 template<
class CloudType>
89 cloudCopyPtr_(nullptr),
90 constProps_(this->particleProperties()),
91 phaseChangeModel_(nullptr)
95 rhoTrans_.setSize(this->composition().carrier().species().
size());
106 const word& specieName = this->composition().carrier().species()[i];
114 this->
name() +
":rhoTrans_" + specieName,
126 if (this->
solution().resetSourcesOnStartup())
133 template<
class CloudType>
141 cloudCopyPtr_(nullptr),
142 constProps_(
c.constProps_),
143 phaseChangeModel_(
c.phaseChangeModel_->
clone()),
144 rhoTrans_(
c.rhoTrans_.size())
148 const word& specieName = this->composition().carrier().species()[i];
156 this->
name() +
":rhoTrans_" + specieName,
170 template<
class CloudType>
179 cloudCopyPtr_(nullptr),
181 phaseChangeModel_(nullptr),
188 template<
class CloudType>
195 template<
class CloudType>
201 CloudType::setParcelThermoProperties(parcel);
203 parcel.Y() = this->composition().YMixture0();
207 template<
class CloudType>
211 const label injectori
214 CloudType::checkParcelProperties(parcel, injectori);
216 if (injectori != -1 && this->injectors()[injectori].fullyDescribed())
218 checkSuppliedComposition
221 this->composition().YMixture0(),
228 template<
class CloudType>
241 template<
class CloudType>
244 cloudReset(cloudCopyPtr_());
245 cloudCopyPtr_.clear();
249 template<
class CloudType>
252 CloudType::resetSourceTerms();
255 rhoTrans_[i].primitiveFieldRef() = 0.0;
260 template<
class CloudType>
266 CloudType::relaxSources(cloudOldTime);
272 dsfType& rhoT = rhoTrans_[fieldi];
273 const dsfType& rhoT0 = cloudOldTime.
rhoTrans()[fieldi];
274 this->
relax(rhoT, rhoT0,
"rho");
279 template<
class CloudType>
282 CloudType::scaleSources();
288 dsfType& rhoT = rhoTrans_[fieldi];
289 this->scale(rhoT,
"rho");
294 template<
class CloudType>
299 typename parcelType::trackingData td(*
this);
301 this->
solve(*
this, td);
306 template<
class CloudType>
311 this->phaseChange().info(
Info);
#define forAll(list, i)
Loop across all elements in list.
void deleteLostParticles()
Remove lost particles from cloud and delete.
label size() const
Return the number of particles in the cloud.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
const fvMesh & mesh() const
Return references to the mesh.
void info() const
Print cloud information.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const objectRegistry & db() const
Return the local objectRegistry.
const word & name() const
Return name.
void size(const label)
Override size to be inconsistent with allocated storage.
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.
virtual ~ReactingCloud()
Destructor.
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.
ReactingCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const fluidThermo &carrierThermo, const bool readFields=true)
Construct given carrier fields and thermo.
volScalarField::Internal & rhoTrans(const label i)
Mass.
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 checkParcelProperties(parcelType &parcel, const label injectori)
Check parcel properties.
void resetSourceTerms()
Reset the cloud source terms.
void setParcelThermoProperties(parcelType &parcel)
Set parcel thermo properties.
Base-class for fluid thermodynamic properties.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return time.
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
word name(const bool)
Return a word representation of a bool.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
errorManip< error > abort(error &err)
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
PtrList< volScalarField > & Y
const word cloudName(propsDict.lookup("cloudName"))