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-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 #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 
41  (
42  diameterModel,
43  velocityGroup,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
51 
52 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::dsm() const
53 {
54  tmp<volScalarField> tInvDsm
55  (
56  new volScalarField
57  (
58  IOobject
59  (
60  "invDsm",
61  phase_.time().timeName(),
62  phase_.mesh()
63  ),
64  phase_.mesh(),
65  dimensionedScalar("invDsm", inv(dimLength), Zero)
66  )
67  );
68 
69  volScalarField& invDsm = tInvDsm.ref();
70 
71  forAll(sizeGroups_, i)
72  {
73  const sizeGroup& fi = sizeGroups_[i];
74 
75  invDsm += fi/fi.d();
76  }
77 
78  return 1.0/tInvDsm;
79 }
80 
81 
83 Foam::diameterModels::velocityGroup::fSum() const
84 {
85  tmp<volScalarField> tsumSizeGroups
86  (
87  new volScalarField
88  (
89  IOobject
90  (
91  "sumSizeGroups",
92  phase_.time().timeName(),
93  phase_.mesh()
94  ),
95  phase_.mesh(),
96  dimensionedScalar("sumSizeGroups", dimless, 0)
97  )
98  );
99 
100  volScalarField& sumSizeGroups = tsumSizeGroups.ref();
101 
102  forAll(sizeGroups_, i)
103  {
104  sumSizeGroups += sizeGroups_[i];
105  }
106 
107  return tsumSizeGroups;
108 }
109 
110 
111 void Foam::diameterModels::velocityGroup::renormalize()
112 {
113  Info<< phase_.name()
114  << " renormalizing sizeGroups"
115  << endl;
116 
117  // Set negative values to zero
118  forAll(sizeGroups_, i)
119  {
120  sizeGroups_[i] *= pos(sizeGroups_[i]);
121  };
122 
123  forAll(sizeGroups_, i)
124  {
125  sizeGroups_[i] /= fSum_;
126  };
127 }
128 
129 
131 Foam::diameterModels::velocityGroup::mvconvection() const
132 {
133  tmp<fv::convectionScheme<Foam::scalar> > mvConvection
134  (
136  (
137  phase_.mesh(),
138  fields_,
139  phase_.alphaRhoPhi(),
140  phase_.mesh().divScheme
141  (
142  "div(" + phase_.alphaRhoPhi()().name() + ",f)"
143  )
144  )
145  );
146 
147  return mvConvection;
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
152 
154 (
155  const dictionary& diameterProperties,
156  const phaseModel& phase
157 )
158 :
159  diameterModel(diameterProperties, phase),
160  popBalName_(diameterProperties.lookup("populationBalance")),
161  f_
162  (
163  IOobject
164  (
165  IOobject::groupName
166  (
167  "f",
168  IOobject::groupName
169  (
170  phase.name(),
171  popBalName_
172  )
173  ),
174  phase.time().timeName(),
175  phase.mesh(),
176  IOobject::MUST_READ,
177  IOobject::AUTO_WRITE
178  ),
179  phase.mesh()
180  ),
181  formFactor_("formFactor", dimless, diameterProperties),
182  sizeGroups_
183  (
184  diameterProperties.lookup("sizeGroups"),
185  sizeGroup::iNew(phase, *this)
186  ),
187  fSum_
188  (
189  IOobject
190  (
191  IOobject::groupName
192  (
193  "fsum",
194  IOobject::groupName
195  (
196  phase.name(),
197  popBalName_
198  )
199  ),
200  phase.time().timeName(),
201  phase.mesh()
202  ),
203  fSum()
204  ),
205  d_
206  (
207  IOobject
208  (
209  IOobject::groupName("d", phase.name()),
210  phase.time().timeName(),
211  phase.mesh(),
212  IOobject::NO_READ,
213  IOobject::AUTO_WRITE
214  ),
215  phase.mesh(),
217  ),
218  dmdt_
219  (
220  IOobject
221  (
222  IOobject::groupName("source", phase.name()),
223  phase.time().timeName(),
224  phase.mesh()
225  ),
226  phase.mesh(),
228  )
229 {
230  if
231  (
232  phase_.mesh().solverDict(popBalName_).lookupOrDefault<Switch>
233  (
234  "renormalizeAtRestart",
235  false
236  )
237  ||
238  phase_.mesh().solverDict(popBalName_).lookupOrDefault<Switch>
239  (
240  "renormalize",
241  false
242  )
243  )
244  {
245  renormalize();
246  }
247 
248  fSum_ = fSum();
249 
250  if
251  (
252  mag(1 - fSum_.weightedAverage(fSum_.mesh().V()).value()) >= 1e-5
253  || mag(1 - max(fSum_).value()) >= 1e-5
254  || mag(1 - min(fSum_).value()) >= 1e-5
255  )
256  {
258  << " Initial values of the sizeGroups belonging to velocityGroup "
259  << this->phase().name()
260  << " must add to" << nl << " unity. This condition might be"
261  << " violated due to wrong entries in the" << nl
262  << " velocityGroupCoeffs subdictionary or bad initial conditions in"
263  << " the startTime" << nl
264  << " directory. The sizeGroups can be renormalized at every"
265  << " timestep or at restart" << nl
266  << " only by setting the corresponding switch renormalize or"
267  << " renormalizeAtRestart" << nl
268  << " in the fvSolution subdictionary " << popBalName_ << "."
269  << " Note that boundary conditions are not" << nl << "renormalized."
270  << exit(FatalError);
271  }
272 
273  forAll(sizeGroups_, i)
274  {
275  fields_.add(sizeGroups_[i]);
276  }
277 
278  d_ = dsm();
279 }
280 
281 
282 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
283 
285 {}
286 
287 
288 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
289 
290 
292 {
293  mvConvection_ = mvconvection();
294 }
295 
296 
298 {
299  d_ = dsm();
300 
301  Info<< this->phase().name() << " Sauter mean diameter, min, max = "
302  << d_.weightedAverage(d_.mesh().V()).value()
303  << ' ' << min(d_).value()
304  << ' ' << max(d_).value()
305  << endl;
306 
307  fSum_ = fSum();
308 
309  Info<< phase_.name() << " sizeGroups-sum volume fraction, min, max = "
310  << fSum_.weightedAverage(phase_.mesh().V()).value()
311  << ' ' << min(fSum_).value()
312  << ' ' << max(fSum_).value()
313  << endl;
314 
315  if
316  (
317  phase_.mesh().solverDict(popBalName_).lookupOrDefault<Switch>
318  (
319  "renormalize",
320  false
321  )
322  )
323  {
324  renormalize();
325  }
326 }
327 
328 
330 read(const dictionary& phaseProperties)
331 {
332  diameterModel::read(phaseProperties);
333 
334  return true;
335 }
336 
337 
340 {
341  return d_;
342 }
343 
344 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const word & name() const
Definition: phaseModel.H:109
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:52
dimensionedScalar pos(const dimensionedScalar &ds)
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
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)
static const char nl
Definition: Ostream.H:265
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
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
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.