All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
velocityGroup.C
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) 2017-2020 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 #include "velocityGroup.H"
27 #include "sizeGroup.H"
28 #include "populationBalanceModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38  defineTypeNameAndDebug(velocityGroup, 0);
39  addToRunTimeSelectionTable(diameterModel, velocityGroup, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 
46 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::dsm() const
47 {
48  tmp<volScalarField> tInvDsm
49  (
51  (
52  "invDsm",
53  phase().mesh(),
55  )
56  );
57 
58  volScalarField& invDsm = tInvDsm.ref();
59 
60  forAll(sizeGroups_, i)
61  {
62  const sizeGroup& fi = sizeGroups_[i];
63 
64  invDsm += fi.a()*fi/fi.x();
65  }
66 
67  return 6.0/tInvDsm;
68 }
69 
70 
71 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::N() const
72 {
73  tmp<volScalarField> tN
74  (
76  (
77  "N",
78  phase().mesh(),
80  )
81  );
82 
83  volScalarField& N = tN.ref();
84 
85  forAll(sizeGroups_, i)
86  {
87  N += phase()*sizeGroups_[i]/sizeGroups_[i].x();
88  }
89 
90  return tN;
91 }
92 
93 
95 Foam::diameterModels::velocityGroup::fSum() const
96 {
97  tmp<volScalarField> tsumSizeGroups
98  (
100  (
101  "sumSizeGroups",
102  phase().mesh(),
104  )
105  );
106 
107  volScalarField& sumSizeGroups = tsumSizeGroups.ref();
108 
109  forAll(sizeGroups_, i)
110  {
111  sumSizeGroups += sizeGroups_[i];
112  }
113 
114  return tsumSizeGroups;
115 }
116 
117 
118 void Foam::diameterModels::velocityGroup::scale()
119 {
120  Info<< "Scaling sizeGroups for velocityGroup " << phase().name() << endl;
121 
122  forAll(sizeGroups_, i)
123  {
124  sizeGroups_[i].max(0);
125  };
126 
127  f_ = fSum();
128 
129  forAll(sizeGroups_, i)
130  {
131  sizeGroups_[i] /= f_;
132 
133  sizeGroups_[i].correctBoundaryConditions();
134  };
135 }
136 
137 
139 Foam::diameterModels::velocityGroup::mvconvection() const
140 {
141  tmp<fv::convectionScheme<Foam::scalar>> mvConvection
142  (
144  (
145  phase().mesh(),
146  fields_,
147  phase().alphaRhoPhi(),
148  phase().mesh().divScheme
149  (
150  "div(" + phase().alphaRhoPhi()().name() + ",f)"
151  )
152  )
153  );
154 
155  return mvConvection;
156 }
157 
158 
159 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
160 
163 {
164  return d_;
165 }
166 
167 
170 {
171  tmp<volScalarField> tA
172  (
174  (
175  "a",
176  phase().mesh(),
178  )
179  );
180 
181  volScalarField& a = tA.ref();
182 
183  forAll(sizeGroups_, i)
184  {
185  const sizeGroup& fi = sizeGroups_[i];
186 
187  a += fi.a()*fi/fi.x();
188  }
189 
190  return phase()*a;
191 }
192 
193 
194 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
195 
197 (
198  const dictionary& diameterProperties,
199  const phaseModel& phase
200 )
201 :
202  diameterModel(diameterProperties, phase),
203  popBalName_(diameterProperties.lookup("populationBalance")),
204  f_
205  (
206  IOobject
207  (
208  IOobject::groupName
209  (
210  "f",
211  IOobject::groupName
212  (
213  phase.name(),
214  popBalName_
215  )
216  ),
217  phase.time().timeName(),
218  phase.mesh(),
219  IOobject::MUST_READ,
220  IOobject::AUTO_WRITE
221  ),
222  phase.mesh()
223  ),
224  sizeGroups_
225  (
226  diameterProperties.lookup("sizeGroups"),
227  sizeGroup::iNew(phase, *this)
228  ),
229  d_(dRef()),
230  dmdt_
231  (
232  IOobject
233  (
234  IOobject::groupName("source", phase.name()),
235  phase.time().timeName(),
236  phase.mesh()
237  ),
238  phase.mesh(),
240  )
241 {
242  forAll(sizeGroups_, i)
243  {
244  fields_.add(sizeGroups_[i]);
245  }
246 
247  d_ = dsm();
248 }
249 
250 
251 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
252 
254 {}
255 
256 
257 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
258 
260 {
261  mvConvection_ = mvconvection();
262 }
263 
264 
266 {
267  if
268  (
269  phase().mesh().solverDict(popBalName_).lookupOrDefault<Switch>
270  (
271  "scale",
272  true
273  )
274  )
275  {
276  scale();
277  }
278 
279  f_ = fSum();
280 
281  f_.correctBoundaryConditions();
282 
283  Info<< phase().name() << " sizeGroups-sum volume fraction, min, max = "
284  << f_.weightedAverage(phase().mesh().V()).value()
285  << ' ' << min(f_).value()
286  << ' ' << max(f_).value()
287  << endl;
288 
289  d_ = dsm();
290 
291  Info<< this->phase().name() << " Sauter mean diameter, min, max = "
292  << d_.weightedAverage(d_.mesh().V()).value()
293  << ' ' << min(d_).value()
294  << ' ' << max(d_).value()
295  << endl;
296 }
297 
298 
300 (
301  const dictionary& phaseProperties
302 )
303 {
304  diameterModel::read(phaseProperties);
305 
306  return true;
307 }
308 
309 
310 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const word & name() const
Definition: phaseModel.H:109
virtual tmp< volScalarField > calcD() const
Get the diameter field.
velocityGroup(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
virtual ~velocityGroup()
Destructor.
word timeName
Definition: getTimeIndex.H:3
void min(const dimensioned< Type > &)
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > calcA() const
Get the surface area per unit volume field.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimDensity
Internal & ref()
Return a reference to the dimensioned internal field.
void preSolve()
Corrections before populationBalanceModel::solve()
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void max(const dimensioned< Type > &)
void postSolve()
Corrections after populationBalanceModel::solve()
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
const autoPtr< label > & i() const
Return label of the sizeGroup within the population balance.
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.