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-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 "fractal.H"
28 #include "sinteringModel.H"
29 #include "fvmDdt.H"
30 #include "fvmDiv.H"
31 #include "fvmSup.H"
32 #include "mixedFvPatchField.H"
33 
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace diameterModels
41 {
42 namespace shapeModels
43 {
46 }
47 }
48 }
49 
50 template<>
51 const char*
53 <
55  3
56 >::names[] = {"hardSphere", "ParkRogak", "conserved"};
57 
58 const Foam::NamedEnum
59 <
61  3
63 
64 
65 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
66 
68 Foam::diameterModels::shapeModels::fractal::dColl() const
69 {
70  tmp<volScalarField> tDColl
71  (
73  (
74  "dColl",
75  group().mesh(),
77  )
78  );
79 
80  volScalarField& dColl = tDColl.ref();
81 
82  dColl =
83  6/kappa_
84  *pow(group().x()*pow3(kappa_)/(36*pi*alphaC_), 1/Df_);
85 
86  return tDColl;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
93 (
94  const dictionary& dict,
95  const sizeGroup& group,
96  const dictionary& groupDict
97 )
98 :
100  kappa_
101  (
102  sizeGroup::fieldIo("kappa", group.i(), group.group()),
103  sizeGroup::field("kappa", group.i(), group.group())
104  ),
105  Df_("Df", dimless, groupDict),
106  alphaC_("alphaC", dimless, groupDict),
107  dColl_
108  (
109  IOobject
110  (
111  "dColl" + group.name().substr(1),
112  group.mesh().time().name(),
113  group.mesh()
114  ),
115  this->dColl()
116  ),
117  Su_
118  (
119  IOobject
120  (
121  IOobject::groupName("Su", kappa_.name()),
122  group.mesh().time().name(),
123  group.mesh()
124  ),
125  group.mesh(),
127  ),
128  sinteringModel_(sinteringModel::New(dict.subDict(type() + "Coeffs"), *this))
129 {
130  // Check and filter for old syntax (remove in due course)
131  if (groupDict.found("kappa"))
132  {
134  << "A 'kappa' entry should not be specified for size-group #"
135  << group.i() << " of population balance "
136  << group.group().popBalName()
137  << ". Instead, the value should be initialised within the field, "
138  << this->name() << " (or the default field, "
139  << IOobject::groupName("kappaDefault", group.phase().name())
140  << ", as appropriate)."
141  << exit(FatalError);
142  }
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
147 
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
156 {
157  return kappa_;
158 }
159 
160 
163 {
164  return Su_;
165 }
166 
167 
169 {
170  const sizeGroup& fi = group();
171  const phaseModel& phase = fi.phase();
172  const volScalarField& alpha = phase;
173 
174  const populationBalanceModel& popBal =
175  group().mesh().lookupObject<populationBalanceModel>
176  (
177  group().group().popBalName()
178  );
179 
180  surfaceScalarField fAlphaPhi
181  (
182  "fAlphaPhi",
183  max(fvc::interpolate(fi, "fi"), small)*phase.alphaPhi()
184  );
185 
186  fvScalarMatrix kappaEqn
187  (
188  fvm::ddt(alpha, fi, kappa_)
189  + fvm::div(fAlphaPhi, kappa_)
190  ==
191  - sinteringModel_->R()
192  + Su_
193  - fvm::Sp(popBal.Sp(fi.i())*fi, kappa_)
194  - correction
195  (
196  fvm::Sp
197  (
198  max(phase.residualAlpha() - alpha*fi, scalar(0))
199  /group().mesh().time().deltaT(),
200  kappa_
201  )
202  )
203  );
204 
205  kappaEqn.relax();
206 
207  kappaEqn.solve();
208 
209  // Bounding of kappa assuming first sizeGroup to represent one primary
210  // particle
211  kappa_ =
212  min
213  (
214  max(kappa_, 6/group().dSph()),
215  6/popBal.sizeGroups().first().dSph()
216  );
217 
218  kappa_.correctBoundaryConditions();
219 
220  // Update the collisional diameter
221  dColl_ = dColl();
222 }
223 
224 
227 {
228  return kappa_*group().x();
229 }
230 
231 
233 (
234  const volScalarField &Su,
235  const sizeGroup &fu,
236  const driftModel &model
237 )
238 {
239  surfaceGrowthTypes sgType
240  (
241  sgTypeNames_.read
242  (
243  model.dict().lookup("surfaceGrowthType")
244  )
245  );
246 
247  const volScalarField& sourceKappa =
248  SecondaryPropertyModelTable()[SecondaryPropertyName(fu)]->fld();
249 
250  switch (sgType)
251  {
252  case sgHardSphere:
253  {
254  Su_ += sourceKappa*fu.dSph()/group().dSph()*Su;
255 
256  break;
257  }
258 
259  case sgParkRogak:
260  {
261  const fractal& sourceShape =
262  refCast<const fractal>
263  (
264  fu.shapeModelPtr()()
265  );
266 
267  volScalarField dp(6/sourceKappa);
268  const volScalarField a(sourceKappa*fu.x());
269  const dimensionedScalar dv(group().x() - fu.x());
270 
271  const volScalarField da1
272  (
273  (2.0/3.0)*dv
274  *(
275  sourceKappa
276  + sourceShape.Df_*(1/sourceShape.d() - 1/dp)
277  )
278  );
279 
280  dp += 6*(dv*a - fu.x()*da1)/sqr(a);
281 
282  const volScalarField np(6*group().x()/pi/pow3(dp));
283  const volScalarField dc(dp*pow(np/alphaC_, 1/Df_));
284 
285  const volScalarField da2
286  (
287  dv*(4/dp + 2*Df_/3*(1/dc - 1/dp))
288  );
289 
290  Su_ += (a + 0.5*da1 + 0.5*da2)/group().x()*Su;
291 
292  break;
293  }
294 
295  case sgConserved:
296  {
298 
299  break;
300  }
301 
302  default:
303  {
305  << "Unknown surface growth type. Valid types are:"
306  << sgTypeNames_ << nl << exit(FatalError);
307  }
308  }
309 }
310 
311 
312 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
Base class for modeling evolution of secondary representative properties of a size class....
const sizeGroup & group() const
Access the sizeGroup.
virtual void addDrift(const volScalarField &Su, const sizeGroup &fu, const driftModel &model)
Add drift contribution to secondary property source.
Base class for drift models.
Definition: driftModel.H:54
const dictionary & dict() const
Return reference to model dictionary.
Definition: driftModel.H:148
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.
const volScalarField & Sp(const label i) const
Return implicit source terms.
Base class for modelling the shape of the particles belonging to a size class through alternative dia...
Definition: shapeModel.H:61
Class for modelling the shape of particle aggregates using the concept of fractal geometry....
Definition: fractal.H:106
virtual void correct()
Correct the collisional diameter.
Definition: fractal.C:168
fractal(const dictionary &dict, const sizeGroup &group, const dictionary &groupDict)
Construct from dictionaries and sizeGroup.
Definition: fractal.C:93
surfaceGrowthTypes
Surface growth type enumeration.
Definition: fractal.H:113
virtual const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: fractal.C:226
virtual const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: fractal.H:190
virtual void addDrift(const volScalarField &Su, const sizeGroup &fu, const driftModel &model)
Add drift contribution to secondary property source.
Definition: fractal.C:233
virtual const volScalarField & fld() const
Return reference to secondary property field.
Definition: fractal.C:155
static const NamedEnum< surfaceGrowthTypes, 3 > sgTypeNames_
Surface growth type names.
Definition: fractal.H:120
virtual volScalarField & src()
Access to secondary property source.
Definition: fractal.C:162
Abstract base class for modelling sintering of primary particles in fractal aggregates.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
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
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
autoPtr< shapeModel > & shapeModelPtr()
Return reference to diameterModel of the phase.
Definition: sizeGroupI.H:65
label i() const
Return index of the size group within the population balance.
Definition: sizeGroupI.H:30
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
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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
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 tmp< surfaceScalarField > alphaPhi() const =0
Return the volumetric flux of the phase.
A class for managing temporary objects.
Definition: tmp.H:55
#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.
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< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
error FatalError
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict