sizeGroup.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-2023 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 "sizeGroup.H"
27 #include "mixedFvPatchField.H"
28 #include "shapeModel.H"
29 
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
35 (
36  const word& name,
37  const label i,
38  const velocityGroup& group,
39  const IOobject::readOption r,
40  const bool registerObject
41 )
42 {
43  return
44  IOobject
45  (
47  (
48  name + (i == -1 ? "Default" : Foam::name(i)),
49  group.phase().name()
50  ),
51  group.phase().mesh().time().name(),
52  group.phase().mesh(),
53  r,
56  );
57 }
58 
59 
61 (
62  const word& name,
63  const label i,
64  const velocityGroup& group
65 )
66 {
68  (
69  fieldIo(name, i, group, IOobject::MUST_READ, false)
70  );
71 
72  return
74  (
75  new volScalarField
76  (
77  io.headerOk()
78  ? io
79  : fieldIo(name, -1, group, IOobject::MUST_READ, false),
80  group.phase().mesh()
81  )
82  );
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const label i,
91  const dictionary& dict,
92  const velocityGroup& group
93 )
94 :
95  volScalarField(fieldIo("f", i, group), field("f", i, group)),
96  i_(i),
97  group_(group),
98  dSph_("dSph", dimLength, dict),
99  x_("x", pi/6*pow3(dSph_)),
100  shapeModel_(shapeModel::New(group_.diameterProperties(), *this, dict))
101 {
102  // Check and filter for old syntax (remove in due course)
103  if (dict.found("value"))
104  {
106  << "A 'value' entry should not be specified for size-group #"
107  << i << " of population balance "
108  << group.popBalName()
109  << ". Instead, the value should be initialised within the field, "
110  << this->name() << " (or the default field, "
111  << IOobject::groupName("fDefault", group.phase().name())
112  << ", as appropriate)."
113  << exit(FatalError);
114  }
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
119 
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
128 {
130  return autoPtr<sizeGroup>(nullptr);
131 }
132 
133 
136 {
137  return shapeModel_->a();
138 }
139 
140 
143 {
144  return shapeModel_->d();
145 }
146 
147 
149 {
150  shapeModel_->correct();
151 }
152 
153 
154 // ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
bool & registerObject()
Register object created from this IOobject with registry if true.
Definition: IOobject.H:346
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:167
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const phaseModel & phase() const
Return the phase.
Base class for modelling the shape of the particles belonging to a size class through alternative dia...
Definition: shapeModel.H:61
void correct()
Correct secondary properties.
Definition: sizeGroup.C:148
static tmp< volScalarField > field(const word &name, const label i, const velocityGroup &group)
Return IO a size-group field. Construction helper.
Definition: sizeGroup.C:61
const velocityGroup & group() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:44
static IOobject fieldIo(const word &name, const label i, const velocityGroup &group, const IOobject::readOption r=IOobject::NO_READ, const bool registerObject=true)
Return IO for a size-group field. Construction helper.
Definition: sizeGroup.C:35
virtual ~sizeGroup()
Destructor.
Definition: sizeGroup.C:120
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: sizeGroup.C:135
autoPtr< sizeGroup > clone() const
Return clone.
Definition: sizeGroup.C:127
sizeGroup(const label i, const dictionary &dict, const velocityGroup &group)
Construct from index, dictionary and velocity group.
Definition: sizeGroup.C:89
label i() const
Return index of the size group within the population balance.
Definition: sizeGroupI.H:30
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroup.C:142
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
const word & popBalName() const
Return name of populationBalance this velocityGroup belongs to.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const char *const group
Group name for atomic constants.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
error FatalError
dictionary dict