Loading...
Searching...
No Matches
ReactingCloudI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30
31template<class CloudType>
35 return *cloudCopyPtr_;
36}
37
38
39template<class CloudType>
40inline const typename CloudType::particleType::constantProperties&
43 return constProps_;
44}
45
46
47template<class CloudType>
48inline typename CloudType::particleType::constantProperties&
51 return constProps_;
52}
53
54
55template<class CloudType>
59 return *compositionModel_;
60}
61
62
63template<class CloudType>
67 return *phaseChangeModel_;
68}
69
70
71template<class CloudType>
74{
75 return *phaseChangeModel_;
76}
77
78
79template<class CloudType>
81(
82 const parcelType& p,
83 const typename parcelType::trackingData& td
84)
85{
86 const auto& comp = this->composition();
87
88 const label celli = p.cell();
89
90 const scalar m = p.nParticle()*p.mass();
91
92 this->rhokTrans()[celli] += m;
93
94 this->UTrans()[celli] += m*p.U();
95
96 const scalar pc = td.pc();
97 const scalar T = p.T();
98 const auto& Y = p.Y();
99
100 forAll(Y, i)
101 {
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);
105
106 this->rhoTrans(gid)[celli] += dm;
107 this->hsTrans()[celli] += dm*hs;
108 }
109}
110
111
112template<class CloudType>
115{
116 return rhoTrans_[i];
117}
118
119
120template<class CloudType>
121inline
125 return rhoTrans_;
126}
127
128
129template<class CloudType>
132{
133 return rhoTrans_;
134}
135
136
137template<class CloudType>
139(
140 const label i,
142) const
143{
144 if (this->solution().coupled())
145 {
146 if (this->solution().semiImplicit("Yi"))
147 {
148 auto trhoTrans = volScalarField::New
149 (
150 IOobject::scopedName(this->name(), "rhoTrans"),
152 this->mesh(),
154 );
155 auto& sourceField = trhoTrans.ref();
156
157 sourceField.primitiveFieldRef() =
158 rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
159
160 const dimensionedScalar YiSMALL("YiSMALL", dimless, SMALL);
161
162 return
163 fvm::Sp(neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
164 + pos0(sourceField)*sourceField;
165 }
166 else
167 {
168 auto tfvm = tmp<fvScalarMatrix>::New(Yi, dimMass/dimTime);
169 auto& fvm = tfvm.ref();
170
171 fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
172
173 return tfvm;
174 }
175 }
178}
179
180
181template<class CloudType>
183Foam::ReactingCloud<CloudType>::Srho(const label i) const
184{
186 (
187 IOobject::scopedName(this->name(), "rhoTrans"),
189 this->mesh(),
191 (
192 rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
193 )
194 );
195 scalarField& rhoi = tRhoi.ref();
196
197 if (this->solution().coupled())
198 {
199 rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
200 }
202 return tRhoi;
203}
204
205
206template<class CloudType>
209{
210 auto trhoTrans = volScalarField::Internal::New
211 (
212 IOobject::scopedName(this->name(), "rhoTrans"),
214 this->mesh(),
216 (
217 rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
218 )
219 );
220 scalarField& sourceField = trhoTrans.ref();
221
222 if (this->solution().coupled())
223 {
224 forAll(rhoTrans_, i)
225 {
226 sourceField += rhoTrans_[i];
227 }
228
229 sourceField /= this->db().time().deltaTValue()*this->mesh().V();
230 }
232 return trhoTrans;
233}
234
235
236template<class CloudType>
239{
240 if (this->solution().coupled())
241 {
242 auto trhoTrans = volScalarField::New
243 (
244 IOobject::scopedName(this->name(), "rhoTrans"),
246 this->mesh(),
248 );
249 scalarField& sourceField = trhoTrans.ref();
251 if (this->solution().semiImplicit("rho"))
252 {
253 forAll(rhoTrans_, i)
254 {
255 sourceField += rhoTrans_[i];
256 }
257 sourceField /= this->db().time().deltaTValue()*this->mesh().V();
258
259 return fvm::SuSp(trhoTrans()/rho, rho);
260 }
261 else
262 {
264 auto& fvm = tfvm.ref();
265
266 forAll(rhoTrans_, i)
267 {
268 sourceField += rhoTrans_[i];
269 }
271 fvm.source() = -trhoTrans()/this->db().time().deltaT();
272
273 return tfvm;
274 }
275 }
276
277 return tmp<fvScalarMatrix>::New(rho, dimMass/dimTime);
278}
279
280
281// ************************************************************************* //
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.
Definition IOobject.C:450
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Templated phase change model class.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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.
Definition solution.H:95
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
basicSpecieMixture & composition
volScalarField & p
PtrList< volScalarField > & Y
bool coupled
dynamicFvMesh & mesh
auto & name
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).
Definition zero.H:127
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.
Definition stdFoam.H:299