30template<
class CloudType>
34 return *cloudCopyPtr_;
38template<
class CloudType>
46template<
class CloudType>
54template<
class CloudType>
58 return *compositionModel_;
62template<
class CloudType>
66 return *phaseChangeModel_;
70template<
class CloudType>
74 return *phaseChangeModel_;
78template<
class CloudType>
86template<
class CloudType>
95template<
class CloudType>
103template<
class CloudType>
112 if (this->
solution().semiImplicit(
"Yi"))
120 this->
name() +
":rhoTrans",
135 rhoTrans_[i]/(this->db().time().deltaTValue()*this->
mesh().V());
140 fvm::Sp(
neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
141 +
pos0(sourceField)*sourceField;
148 fvm.
source() = -rhoTrans_[i]/this->db().time().deltaTValue();
158template<
class CloudType>
168 this->
name() +
":rhoTrans",
186 rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->
mesh().V());
193template<
class CloudType>
203 this->
name() +
":rhoTrans",
213 rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
223 sourceField += rhoTrans_[i];
226 sourceField /= this->db().time().deltaTValue()*this->
mesh().V();
233template<
class CloudType>
245 this->
name() +
":rhoTrans",
259 if (this->
solution().semiImplicit(
"rho"))
264 sourceField += rhoTrans_[i];
266 sourceField /= this->db().time().deltaTValue()*this->
mesh().V();
268 return fvm::SuSp(trhoTrans()/
rho,
rho);
272 tmp<fvScalarMatrix> tfvm(
new fvScalarMatrix(
rho, dimMass/dimTime));
273 fvScalarMatrix& fvm = tfvm.ref();
277 sourceField += rhoTrans_[i];
280 fvm.source() = -trhoTrans()/this->db().time().deltaT();
286 return tmp<fvScalarMatrix>::New(
rho, dimMass/dimTime);
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Class to hold DSMC particle constant properties.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated phase change model class.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Templated base class for reacting cloud.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Field< Type > & source() noexcept
Selector class for relaxation factors, solver type and solution.
A class for managing temporary objects.
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
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)
#define forAll(list, i)
Loop across all elements in list.