Go to the documentation of this file.
31 template<
class ParcelType>
36 pMin_(this->dict_, 0.0),
37 constantVolume_(this->dict_, false),
38 volUpdateType_(this->dict_, mUndefined)
42 template<
class ParcelType>
48 ParcelType::constantProperties(
cp),
50 constantVolume_(
cp.constantVolume_),
51 volUpdateType_(
cp.volUpdateType_)
55 template<
class ParcelType>
61 ParcelType::constantProperties(parentDict),
62 pMin_(this->dict_,
"pMin", 1000.0),
63 constantVolume_(this->dict_,
"constantVolume",
false),
64 volUpdateType_(this->dict_,
"volumeUpdateMethod")
69 if (this->dict_.found(
"constantVolume"))
71 volUpdateType_.setValue(mUndefined);
73 else if (this->dict_.readIfPresent(
"volumeUpdateMethod", updateName))
76 if (updateName ==
"constantRho")
78 volUpdateType_.setValue(mConstRho);
80 else if (updateName ==
"constantVolume")
82 volUpdateType_.setValue(mConstVol);
84 else if (updateName ==
"updateRhoAndVol")
86 volUpdateType_.setValue(mUpdateRhoAndVol);
93 <<
"Unknown volumeUpdateMethod type " << updateName
94 <<
"\n\nValid volumeUpdateMethod types :\n"
95 <<
"(constantRho constantVolume updateRhoAndVol)" <<
nl
101 constantVolume_.setValue(
false);
106 template<
class ParcelType>
112 const label tetFacei,
122 template<
class ParcelType>
130 ParcelType(
mesh, position, celli),
136 template<
class ParcelType>
142 const label tetFacei,
145 const scalar nParticle0,
147 const scalar dTarget0,
150 const vector& angularMomentum0,
183 template<
class ParcelType>
187 return pMin_.value();
191 template<
class ParcelType>
195 return constantVolume_.value();
199 template<
class ParcelType>
203 return volUpdateType_.value();
209 template<
class ParcelType>
216 template<
class ParcelType>
223 template<
class ParcelType>
230 template<
class ParcelType>
238 template<
class ParcelType>
246 template<
class ParcelType>
253 template<
class ParcelType>
scalar mass0() const
Return const access to initial mass [kg].
const scalarField & YSolid() const
Return const access to mass fractions of solids.
A class for handling words, derived from Foam::string.
const scalarField & Y() const
Return const access to mass fractions of mixture [].
Mesh consisting of general polyhedral cells.
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
PtrList< coordinateSystem > coordinates(solidRegions.size())
label volUpdateType() const
Return const access to the constant volume flag.
const scalarField & YGas() const
Return const access to mass fractions of gases.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar mass0_
Initial mass [kg].
errorManipArg< error, int > exit(error &err, const int errNo=1)
bool constantVolume() const
Return const access to the constant volume flag.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Class to hold reacting parcel constant properties.
constantProperties()
Null constructor.
scalar pMin() const
Return const access to the minimum pressure.
scalarField Y_
Mass fractions of mixture [].
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
scalarList Y0(nSpecie, Zero)