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 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
33 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const word& name,
40  const dictionary& dict,
41  const phaseModel& phase,
43  const fvMesh& mesh
44 )
45 :
47  (
48  IOobject
49  (
50  IOobject::groupName
51  (
52  name,
53  velocityGroup.phase().name()
54  ),
55  mesh.time().name(),
56  mesh,
57  IOobject::READ_IF_PRESENT,
58  IOobject::AUTO_WRITE
59  ),
60  mesh,
61  dimensionedScalar(name, dimless, dict.lookup<scalar>("value")),
62  velocityGroup.f().boundaryField().types()
63  ),
64  dict_(dict),
65  phase_(phase),
66  velocityGroup_(velocityGroup),
67  dSph_("dSph", dimLength, dict),
68  x_("x", pi/6*pow3(dSph_)),
69  value_(dict.lookup<scalar>("value"))
70 {
71  // Adjust refValue at mixedFvPatchField boundaries
72  forAll(this->boundaryField(), patchi)
73  {
74  typedef mixedFvPatchField<scalar> mixedFvPatchScalarField;
75 
76  if
77  (
78  isA<const mixedFvPatchScalarField>(this->boundaryField()[patchi])
79  )
80  {
81  mixedFvPatchScalarField& f =
82  refCast<mixedFvPatchScalarField>
83  (
84  this->boundaryFieldRef()[patchi]
85  );
86 
87  f.refValue() = value_;
88  }
89  }
90 
91  shapeModel_ = shapeModel::New(velocityGroup_.diameterProperties(), *this);
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
105 {
106  notImplemented("sizeGroup::clone() const");
107  return autoPtr<sizeGroup>(nullptr);
108 }
109 
110 
112 {
113  if (!i_.valid())
114  {
115  const populationBalanceModel& popBal =
116  this->mesh().lookupObject<populationBalanceModel>
117  (
118  velocityGroup_.popBalName()
119  );
120 
121  forAll(popBal.sizeGroups(), j)
122  {
123  if (&popBal.sizeGroups()[j] == &*this)
124  {
125  i_.set(new label(j));
126  }
127  }
128  }
129 
130  return i_();
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 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
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 dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
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.
static autoPtr< shapeModel > New(const dictionary &dict, const sizeGroup &group)
Definition: shapeModel.C:56
void correct()
Correct secondary properties.
Definition: sizeGroup.C:148
const label & i() const
Return index of the size group within the population balance.
Definition: sizeGroup.C:111
sizeGroup(const word &name, const dictionary &dict, const phaseModel &phase, const velocityGroup &velocityGroup, const fvMesh &mesh)
Definition: sizeGroup.C:38
virtual ~sizeGroup()
Destructor.
Definition: sizeGroup.C:97
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:104
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:86
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:343
label patchi
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)
const dimensionSet dimless
const dimensionSet dimLength
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList f(nPoints)
dictionary dict