populationBalanceSetSizeDistribution.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) 2024-2025 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 
27 #include "distribution.H"
29 #include "velocityGroup.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
39  (
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
52 (
53  const word& name,
54  const Time& runTime,
55  const dictionary& dict
56 )
57 :
58  fvMeshFunctionObject(name, runTime, dict),
59  popBalName_(dict.lookupOrDefault("populationBalance", word::null)),
60  phaseName_(dict.lookupOrDefault("phase", word::null)),
61  distribution_
62  (
63  distribution::New(dimLength, dict.subDict("distribution"), 3, -1)
64  )
65 {
66  const bool havePopBal = popBalName_ != word::null;
67  const bool havePhase = phaseName_ != word::null;
68  if (havePopBal == havePhase)
69  {
71  << (havePopBal ? "both" : "neither") << " of keywords "
72  << "populationBalance " << (havePopBal ? "and" : "or")
73  << " phase defined in dictionary " << dict.name()
74  << exit(FatalIOError);
75  }
76 
77  read(dict);
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
82 
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  return true;
93 }
94 
95 
97 {
99  {
101  << "The " << typeName << " function cannot be executed at run-time"
102  << endl;
103 
104  return false;
105  }
106 
107  Info<< nl << indent << name() << token::COLON << endl << incrIndent;
108 
109  if (popBalName_ != word::null)
110  {
112  obr_.lookupObjectRef<diameterModels::populationBalanceModel>
113  (
114  popBalName_
115  );
116 
118 
119  // Set the size-group fractions for the population balance
120  forAll(popBal.sizeGroups(), sizeGroupi)
121  {
122  diameterModels::sizeGroup& fi = popBal.sizeGroups()[sizeGroupi];
123 
124  fi == popBal.etaV(sizeGroupi, distribution_());
125 
126  Info<< indent << "Writing " << fi.name() << endl;
127 
128  fi.write();
129 
130  velocityGroupPtrs.insert(fi.phase().name(), &fi.group());
131  }
132 
133  // Set the volume fractions if there are multiple velocity groups
134  if (velocityGroupPtrs.size() > 1)
135  {
136  tmp<volScalarField> talphas =
138  (
139  IOobject::groupName("alpha", popBal.name()),
140  mesh(),
141  dimensionedScalar(dimless, scalar(0))
142  );
143  volScalarField& alphas = talphas.ref();
144 
146  (
148  velocityGroupPtrs,
149  iter
150  )
151  {
152  alphas += iter()->phase();
153  }
154 
156  (
158  velocityGroupPtrs,
159  iter
160  )
161  {
163  const_cast<phaseModel&>(iter()->phase());
164 
165  const label i0 = iter()->sizeGroups().first().i();
166  const label i1 = iter()->sizeGroups().last().i();
167 
168  alpha == popBal.etaV(labelPair(i0, i1), distribution_)*alphas;
169 
170  Info<< indent << "Writing " << alpha.name() << endl;
171 
172  alpha.write();
173  }
174  }
175  }
176  else
177  {
179  refCast<diameterModels::velocityGroup>
180  (
181  obr_.lookupObjectRef<phaseSystem>
182  (
184  ).phases()[phaseName_].diameter()
185  );
186 
187  // Set the size-group fractions for this velocity group
188  forAll(velGrp.sizeGroups(), i)
189  {
190  diameterModels::sizeGroup& fi = velGrp.sizeGroups()[i];
191 
192  const label sizeGroupi = velGrp.sizeGroups().first().i() + i;
193 
194  fi == velGrp.popBal().etaV(sizeGroupi, distribution_());
195 
196  Info<< indent << "Writing " << fi.name() << endl;
197 
198  fi.write();
199  }
200  }
201 
202  Info<< decrIndent;
203 
204  return true;
205 }
206 
207 
208 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
An STL-conforming hash table.
Definition: HashTable.H:127
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
dimensionedScalar etaV(const label i, const dimensionedScalar &v) const
Return the volume allocation coefficient for a single volume.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const velocityGroup & group() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:44
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
const populationBalanceModel & popBal() const
Return the populationBalance this velocityGroup belongs to.
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
Abstract base-class for Time/database functionObjects.
static bool postProcess
Global post-processing mode switch.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool read(const dictionary &)
Read optional controls.
Sets the population balance size distribution by overwriting the values in the size-group fraction fi...
populationBalanceSetSizeDistribution(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool write()
Calculate and write the size group fraction fields.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
Class to represent a system of phases.
Definition: phaseSystem.H:74
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:252
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict