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-2022 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"
28 #include "populationBalanceModel.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36  defineTypeNameAndDebug(velocityGroup, 0);
37  addToRunTimeSelectionTable(diameterModel, velocityGroup, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
43 
44 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::dsm() const
45 {
46  tmp<volScalarField> tInvDsm
47  (
49  (
50  "invDsm",
51  phase().mesh(),
53  )
54  );
55 
56  volScalarField& invDsm = tInvDsm.ref();
57 
58  forAll(sizeGroups_, i)
59  {
60  const sizeGroup& fi = sizeGroups_[i];
61 
62  invDsm += fi.a()*fi/fi.x();
63  }
64 
65  return 6/tInvDsm;
66 }
67 
68 
69 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::N() const
70 {
71  tmp<volScalarField> tN
72  (
74  (
75  "N",
76  phase().mesh(),
78  )
79  );
80 
81  volScalarField& N = tN.ref();
82 
83  forAll(sizeGroups_, i)
84  {
85  N += phase()*sizeGroups_[i]/sizeGroups_[i].x();
86  }
87 
88  return tN;
89 }
90 
91 
93 Foam::diameterModels::velocityGroup::fSum() const
94 {
95  tmp<volScalarField> tsumSizeGroups
96  (
98  (
99  "sumSizeGroups",
100  phase().mesh(),
102  )
103  );
104 
105  volScalarField& sumSizeGroups = tsumSizeGroups.ref();
106 
107  forAll(sizeGroups_, i)
108  {
109  sumSizeGroups += sizeGroups_[i];
110  }
111 
112  return tsumSizeGroups;
113 }
114 
115 
116 void Foam::diameterModels::velocityGroup::scale()
117 {
118  Info<< "Scaling sizeGroups for velocityGroup " << phase().name() << endl;
119 
120  forAll(sizeGroups_, i)
121  {
122  sizeGroups_[i].max(0);
123  };
124 
125  f_ = fSum();
126 
127  forAll(sizeGroups_, i)
128  {
129  sizeGroups_[i] /= f_;
130 
131  sizeGroups_[i].correctBoundaryConditions();
132  };
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137 
139 (
140  const dictionary& diameterProperties,
141  const phaseModel& phase
142 )
143 :
144  diameterModel(diameterProperties, phase),
145  popBalName_(diameterProperties.lookup("populationBalance")),
146  f_
147  (
148  IOobject
149  (
150  IOobject::groupName
151  (
152  "f",
153  phase.name()
154  ),
155  phase.time().timeName(),
156  phase.mesh(),
157  IOobject::MUST_READ,
158  IOobject::AUTO_WRITE
159  ),
160  phase.mesh()
161  ),
162  sizeGroups_
163  (
164  diameterProperties.lookup("sizeGroups"),
165  sizeGroup::iNew(phase, *this)
166  ),
167  d_(IOobject::groupName("d", phase.name()), dsm())
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
172 
174 {}
175 
176 
177 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178 
180 {
181  return d_;
182 }
183 
184 
186 {
187  tmp<volScalarField> tA
188  (
190  (
191  "a",
192  phase().mesh(),
194  )
195  );
196 
197  volScalarField& a = tA.ref();
198 
199  forAll(sizeGroups_, i)
200  {
201  const sizeGroup& fi = sizeGroups_[i];
202 
203  a += fi.a()*fi/fi.x();
204  }
205 
206  return phase()*a;
207 }
208 
209 
211 {
212  const populationBalanceModel& popBal =
213  phase().mesh().lookupObject<populationBalanceModel>(popBalName_);
214 
215  const pimpleControl& pimple =
216  phase().mesh().lookupObject<pimpleControl>("solutionControl");
217 
218  if (!popBal.solveOnFinalIterOnly() || pimple.finalPimpleIter())
219  {
220  forAll(sizeGroups_, i)
221  {
222  sizeGroups_[i].correct();
223  }
224 
225  if
226  (
227  phase().mesh().solution().solverDict(popBalName_)
228  .lookupOrDefault<Switch>
229  (
230  "scale",
231  true
232  )
233  )
234  {
235  scale();
236  }
237 
238  f_ = fSum();
239 
240  f_.correctBoundaryConditions();
241 
242  Info<< phase().name() << " sizeGroups-sum volume fraction, min, max = "
243  << f_.weightedAverage(phase().mesh().V()).value()
244  << ' ' << min(f_).value()
245  << ' ' << max(f_).value()
246  << endl;
247 
248  d_ = dsm();
249 
250  Info<< this->phase().name() << " Sauter mean diameter, min, max = "
251  << d_.weightedAverage(d_.mesh().V()).value()
252  << ' ' << min(d_).value()
253  << ' ' << max(d_).value()
254  << endl;
255  }
256 }
257 
258 
260 (
261  const dictionary& phaseProperties
262 )
263 {
264  diameterModel::read(phaseProperties);
265 
266  return true;
267 }
268 
269 
270 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volScalarField > d() const
Get the diameter field.
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
sizeGroup(const word &name, const dictionary &dict, const phaseModel &phase, const velocityGroup &velocityGroup, const fvMesh &mesh)
const word & name() const
Definition: phaseModel.H:110
const dimensionSet dimless
fvMesh & mesh
velocityGroup(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > a() const
Get the surface area per unit volume field.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
stressControl lookup("compactNormalStress") >> compactNormalStress
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)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Internal & ref()
Return a reference to the dimensioned internal field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void max(const dimensioned< Type > &)
const dimensionSet dimVolume
messageStream Info
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
virtual void correct()
Correct the model.
const label & i() const
Return index of the size group within the population balance.