ReactingCloudI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
27 
28 template<class CloudType>
31 {
32  return cloudCopyPtr_();
33 }
34 
35 
36 template<class CloudType>
37 inline const typename CloudType::particleType::constantProperties&
39 {
40  return constProps_;
41 }
42 
43 
44 template<class CloudType>
45 inline typename CloudType::particleType::constantProperties&
47 {
48  return constProps_;
49 }
50 
51 
52 template<class CloudType>
55 {
56  return compositionModel_;
57 }
58 
59 
60 template<class CloudType>
63 {
64  return phaseChangeModel_;
65 }
66 
67 
68 template<class CloudType>
71 {
72  return phaseChangeModel_();
73 }
74 
75 
76 template<class CloudType>
79 {
80  return rhoTrans_[i];
81 }
82 
83 
84 template<class CloudType>
85 inline
88 {
89  return rhoTrans_;
90 }
91 
92 
93 template<class CloudType>
96 {
97  return rhoTrans_;
98 }
99 
100 
101 template<class CloudType>
103 (
104  const label i,
105  volScalarField& Yi
106 ) const
107 {
108  if (this->solution().coupled())
109  {
110  if (this->solution().semiImplicit("Yi"))
111  {
112  tmp<volScalarField> trhoTrans
113  (
114  new volScalarField
115  (
116  IOobject
117  (
118  this->name() + ":rhoTrans",
119  this->db().time().timeName(),
120  this->db(),
121  IOobject::NO_READ,
122  IOobject::NO_WRITE,
123  false
124  ),
125  this->mesh(),
127  )
128  );
129 
130  volScalarField& sourceField = trhoTrans.ref();
131 
132  sourceField.primitiveFieldRef() =
133  rhoTrans_[i]/(this->db().time().deltaTValue()*this->mesh().V());
134 
135  const dimensionedScalar Yismall("Yismall", dimless, small);
136 
137  return
138  fvm::Sp(neg(sourceField)*sourceField/(Yi + Yismall), Yi)
139  + pos0(sourceField)*sourceField;
140  }
141  else
142  {
144  fvScalarMatrix& fvm = tfvm.ref();
145 
146  fvm.source() = -rhoTrans_[i]/this->db().time().deltaTValue();
147 
148  return tfvm;
149  }
150  }
151 
153 }
154 
155 
156 template<class CloudType>
159 {
161  (
163  (
164  IOobject
165  (
166  this->name() + ":rhoTrans",
167  this->db().time().timeName(),
168  this->db(),
169  IOobject::NO_READ,
170  IOobject::NO_WRITE,
171  false
172  ),
173  this->mesh(),
175  (
176  "zero",
177  rhoTrans_[0].dimensions()/dimTime/dimVolume,
178  0.0
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  "zero",
214  rhoTrans_[0].dimensions()/dimTime/dimVolume,
215  0.0
216  )
217  )
218  );
219 
220  if (this->solution().coupled())
221  {
222  scalarField& sourceField = trhoTrans.ref();
223  forAll(rhoTrans_, i)
224  {
225  sourceField += rhoTrans_[i];
226  }
227 
228  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
229  }
230 
231  return trhoTrans;
232 }
233 
234 
235 template<class CloudType>
238 {
239  if (this->solution().coupled())
240  {
241  tmp<volScalarField> trhoTrans
242  (
243  new volScalarField
244  (
245  IOobject
246  (
247  this->name() + ":rhoTrans",
248  this->db().time().timeName(),
249  this->db(),
250  IOobject::NO_READ,
251  IOobject::NO_WRITE,
252  false
253  ),
254  this->mesh(),
256  )
257  );
258 
259  scalarField& sourceField = trhoTrans.ref();
260 
261  if (this->solution().semiImplicit("rho"))
262  {
263 
264  forAll(rhoTrans_, i)
265  {
266  sourceField += rhoTrans_[i];
267  }
268  sourceField /= this->db().time().deltaTValue()*this->mesh().V();
269 
270  return fvm::SuSp(trhoTrans()/rho, rho);
271  }
272  else
273  {
275  fvScalarMatrix& fvm = tfvm.ref();
276 
277  forAll(rhoTrans_, i)
278  {
279  sourceField += rhoTrans_[i];
280  }
281 
282  fvm.source() = -trhoTrans()/this->db().time().deltaT();
283 
284  return tfvm;
285  }
286  }
287 
289 }
290 
291 
292 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
tmp< GeometricField< Type, fvPatchField, volMesh > > SuSp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:102
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Templated phase change model class.
Definition: ReactingCloud.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const PtrList< volScalarField::Internal > & rhoTrans() const
Return const access to mass source fields.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
dimensionedScalar neg(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const ReactingCloud & cloudCopy() const
Return a reference to the cloud copy.
dynamicFvMesh & mesh
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:72
const parcelType::constantProperties & constProps() const
Return the constant properties.
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos0(const dimensionedScalar &ds)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
Field< Type > & source()
Definition: fvMatrix.H:294
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Templated base class for reacting cloud.
Definition: ReactingCloud.H:62
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
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:48
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
Templated reacting parcel composition model class Consists of carrier species (via thermo package)...
Definition: ReactingCloud.H:52
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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