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-2019 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  (
57  (
58  "invDsm",
59  phase_.mesh(),
61  )
62  );
63 
64  volScalarField& invDsm = tInvDsm.ref();
65 
66  forAll(sizeGroups_, i)
67  {
68  const sizeGroup& fi = sizeGroups_[i];
69 
70  invDsm += fi/fi.d();
71  }
72 
73  return 1.0/tInvDsm;
74 }
75 
76 
78 Foam::diameterModels::velocityGroup::fSum() const
79 {
80  tmp<volScalarField> tsumSizeGroups
81  (
83  (
84  "sumSizeGroups",
85  phase_.mesh(),
87  )
88  );
89 
90  volScalarField& sumSizeGroups = tsumSizeGroups.ref();
91 
92  forAll(sizeGroups_, i)
93  {
94  sumSizeGroups += sizeGroups_[i];
95  }
96 
97  return tsumSizeGroups;
98 }
99 
100 
101 void Foam::diameterModels::velocityGroup::renormalize()
102 {
103  Info<< "Renormalizing sizeGroups for velocityGroup "
104  << phase_.name()
105  << endl;
106 
107  // Set negative values to zero
108  forAll(sizeGroups_, i)
109  {
110  sizeGroups_[i] *= pos(sizeGroups_[i]);
111  };
112 
113  forAll(sizeGroups_, i)
114  {
115  sizeGroups_[i] /= fSum_;
116  };
117 }
118 
119 
121 Foam::diameterModels::velocityGroup::mvconvection() const
122 {
123  tmp<fv::convectionScheme<Foam::scalar> > mvConvection
124  (
126  (
127  phase_.mesh(),
128  fields_,
129  phase_.alphaRhoPhi(),
130  phase_.mesh().divScheme
131  (
132  "div(" + phase_.alphaRhoPhi()().name() + ",f)"
133  )
134  )
135  );
136 
137  return mvConvection;
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
144 (
145  const dictionary& diameterProperties,
146  const phaseModel& phase
147 )
148 :
149  diameterModel(diameterProperties, phase),
150  popBalName_(diameterProperties.lookup("populationBalance")),
151  f_
152  (
153  IOobject
154  (
155  IOobject::groupName
156  (
157  "f",
158  IOobject::groupName
159  (
160  phase.name(),
161  popBalName_
162  )
163  ),
164  phase.time().timeName(),
165  phase.mesh(),
166  IOobject::MUST_READ,
167  IOobject::AUTO_WRITE
168  ),
169  phase.mesh()
170  ),
171  formFactor_("formFactor", dimless, diameterProperties),
172  sizeGroups_
173  (
174  diameterProperties.lookup("sizeGroups"),
175  sizeGroup::iNew(phase, *this)
176  ),
177  fSum_
178  (
179  IOobject
180  (
181  IOobject::groupName
182  (
183  "fsum",
184  IOobject::groupName
185  (
186  phase.name(),
187  popBalName_
188  )
189  ),
190  phase.time().timeName(),
191  phase.mesh()
192  ),
193  fSum()
194  ),
195  d_
196  (
197  IOobject
198  (
199  IOobject::groupName("d", phase.name()),
200  phase.time().timeName(),
201  phase.mesh(),
202  IOobject::NO_READ,
203  IOobject::AUTO_WRITE
204  ),
205  phase.mesh(),
207  ),
208  dmdt_
209  (
210  IOobject
211  (
212  IOobject::groupName("source", phase.name()),
213  phase.time().timeName(),
214  phase.mesh()
215  ),
216  phase.mesh(),
218  )
219 {
220  if
221  (
222  phase_.mesh().solverDict(popBalName_).lookupOrDefault<Switch>
223  (
224  "renormalizeAtRestart",
225  false
226  )
227  ||
228  phase_.mesh().solverDict(popBalName_).lookupOrDefault<Switch>
229  (
230  "renormalize",
231  false
232  )
233  )
234  {
235  renormalize();
236  }
237 
238  fSum_ = fSum();
239 
240  if
241  (
242  mag(1 - fSum_.weightedAverage(fSum_.mesh().V()).value()) >= 1e-5
243  || mag(1 - max(fSum_).value()) >= 1e-5
244  || mag(1 - min(fSum_).value()) >= 1e-5
245  )
246  {
248  << " Initial values of the sizeGroups belonging to velocityGroup "
249  << this->phase().name()
250  << " must add to" << nl << " unity. This condition might be"
251  << " violated due to wrong entries in the" << nl
252  << " velocityGroupCoeffs subdictionary or bad initial conditions in"
253  << " the startTime" << nl
254  << " directory. The sizeGroups can be renormalized at every"
255  << " timestep or at restart" << nl
256  << " only by setting the corresponding switch renormalize or"
257  << " renormalizeAtRestart" << nl
258  << " in the fvSolution subdictionary " << popBalName_ << "."
259  << " Note that boundary conditions are not" << nl << "renormalized."
260  << exit(FatalError);
261  }
262 
263  forAll(sizeGroups_, i)
264  {
265  fields_.add(sizeGroups_[i]);
266  }
267 
268  d_ = dsm();
269 }
270 
271 
272 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
273 
275 {}
276 
277 
278 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
279 
280 
282 {
283  mvConvection_ = mvconvection();
284 }
285 
286 
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  fSum_ = fSum();
298 
299  Info<< phase_.name() << " sizeGroups-sum volume fraction, min, max = "
300  << fSum_.weightedAverage(phase_.mesh().V()).value()
301  << ' ' << min(fSum_).value()
302  << ' ' << max(fSum_).value()
303  << endl;
304 
305  if
306  (
307  phase_.mesh().solverDict(popBalName_).lookupOrDefault<Switch>
308  (
309  "renormalize",
310  false
311  )
312  )
313  {
314  renormalize();
315  }
316 }
317 
318 
320 read(const dictionary& phaseProperties)
321 {
322  diameterModel::read(phaseProperties);
323 
324  return true;
325 }
326 
327 
330 {
331  return d_;
332 }
333 
334 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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.