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 -------------------------------------------------------------------------------
11 License
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/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<class CloudType>
34 {
35  return *cloudCopyPtr_;
36 }
37 
38 
39 template<class CloudType>
40 inline const typename CloudType::particleType::constantProperties&
42 {
43  return constProps_;
44 }
45 
46 
47 template<class CloudType>
48 inline typename CloudType::particleType::constantProperties&
50 {
51  return constProps_;
52 }
53 
54 
55 template<class CloudType>
58 {
59  return *compositionModel_;
60 }
61 
62 
63 template<class CloudType>
66 {
67  return *phaseChangeModel_;
68 }
69 
70 
71 template<class CloudType>
74 {
75  return *phaseChangeModel_;
76 }
77 
78 
79 template<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 
112 template<class CloudType>
115 {
116  return rhoTrans_[i];
117 }
118 
119 
120 template<class CloudType>
121 inline
124 {
125  return rhoTrans_;
126 }
127 
128 
129 template<class CloudType>
132 {
133  return rhoTrans_;
134 }
135 
136 
137 template<class CloudType>
139 (
140  const label i,
141  volScalarField& Yi
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"),
151  IOobject::NO_REGISTER,
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 
181 template<class CloudType>
183 Foam::ReactingCloud<CloudType>::Srho(const label i) const
184 {
185  auto tRhoi = volScalarField::Internal::New
186  (
187  IOobject::scopedName(this->name(), "rhoTrans"),
188  IOobject::NO_REGISTER,
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 
206 template<class CloudType>
209 {
210  auto trhoTrans = volScalarField::Internal::New
211  (
212  IOobject::scopedName(this->name(), "rhoTrans"),
213  IOobject::NO_REGISTER,
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 
236 template<class CloudType>
239 {
240  if (this->solution().coupled())
241  {
242  auto trhoTrans = volScalarField::New
243  (
244  IOobject::scopedName(this->name(), "rhoTrans"),
245  IOobject::NO_REGISTER,
246  this->mesh(),
248  );
249  scalarField& sourceField = trhoTrans.ref();
250 
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  }
270 
271  fvm.source() = -trhoTrans()/this->db().time().deltaT();
272 
273  return tfvm;
274  }
275  }
276 
278 }
279 
280 
281 // ************************************************************************* //
Templated phase change model class.
Definition: ReactingCloud.H:57
basicSpecieMixture & composition
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic dimensioned Type class.
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
const dimensionSet dimless
Dimensionless.
zeroField SuSp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
dimensionedScalar neg(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const parcelType::constantProperties & constProps() const
Return the constant properties.
dimensionedScalar pos0(const dimensionedScalar &ds)
const volScalarField & T
Templated base class for reacting cloud.
Definition: ReactingCloud.H:64
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
PtrList< volScalarField > & Y
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
bool coupled
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ReactingCloud.H:54
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127