31template<
class CloudType>
35 return *cloudCopyPtr_;
39template<
class CloudType>
40inline const typename CloudType::particleType::constantProperties&
47template<
class CloudType>
48inline typename CloudType::particleType::constantProperties&
55template<
class CloudType>
63template<
class CloudType>
71template<
class CloudType>
79template<
class CloudType>
83 const typename parcelType::trackingData&
td
88 const label celli =
p.cell();
90 const scalar m =
p.nParticle()*
p.mass();
92 this->rhokTrans()[celli] += m;
94 this->UTrans()[celli] += m*
p.U();
96 const scalar pc =
td.pc();
97 const scalar
T =
p.T();
98 const auto&
Y =
p.Y();
102 const scalar dm = m*
p.Y[i];
103 const label gid = comp.localToCarrierId(0, i);
104 const scalar hs = comp.carrier().Hs(gid, pc,
T);
106 this->rhoTrans(gid)[celli] += dm;
107 this->hsTrans()[celli] += dm*hs;
112template<
class CloudType>
120template<
class CloudType>
129template<
class CloudType>
137template<
class CloudType>
146 if (this->
solution().semiImplicit(
"Yi"))
155 auto& sourceField = trhoTrans.ref();
157 sourceField.primitiveFieldRef() =
158 rhoTrans_[i]/(this->db().time().deltaTValue()*this->
mesh().V());
163 fvm::Sp(
neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
164 +
pos0(sourceField)*sourceField;
168 auto tfvm = tmp<fvScalarMatrix>::New(Yi, dimMass/dimTime);
169 auto& fvm = tfvm.ref();
171 fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
181template<
class CloudType>
199 rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->
mesh().V());
206template<
class CloudType>
226 sourceField += rhoTrans_[i];
229 sourceField /= this->db().time().deltaTValue()*this->
mesh().V();
236template<
class CloudType>
251 if (this->
solution().semiImplicit(
"rho"))
257 sourceField /= this->
db().time().deltaTValue()*this->
mesh().V();
271 fvm.source() = -trhoTrans()/this->
db().time().deltaT();
277 return tmp<fvScalarMatrix>::New(
rho, dimMass/dimTime);
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
const objectRegistry & db() const noexcept
Return the local objectRegistry.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
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.
autoPtr< PhaseChangeModel< ReactingCloud< CloudType > > > phaseChangeModel_
Reacting phase change model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
PtrList< volScalarField::Internal > rhoTrans_
Mass transfer fields - one per carrier phase specie.
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.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
volScalarField::Internal & rhoTrans(const label i)
Return reference to mass source for field i.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
parcelType::constantProperties constProps_
Parcel constant properties.
autoPtr< CompositionModel< ReactingCloud< CloudType > > > compositionModel_
Reacting composition model.
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
Selector class for relaxation factors, solver type and solution.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
basicSpecieMixture & composition
PtrList< volScalarField > & Y
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace of functions to calculate implicit derivatives returning a matrix.
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.