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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class CloudType>
33 {
34  return *cloudCopyPtr_;
35 }
36 
37 
38 template<class CloudType>
39 inline const typename CloudType::particleType::constantProperties&
41 {
42  return constProps_;
43 }
44 
45 
46 template<class CloudType>
47 inline typename CloudType::particleType::constantProperties&
49 {
50  return constProps_;
51 }
52 
53 
54 template<class CloudType>
57 {
58  return *compositionModel_;
59 }
60 
61 
62 template<class CloudType>
65 {
66  return *phaseChangeModel_;
67 }
68 
69 
70 template<class CloudType>
73 {
74  return *phaseChangeModel_;
75 }
76 
77 
78 template<class CloudType>
80 (
81  const parcelType& p,
82  const typename parcelType::trackingData& td
83 )
84 {
85  const auto& comp = this->composition();
86 
87  const label celli = p.cell();
88 
89  const scalar m = p.nParticle()*p.mass();
90 
91  this->rhokTrans()[celli] += m;
92 
93  this->UTrans()[celli] += m*p.U();
94 
95  const scalar pc = td.pc();
96  const scalar T = p.T();
97  const auto& Y = p.Y();
98 
99  forAll(Y, i)
100  {
101  const scalar dm = m*p.Y[i];
102  const label gid = comp.localToCarrierId(0, i);
103  const scalar hs = comp.carrier().Hs(gid, pc, T);
104 
105  this->rhoTrans(gid)[celli] += dm;
106  this->hsTrans()[celli] += dm*hs;
107  }
108 }
109 
110 
111 template<class CloudType>
114 {
115  return rhoTrans_[i];
116 }
117 
118 
119 template<class CloudType>
120 inline
123 {
124  return rhoTrans_;
125 }
126 
127 
128 template<class CloudType>
131 {
132  return rhoTrans_;
133 }
134 
135 
136 template<class CloudType>
138 (
139  const label i,
140  volScalarField& Yi
141 ) const
142 {
143  if (this->solution().coupled())
144  {
145  if (this->solution().semiImplicit("Yi"))
146  {
147  tmp<volScalarField> trhoTrans
148  (
149  new volScalarField
150  (
151  IOobject
152  (
153  this->name() + ":rhoTrans",
154  this->db().time().timeName(),
155  this->db(),
156  IOobject::NO_READ,
157  IOobject::NO_WRITE,
158  IOobject::NO_REGISTER
159  ),
160  this->mesh(),
162  )
163  );
164 
165  volScalarField& sourceField = trhoTrans.ref();
166 
167  sourceField.primitiveFieldRef() =
168  rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
169 
170  const dimensionedScalar YiSMALL("YiSMALL", dimless, SMALL);
171 
172  return
173  fvm::Sp(neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
174  + pos0(sourceField)*sourceField;
175  }
176  else
177  {
178  tmp<fvScalarMatrix> tfvm(new fvScalarMatrix(Yi, dimMass/dimTime));
179  fvScalarMatrix& fvm = tfvm.ref();
180 
181  fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
182 
183  return tfvm;
184  }
185  }
188 }
189 
190 
191 template<class CloudType>
193 Foam::ReactingCloud<CloudType>::Srho(const label i) const
194 {
196  (
198  (
199  IOobject
200  (
201  this->name() + ":rhoTrans",
202  this->db().time().timeName(),
203  this->db(),
204  IOobject::NO_READ,
205  IOobject::NO_WRITE,
206  IOobject::NO_REGISTER
207  ),
208  this->mesh(),
210  (
211  rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
212  )
213  )
214  );
215 
216  if (this->solution().coupled())
217  {
218  scalarField& rhoi = tRhoi.ref();
219  rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
220  }
222  return tRhoi;
223 }
224 
225 
226 template<class CloudType>
229 {
231  (
233  (
234  IOobject
235  (
236  this->name() + ":rhoTrans",
237  this->db().time().timeName(),
238  this->db(),
239  IOobject::NO_READ,
240  IOobject::NO_WRITE,
241  IOobject::NO_REGISTER
242  ),
243  this->mesh(),
245  (
246  rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
247  )
248  )
249  );
250 
251  if (this->solution().coupled())
252  {
253  scalarField& sourceField = trhoTrans.ref();
254  forAll(rhoTrans_, i)
255  {
256  sourceField += rhoTrans_[i];
257  }
258 
259  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
260  }
262  return trhoTrans;
263 }
264 
265 
266 template<class CloudType>
269 {
270  if (this->solution().coupled())
271  {
272  tmp<volScalarField> trhoTrans
273  (
274  new volScalarField
275  (
276  IOobject
277  (
278  this->name() + ":rhoTrans",
279  this->db().time().timeName(),
280  this->db(),
281  IOobject::NO_READ,
282  IOobject::NO_WRITE,
283  IOobject::NO_REGISTER
284  ),
285  this->mesh(),
287  )
288  );
289 
290  scalarField& sourceField = trhoTrans.ref();
291 
292  if (this->solution().semiImplicit("rho"))
293  {
294 
295  forAll(rhoTrans_, i)
296  {
297  sourceField += rhoTrans_[i];
298  }
299  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
300 
301  return fvm::SuSp(trhoTrans()/rho, rho);
302  }
303  else
304  {
306  fvScalarMatrix& fvm = tfvm.ref();
307 
308  forAll(rhoTrans_, i)
309  {
310  sourceField += rhoTrans_[i];
311  }
312 
313  fvm.source() = -trhoTrans()/this->db().time().deltaT();
314 
315  return tfvm;
316  }
317  }
318 
320 }
321 
322 
323 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
Templated phase change model class.
Definition: ReactingCloud.H:57
basicSpecieMixture & composition
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
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
word timeName
Definition: getTimeIndex.H:3
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
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
const parcelType::constantProperties & constProps() const
Return the constant properties.
dimensionedScalar pos0(const dimensionedScalar &ds)
const volScalarField & T
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
Templated base class for reacting cloud.
Definition: ReactingCloud.H:64
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
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...
Definition: areaFieldsFwd.H:42
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.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
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
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
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