fractal.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) 2019-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 
26 #include "fractal.H"
27 #include "populationBalanceModel.H"
29 #include "fvmDdt.H"
30 #include "fvmDiv.H"
31 #include "fvmSup.H"
32 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace diameterModels
40 {
41 namespace shapeModels
42 {
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
51 
53 Foam::diameterModels::shapeModels::fractal::dColl() const
54 {
55  tmp<volScalarField> tDColl
56  (
58  (
59  "dColl",
60  group().mesh(),
62  )
63  );
64 
65  volScalarField& dColl = tDColl.ref();
66 
67  dColl =
68  6/kappa_
69  *pow(group().x()*pow3(kappa_)/(36*pi*alphaC_), 1/Df_);
70 
71  return tDColl;
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
78 (
79  const dictionary& dict,
80  const sizeGroup& group,
81  const dictionary& groupDict
82 )
83 :
85  kappa_
86  (
87  sizeGroup::fieldIo("kappa", group.i(), group.group()),
88  sizeGroup::field("kappa", group.i(), group.group())
89  ),
90  Df_("Df", dimless, groupDict),
91  alphaC_("alphaC", dimless, groupDict),
92  dColl_
93  (
94  IOobject
95  (
96  "dColl" + group.name().substr(1),
97  group.mesh().time().name(),
98  group.mesh()
99  ),
100  this->dColl()
101  ),
102  Su_
103  (
104  IOobject
105  (
106  IOobject::groupName("Su", kappa_.name()),
107  group.mesh().time().name(),
108  group.mesh()
109  ),
110  group.mesh(),
112  )
113 {
114  // Check and filter for old syntax (remove in due course)
115  if (groupDict.found("kappa"))
116  {
118  << "A 'kappa' entry should not be specified for size-group #"
119  << group.i() << " of population balance "
120  << group.group().popBalName()
121  << ". Instead, the value should be initialised within the field, "
122  << kappa_.name() << " (or the default field, "
123  << IOobject::groupName("kappaDefault", group.phase().name())
124  << ", as appropriate)."
125  << exit(FatalError);
126  }
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
131 
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
140 {
141  return kappa_;
142 }
143 
144 
147 {
148  return Su_;
149 }
150 
151 
153 {
154  const sizeGroup& fi = group();
155  const label i = fi.i();
156  const phaseModel& phase = fi.phase();
157  const volScalarField& alpha = phase;
158  const volScalarField& rho = phase.rho();
159 
160  const populationBalanceModel& popBal = fi.group().popBal();
161 
162  const volScalarField alphaFi
163  (
165  (
166  alpha.member() + fi.member().capitalise(),
167  alpha.group()
168  ),
169  alpha*fi
170  );
171 
172  const surfaceScalarField alphaFiPhi
173  (
175  (
176  alpha.member() + fi.member().capitalise() + "Phi",
177  alpha.group()
178  ),
179  max(fvc::interpolate(fi, "fi"), small)*phase.alphaPhi()
180  );
181 
182  // !!! Create pointers to other kappa fields. The shortcut taken here is
183  // that only adjacent fields are set. We know that the transfers never
184  // extend beyond these adjacent groups. Eventually the plan is to store all
185  // the fields in a pointer list anyway, rather than having nested per-group
186  // models, so then it will be possible just to pass the list directly.
187  UPtrList<const volScalarField> flds(popBal.sizeGroups().size());
188  for (label deltai = -1; deltai <= +1; ++ deltai)
189  {
190  const label j = fi.i() + deltai;
191  if (j < 0 || j >= popBal.sizeGroups().size()) continue;
192  const sizeGroup& fj = popBal.sizeGroups()[j];
193  flds.set(fj.i(), &model(fj).fld());
194  }
195 
196  fvScalarMatrix kappaEqn
197  (
198  fvm::ddt(alpha, fi, kappa_)
199  + fvm::div(alphaFiPhi, kappa_)
200  ==
201  Su_ + fvm::Sp(popBal.Sp(i)*fi, kappa_)
202  + popBal.expansionSu(i, flds) + fvm::Sp(popBal.expansionSp(i)*fi, kappa_)
203  + popBal.modelSourceSu(i, flds)
204  + popBal.fluid().fvModels().source(alphaFi, rho, kappa_)/rho
205  - correction
206  (
207  fvm::Sp
208  (
209  max(phase.residualAlpha() - alpha*fi, scalar(0))
210  /kappa_.mesh().time().deltaT(),
211  kappa_
212  )
213  )
214  );
215 
216  kappaEqn.relax();
217 
218  popBal.fluid().fvConstraints().constrain(kappaEqn);
219 
220  kappaEqn.solve();
221 
222  popBal.fluid().fvConstraints().constrain(kappa_);
223 
224  // Bound kappa so that the surface-area-volume ratio is greater than that
225  // of spherical particles of this group, but less than that of the
226  // particles represented by the first size group
227  const sizeGroup& f0 = popBal.sizeGroups().first();
228  kappa_ = min(max(kappa_, 6/fi.dSph()), 6/f0.dSph());
229 
230  kappa_.correctBoundaryConditions();
231 
232  // Update the collisional diameter
233  dColl_ = dColl();
234 }
235 
236 
239 {
240  return kappa_*group().x();
241 }
242 
243 
246 {
247  return dColl_;
248 }
249 
250 
251 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:146
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
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
bool set(const label) const
Is element set.
Definition: UPtrListI.H:87
Base class for modelling evolution of secondary representative properties of a size class....
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
tmp< volScalarField::Internal > Sp(const label i) const
Return the implicit coalescence and breakup source term.
const phaseSystem & fluid() const
Return reference to the phaseSystem.
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
tmp< volScalarField::Internal > modelSourceSu(const label i, const UPtrList< const volScalarField > &flds=UPtrList< const volScalarField >()) const
Return the explicit model source source term.
tmp< volScalarField::Internal > expansionSp(const label i) const
Return the implicit expansion source term.
tmp< volScalarField::Internal > expansionSu(const label i, const UPtrList< const volScalarField > &flds=UPtrList< const volScalarField >()) const
Return the explicit expansion source term.
Base class for modelling the shape of the particles belonging to a size class through alternative dia...
Definition: shapeModel.H:61
const sizeGroup & group() const
Return reference to size group.
Definition: shapeModel.C:94
Class for modelling the shape of particle aggregates using the concept of fractal geometry....
Definition: fractal.H:94
virtual void correct()
Correct the collisional diameter.
Definition: fractal.C:152
fractal(const dictionary &dict, const sizeGroup &group, const dictionary &groupDict)
Construct from dictionaries and sizeGroup.
Definition: fractal.C:78
virtual const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: fractal.C:238
virtual volScalarField::Internal & src()
Access to secondary property source.
Definition: fractal.C:146
virtual const volScalarField & fld() const
Return reference to secondary property field.
Definition: fractal.C:139
virtual const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: fractal.C:245
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
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
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 word & popBalName() const
Return name of populationBalance this velocityGroup belongs to.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
Foam::fvModels & fvModels(fvMesh &mesh)
Access the fvModels.
Definition: phaseSystemI.H:226
Foam::fvConstraints & fvConstraints(fvMesh &mesh)
Access the fvConstraints.
Definition: phaseSystemI.H:238
A class for managing temporary objects.
Definition: tmp.H:55
word capitalise() const
Return the word with the first letter capitalised.
Definition: wordI.H:131
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
addToRunTimeSelectionTable(shapeModel, fractal, dictionary)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
const dimensionSet dimless
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict