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>
81 {
82  return rhoTrans_[i];
83 }
84 
85 
86 template<class CloudType>
87 inline
90 {
91  return rhoTrans_;
92 }
93 
94 
95 template<class CloudType>
98 {
99  return rhoTrans_;
100 }
101 
102 
103 template<class CloudType>
105 (
106  const label i,
107  volScalarField& Yi
108 ) const
109 {
110  if (this->solution().coupled())
111  {
112  if (this->solution().semiImplicit("Yi"))
113  {
114  tmp<volScalarField> trhoTrans
115  (
116  new volScalarField
117  (
118  IOobject
119  (
120  this->name() + ":rhoTrans",
121  this->db().time().timeName(),
122  this->db(),
123  IOobject::NO_READ,
124  IOobject::NO_WRITE,
125  false
126  ),
127  this->mesh(),
129  )
130  );
131 
132  volScalarField& sourceField = trhoTrans.ref();
133 
134  sourceField.primitiveFieldRef() =
135  rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
136 
137  const dimensionedScalar YiSMALL("YiSMALL", dimless, SMALL);
138 
139  return
140  fvm::Sp(neg(sourceField)*sourceField/(Yi + YiSMALL), Yi)
141  + pos0(sourceField)*sourceField;
142  }
143  else
144  {
146  fvScalarMatrix& fvm = tfvm.ref();
147 
148  fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
149 
150  return tfvm;
151  }
152  }
153 
155 }
156 
157 
158 template<class CloudType>
161 {
163  (
165  (
166  IOobject
167  (
168  this->name() + ":rhoTrans",
169  this->db().time().timeName(),
170  this->db(),
171  IOobject::NO_READ,
172  IOobject::NO_WRITE,
173  false
174  ),
175  this->mesh(),
177  (
178  rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
179  )
180  )
181  );
182 
183  if (this->solution().coupled())
184  {
185  scalarField& rhoi = tRhoi.ref();
186  rhoi = rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
187  }
188 
189  return tRhoi;
190 }
191 
192 
193 template<class CloudType>
196 {
198  (
200  (
201  IOobject
202  (
203  this->name() + ":rhoTrans",
204  this->db().time().timeName(),
205  this->db(),
206  IOobject::NO_READ,
207  IOobject::NO_WRITE,
208  false
209  ),
210  this->mesh(),
212  (
213  rhoTrans_[0].dimensions()/dimTime/dimVolume, Zero
214  )
215  )
216  );
217 
218  if (this->solution().coupled())
219  {
220  scalarField& sourceField = trhoTrans.ref();
221  forAll(rhoTrans_, i)
222  {
223  sourceField += rhoTrans_[i];
224  }
225 
226  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
227  }
228 
229  return trhoTrans;
230 }
231 
232 
233 template<class CloudType>
236 {
237  if (this->solution().coupled())
238  {
239  tmp<volScalarField> trhoTrans
240  (
241  new volScalarField
242  (
243  IOobject
244  (
245  this->name() + ":rhoTrans",
246  this->db().time().timeName(),
247  this->db(),
248  IOobject::NO_READ,
249  IOobject::NO_WRITE,
250  false
251  ),
252  this->mesh(),
254  )
255  );
256 
257  scalarField& sourceField = trhoTrans.ref();
258 
259  if (this->solution().semiImplicit("rho"))
260  {
261 
262  forAll(rhoTrans_, i)
263  {
264  sourceField += rhoTrans_[i];
265  }
266  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
267 
268  return fvm::SuSp(trhoTrans()/rho, rho);
269  }
270  else
271  {
273  fvScalarMatrix& fvm = tfvm.ref();
274 
275  forAll(rhoTrans_, i)
276  {
277  sourceField += rhoTrans_[i];
278  }
279 
280  fvm.source() = -trhoTrans()/this->db().time().deltaT();
281 
282  return tfvm;
283  }
284  }
285 
287 }
288 
289 
290 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::PhaseChangeModel
Templated phase change model class.
Definition: ReactingCloud.H:61
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:55
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::ReactingCloud::SYi
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
Definition: ReactingCloudI.H:105
Foam::ReactingCloud::constProps
const parcelType::constantProperties & constProps() const
Return the constant properties.
Definition: ReactingCloudI.H:40
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fam::SuSp
tmp< faMatrix< Type > > SuSp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:151
Foam::CompositionModel
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
Definition: ReactingCloud.H:58
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::Field< scalar >
Foam::ReactingCloud::composition
const CompositionModel< ReactingCloud< CloudType > > & composition() const
Return const access to reacting composition model.
Definition: ReactingCloudI.H:56
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
timeName
word timeName
Definition: getTimeIndex.H:3
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::ReactingCloud::cloudCopy
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
Definition: ReactingCloudI.H:32
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:445
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::ReactingCloud::phaseChange
const PhaseChangeModel< ReactingCloud< CloudType > > & phaseChange() const
Return const access to reacting phase change model.
Definition: ReactingCloudI.H:64
Foam::ReactingCloud::Srho
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
Definition: ReactingCloudI.H:195
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::ReactingCloud::rhoTrans
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
Definition: ReactingCloudI.H:89
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::ReactingCloud
Templated base class for reacting cloud.
Definition: ReactingCloud.H:68