velocityGroup.C
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) 2017-2019 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "velocityGroup.H"
30#include "sizeGroup.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace diameterModels
40{
42
44 (
48 );
49}
50}
51
52
53// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54
55Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::dsm() const
56{
57 tmp<volScalarField> tInvDsm
58 (
60 (
61 "invDsm",
62 phase_.mesh(),
64 )
65 );
66
67 volScalarField& invDsm = tInvDsm.ref();
68
69 forAll(sizeGroups_, i)
70 {
71 const sizeGroup& fi = sizeGroups_[i];
72
73 invDsm += fi/fi.d();
74 }
75
76 return 1.0/tInvDsm;
77}
78
79
81Foam::diameterModels::velocityGroup::fSum() const
82{
83 tmp<volScalarField> tsumSizeGroups
84 (
86 (
87 "sumSizeGroups",
88 phase_.mesh(),
89 dimensionedScalar("zero", dimless, 0)
90 )
91 );
92
93 volScalarField& sumSizeGroups = tsumSizeGroups.ref();
94
95 forAll(sizeGroups_, i)
96 {
97 sumSizeGroups += sizeGroups_[i];
98 }
99
100 return tsumSizeGroups;
101}
102
103
104void Foam::diameterModels::velocityGroup::renormalize()
105{
106 Info<< "Renormalizing sizeGroups for velocityGroup "
107 << phase_.name()
108 << endl;
109
110 // Set negative values to zero
111 forAll(sizeGroups_, i)
112 {
113 sizeGroups_[i] *= pos(sizeGroups_[i]);
114 };
115
116 forAll(sizeGroups_, i)
117 {
118 sizeGroups_[i] /= fSum_;
119 };
120}
121
122
124Foam::diameterModels::velocityGroup::mvconvection() const
125{
126 tmp<fv::convectionScheme<Foam::scalar>> mvConvection
127 (
129 (
130 phase_.mesh(),
131 fields_,
132 phase_.alphaRhoPhi(),
133 phase_.mesh().divScheme
134 (
135 "div(" + phase_.alphaRhoPhi()().name() + ",f)"
136 )
137 )
138 );
139
140 return mvConvection;
141}
142
143
144// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
145
147(
148 const dictionary& diameterProperties,
149 const phaseModel& phase
150)
151:
152 diameterModel(diameterProperties, phase),
153 popBalName_(diameterProperties.lookup("populationBalance")),
154 f_
155 (
157 (
158 IOobject::groupName
159 (
160 "f",
161 IOobject::groupName
162 (
163 phase.name(),
164 popBalName_
165 )
166 ),
167 phase.time().timeName(),
168 phase.mesh(),
169 IOobject::MUST_READ,
170 IOobject::AUTO_WRITE
171 ),
172 phase.mesh()
173 ),
174 formFactor_("formFactor", dimless, diameterProperties),
175 sizeGroups_
176 (
177 diameterProperties.lookup("sizeGroups"),
178 sizeGroup::iNew(phase, *this)
179 ),
180 fSum_
181 (
183 (
184 IOobject::groupName
185 (
186 "fsum",
187 IOobject::groupName
188 (
189 phase.name(),
190 popBalName_
191 )
192 ),
193 phase.time().timeName(),
194 phase.mesh()
195 ),
196 fSum()
197 ),
198 d_
199 (
201 (
202 IOobject::groupName("d", phase.name()),
203 phase.time().timeName(),
204 phase.mesh(),
205 IOobject::NO_READ,
206 IOobject::AUTO_WRITE
207 ),
208 phase.mesh(),
210 ),
211 dmdt_
212 (
214 (
215 IOobject::groupName("source", phase.name()),
216 phase.time().timeName(),
217 phase.mesh()
218 ),
219 phase.mesh(),
221 )
222{
223 if
224 (
226 (
227 "renormalizeAtRestart",
228 false
229 )
230 ||
232 (
233 "renormalize",
234 false
235 )
236 )
237 {
238 renormalize();
239 }
240
241 fSum_ = fSum();
242
243 if
244 (
245 mag(1 - fSum_.weightedAverage(fSum_.mesh().V()).value()) >= 1e-5
246 || mag(1 - max(fSum_).value()) >= 1e-5
247 || mag(1 - min(fSum_).value()) >= 1e-5
248 )
249 {
251 << " Initial values of the sizeGroups belonging to velocityGroup "
252 << this->phase().name()
253 << " must add to" << nl << " unity. This condition might be"
254 << " violated due to wrong entries in the" << nl
255 << " velocityGroupCoeffs subdictionary or bad initial conditions in"
256 << " the startTime" << nl
257 << " directory. The sizeGroups can be renormalized at every"
258 << " timestep or at restart" << nl
259 << " only by setting the corresponding switch renormalize or"
260 << " renormalizeAtRestart" << nl
261 << " in the fvSolution subdictionary " << popBalName_ << "."
262 << " Note that boundary conditions are not" << nl << "renormalized."
263 << exit(FatalError);
264 }
265
266 forAll(sizeGroups_, i)
267 {
268 fields_.add(sizeGroups_[i]);
269 }
270
271 d_ = dsm();
272}
273
274
275// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
276
278{}
279
280
281// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
282
283
285{
286 mvConvection_ = mvconvection();
287}
288
289
291{
292 d_ = dsm();
293
294 Info<< this->phase().name() << " Sauter mean diameter, min, max = "
295 << d_.weightedAverage(d_.mesh().V()).value()
296 << ' ' << min(d_).value()
297 << ' ' << max(d_).value()
298 << endl;
299
300 fSum_ = fSum();
301
302 Info<< phase_.name() << " sizeGroups-sum volume fraction, min, max = "
303 << fSum_.weightedAverage(phase_.mesh().V()).value()
304 << ' ' << min(fSum_).value()
305 << ' ' << max(fSum_).value()
306 << endl;
307
308 if
309 (
310 phase_.mesh().solverDict(popBalName_).getOrDefault<Switch>
311 (
312 "renormalize",
313 false
314 )
315 )
316 {
317 renormalize();
318 }
319}
320
321
324{
326
327 return true;
328}
329
330
333{
334 return d_;
335}
336
337// ************************************************************************* //
tmp< fv::convectionScheme< scalar > > mvConvection(fv::convectionScheme< scalar >::New(mesh, fields, phi, mesh.divScheme("div(phi,Yi_h)")))
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &weightField) const
Calculate and return weighted average.
const Mesh & mesh() const
Return mesh.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual bool read()
Re-read model coefficients if they have changed.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:54
const phaseModel & phase() const
Return the phase.
const phaseModel & phase_
Definition: diameterModel.H:61
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:99
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
void postSolve()
Corrections after populationBalanceModel::solve()
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
virtual ~velocityGroup()
Destructor.
void preSolve()
Corrections before populationBalanceModel::solve()
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const word & name() const
Definition: phaseModel.H:144
Helper class to manage multi-specie phase properties.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
const word & name() const
Definition: phase.H:111
Lookup type of boundary radiation properties.
Definition: lookup.H:66
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:366
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333