growthSizeGroupFvScalarFieldSource.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) 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 "populationBalanceModel.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
32 Foam::growthSizeGroupFvScalarFieldSource::w
33 (
34  const fvSource& model,
35  const diameterModels::sizeGroup& fi
36 ) const
37 {
38  // Get the moment
39  const label q = this->q();
40 
41  // Name of the weight normalisation field
42  const word wName = fi.group().phase().name() + ":" + model.name() + ":w";
43 
44  // Quick return for volume moments that do not need to compute a weight
45  if (q == 3)
46  {
47  return DimensionedField<scalar, volMesh>::New(wName, fi.mesh(), fi.x());
48  }
49 
50  // Create the weight normalisation field if it does not yet exist
51  const bool haveW = db().foundObject<volScalarField::Internal>(wName);
52  if (!haveW)
53  {
56  (
57  IOobject
58  (
59  wName,
60  internalField().mesh().time().name(),
61  internalField().mesh()
62  ),
63  internalField().mesh(),
64  pow(dimVolume, scalar(q)/3 - 1)
65  );
66 
67  wPtr->store();
68  }
69 
70  // Update the weight normalisation field if it is out of date
72  db().lookupObjectRef<volScalarField::Internal>(wName);
73  if (!haveW || !w.hasStoredOldTimes())
74  {
75  const UPtrList<diameterModels::sizeGroup>& velGrpFis =
76  fi.group().sizeGroups();
77 
78  w.primitiveFieldRef() = scalar(0);
79 
80  forAll(velGrpFis, i)
81  {
82  w.primitiveFieldRef() +=
83  velGrpFis[i]*pow(velGrpFis[i].x(), scalar(q)/3 - 1);
84  }
85  }
86 
87  // Return the normalised weight for this size-group
88  return pow(fi.x(), scalar(q)/3)/w;
89 }
90 
91 
92 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
93 
96 (
97  const fvSource& model,
99 ) const
100 {
101  const diameterModels::sizeGroup& fi =
102  refCast<const diameterModels::sizeGroup>(this->internalField());
103  const UPtrList<diameterModels::sizeGroup>& popBalFis =
104  fi.group().popBal().sizeGroups();
105 
106  const dimensionedScalar& xi = fi.x();
107  const DimensionedField<scalar, volMesh> wi(w(model, fi));
108 
110 
111  if (fi.i() == 0)
112  {
113  tinternalCoeff = neg(source)*wi/xi;
114  }
115  else
116  {
117  const dimensionedScalar& xiMinus1 = popBalFis[fi.i() - 1].x();
118  tinternalCoeff = neg(source)*wi/(xi - xiMinus1);
119  }
120 
121  if (fi.i() != popBalFis.size() - 1)
122  {
123  const dimensionedScalar& xiPlus1 = popBalFis[fi.i() + 1].x();
124  tinternalCoeff.ref() -= pos(source)*wi/(xiPlus1 - xi);
125  }
126  else
127  {
128  tinternalCoeff.ref() += pos(source)*wi/xi;
129  }
130 
131  return tinternalCoeff;
132 }
133 
134 
137 (
138  const fvSource& model
139 ) const
140 {
141  const diameterModels::sizeGroup& fi =
142  refCast<const diameterModels::sizeGroup>(this->internalField());
143  const UPtrList<diameterModels::sizeGroup>& popBalFis =
144  fi.group().popBal().sizeGroups();
145 
146  const dimensionedScalar& xi = fi.x();
147 
149 
150  if (fi.i() != 0)
151  {
152  const diameterModels::sizeGroup& fiMinus1 = popBalFis[fi.i() - 1];
153  const dimensionedScalar& xiMinus1 = fiMinus1.x();
154  tsourceCoeffs.first() =
155  fiMinus1()*w(model, fiMinus1)*(xi/xiMinus1)/(xi - xiMinus1);
156  }
157 
158  if (fi.i() != popBalFis.size() - 1)
159  {
160  const diameterModels::sizeGroup& fiPlus1 = popBalFis[fi.i() + 1];
161  const dimensionedScalar& xiPlus1 = fiPlus1.x();
162  tsourceCoeffs.second() =
163  - fiPlus1()*w(model, fiPlus1)*(xi/xiPlus1)/(xiPlus1 - xi);
164  }
165 
166  return tsourceCoeffs;
167 }
168 
169 
172 (
173  const fvSource& model,
175 ) const
176 {
177  const diameterModels::sizeGroup& fi =
178  refCast<const diameterModels::sizeGroup>(this->internalField());
179  const UPtrList<diameterModels::sizeGroup>& velGrpFis =
180  fi.group().sizeGroups();
181 
183  sourceCoeffs(model);
184 
185  return
186  fi.i() == velGrpFis.first().i()
187  ? neg(source)*tsourceCoeffs.second()
188  : fi.i() == velGrpFis.last().i()
189  ? pos(source)*tsourceCoeffs.first()
190  : pos(source)*tsourceCoeffs.first()
191  + neg(source)*tsourceCoeffs.second();
192 }
193 
194 
195 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
const Type & second() const
Return second.
Definition: PairI.H:121
const Type & first() const
Return first.
Definition: PairI.H:107
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
T & first()
Return reference to the first element of the list.
Definition: UPtrListI.H:43
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
T & last()
Return reference to the last element of the list.
Definition: UPtrListI.H:57
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const velocityGroup & group() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:44
label i() const
Return index of the size group within the population balance.
Definition: sizeGroupI.H:30
const populationBalanceModel & popBal() const
Return the populationBalance this velocityGroup belongs to.
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
Base class for finite volume sources.
Definition: fvSource.H:52
virtual tmp< DimensionedField< scalar, volMesh > > sourceCoeff(const fvSource &model, const DimensionedField< scalar, volMesh > &source) const
Return the combined source coefficient.
virtual tmp< DimensionedField< scalar, volMesh > > internalCoeff(const fvSource &model, const DimensionedField< scalar, volMesh > &source) const
Return the internal coefficient.
virtual Pair< tmp< DimensionedField< scalar, volMesh > > > sourceCoeffs(const fvSource &model) const
Return the source coefficients for exchange with the groups below.
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
dimensionedScalar pos(const dimensionedScalar &ds)
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
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVolume
dimensionedScalar neg(const dimensionedScalar &ds)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.