All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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 "populationBalanceModel.H"
28 #include "mixedFvPatchField.H"
29 #include "shapeModel.H"
30 #include "mathematicalConstants.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const word& name,
42  const dictionary& dict,
43  const phaseModel& phase,
44  const velocityGroup& velocityGroup,
45  const fvMesh& mesh
46 )
47 :
49  (
50  IOobject
51  (
52  IOobject::groupName
53  (
54  name,
55  IOobject::groupName
56  (
57  velocityGroup.phase().name(),
58  velocityGroup.popBalName()
59  )
60  ),
61  mesh.time().timeName(),
62  mesh,
63  IOobject::READ_IF_PRESENT,
64  IOobject::AUTO_WRITE
65  ),
66  mesh,
67  dimensionedScalar(name, dimless, dict.lookup<scalar>("value")),
68  velocityGroup.f().boundaryField().types()
69  ),
70  dict_(dict),
71  phase_(phase),
72  velocityGroup_(velocityGroup),
73  dSph_("dSph", dimLength, dict),
74  x_("x", pi/6.0*pow3(dSph_)),
75  value_(dict.lookup<scalar>("value"))
76 {
77  // Adjust refValue at mixedFvPatchField boundaries
78  forAll(this->boundaryField(), patchi)
79  {
80  typedef mixedFvPatchField<scalar> mixedFvPatchScalarField;
81 
82  if
83  (
84  isA<mixedFvPatchScalarField>(this->boundaryFieldRef()[patchi])
85  )
86  {
87  mixedFvPatchScalarField& f =
88  refCast<mixedFvPatchScalarField>
89  (
90  this->boundaryFieldRef()[patchi]
91  );
92 
93  f.refValue() = value_;
94  }
95  }
96 
97  shapeModel_ = shapeModel::New(velocityGroup_.diameterProperties(), *this);
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102 
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
111 {
112  notImplemented("sizeGroup::clone() const");
113  return autoPtr<sizeGroup>(nullptr);
114 }
115 
116 
118 {
119  if (!i_.valid())
120  {
121  const populationBalanceModel& popBal =
122  this->mesh().lookupObject<populationBalanceModel>
123  (
124  velocityGroup_.popBalName()
125  );
126 
127  forAll(popBal.sizeGroups(), j)
128  {
129  if (&popBal.sizeGroups()[j] == &*this)
130  {
131  i_.set(new label(j));
132  }
133  }
134  }
135 
136  return i_;
137 }
138 
139 
142 {
143  return shapeModel_->a();
144 }
145 
146 
149 {
150  return shapeModel_->d();
151 }
152 
153 
155 {
156  shapeModel_->correct();
157 }
158 
159 
160 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
sizeGroup(const word &name, const dictionary &dict, const phaseModel &phase, const velocityGroup &velocityGroup, const fvMesh &mesh)
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
void correct()
Correct secondary properties.
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
virtual ~sizeGroup()
Destructor.
autoPtr< sizeGroup > clone() const
Return clone.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
labelList f(nPoints)
static autoPtr< shapeModel > New(const dictionary &dict, const sizeGroup &group)
dimensionedScalar pow3(const dimensionedScalar &ds)
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:356
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const autoPtr< label > & i() const
Return label of the sizeGroup within the population balance.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53