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-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 "fractal.H"
27 #include "sinteringModel.H"
28 #include "fvmDdt.H"
29 #include "fvmDiv.H"
30 #include "fvcDiv.H"
31 #include "fvmSup.H"
32 #include "fvcSup.H"
33 #include "fvcDdt.H"
34 #include "mixedFvPatchField.H"
35 #include "mathematicalConstants.H"
36 #include "populationBalanceModel.H"
38 #include "phaseSystem.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace diameterModels
45 {
46 namespace shapeModels
47 {
48  defineTypeNameAndDebug(fractal, 0);
49  addToRunTimeSelectionTable(shapeModel, fractal, dictionary);
50 }
51 }
52 }
53 
54 template<>
55 const char*
57 <
59  3
60 >::names[] = {"hardSphere", "ParkRogak", "conserved"};
61 
62 const Foam::NamedEnum
63 <
65  3
67 
68 
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
75 (
76  const dictionary& dict,
77  const sizeGroup& group
78 )
79 :
80  SecondaryPropertyModel<shapeModel>(dict, group),
81  kappa_
82  (
83  IOobject
84  (
85  "kappa" + group.name().substr(1),
86  group.mesh().time().timeName(),
87  group.mesh(),
88  IOobject::READ_IF_PRESENT,
89  IOobject::AUTO_WRITE
90  ),
91  group.mesh(),
93  (
94  "kappa",
95  inv(dimLength),
96  group.dict()
97  ),
98  group.VelocityGroup().f().boundaryField().types()
99  ),
100  Df_("Df", dimless, group.dict()),
101  alphaC_("alphaC", dimless, group.dict()),
102  dColl_
103  (
104  IOobject
105  (
106  "dColl" + group.name().substr(1),
107  group.mesh().time().timeName(),
108  group.mesh()
109  ),
110  this->dColl()
111  ),
112  Su_
113  (
114  IOobject
115  (
116  IOobject::groupName("Su", kappa_.name()),
117  group.mesh().time().timeName(),
118  group.mesh()
119  ),
120  group.mesh(),
121  dimensionedScalar(kappa_.dimensions()/dimTime, Zero)
122  )
123 {
124  // Adjust refValue at mixedFvPatchField boundaries
125  forAll(kappa_.boundaryField(), patchi)
126  {
127  typedef mixedFvPatchField<scalar> mixedFvPatchScalarField;
128 
129  if
130  (
131  isA<mixedFvPatchScalarField>(kappa_.boundaryFieldRef()[patchi])
132  )
133  {
134  mixedFvPatchScalarField& kappa =
135  refCast<mixedFvPatchScalarField>
136  (
137  kappa_.boundaryFieldRef()[patchi]
138  );
139 
140  kappa.refValue() = sizeGroup_.dict().lookup<scalar>("kappa");
141  }
142  }
143 
144  sinteringModel_ =
145  sinteringModel::New(dict.subDict(type() + "Coeffs"), *this);
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
159 {
160  return kappa_;
161 }
162 
163 
166 {
167  return Su_;
168 }
169 
170 
172 Foam::diameterModels::shapeModels::fractal::dColl() const
173 {
174  tmp<volScalarField> tDColl
175  (
177  (
178  "dColl",
179  sizeGroup_.mesh(),
181  )
182  );
183 
184  volScalarField& dColl = tDColl.ref();
185 
186  dColl =
187  6.0/kappa_
188  *pow(sizeGroup_.x()*pow3(kappa_)/(36.0*pi*alphaC_), 1.0/Df_);
189 
190  return tDColl;
191 }
192 
193 
195 {
196  const sizeGroup& fi = sizeGroup_;
197  const phaseModel& phase = fi.phase();
198  const volScalarField& alpha = phase;
199  const volScalarField& rho = phase.thermo().rho();
200 
201  const populationBalanceModel& popBal =
202  sizeGroup_.mesh().lookupObject<populationBalanceModel>
203  (
204  sizeGroup_.VelocityGroup().popBalName()
205  );
206 
207  surfaceScalarField fAlphaRhoPhi
208  (
209  "fAlphaRhoPhi",
210  max(fvc::interpolate(fi, "fi"), SMALL)*phase.alphaRhoPhi()
211  );
212 
213  fvScalarMatrix kappaEqn
214  (
215  fvc::ddt(alpha, rho, fi)*kappa_.oldTime()
216  + alpha*rho*fi*fvm::ddt(kappa_)
217  + fvm::div(fAlphaRhoPhi, kappa_)
218  + fvm::SuSp
219  (
220  fi
221  *(
222  fi.VelocityGroup().dmdt()
223  - (fvc::ddt(alpha, rho) + fvc::div(phase.alphaRhoPhi()))
224  ),
225  kappa_
226  )
227  ==
228  - sinteringModel_->R()
229  + fvc::Su(Su_*rho, kappa_)
230  - fvm::SuSp(popBal.SuSp(fi.i()())*fi*rho, kappa_)
231  + fvc::ddt(fi.phase().residualAlpha()*rho, kappa_)
232  - fvm::ddt(fi.phase().residualAlpha()*rho, kappa_)
233  );
234 
235  kappaEqn.relax();
236 
237  kappaEqn.solve();
238 
239  // Bounding of kappa assuming first sizeGroup to represent one primary
240  // particle
241  kappa_ =
242  min
243  (
244  max(kappa_, 6.0/sizeGroup_.dSph()),
245  6.0/popBal.sizeGroups().first().dSph()
246  );
247 
248  kappa_.correctBoundaryConditions();
249 
250  // Update the collisional diameter
251  dColl_ = dColl();
252 }
253 
254 
257 {
258  return kappa_*sizeGroup_.x();
259 }
260 
261 
263 (
264  const volScalarField &Su,
265  const sizeGroup &fu,
266  const driftModel &model
267 )
268 {
269  surfaceGrowthTypes sgType
270  (
271  sgTypeNames_.read
272  (
273  model.dict().lookup("surfaceGrowthType")
274  )
275  );
276 
277  const volScalarField& sourceKappa =
278  SecondaryPropertyModelTable()[SecondaryPropertyName(fu)]->fld();
279 
280  switch (sgType)
281  {
282  case sgHardSphere:
283  {
284  Su_ += sourceKappa*fu.dSph()/sizeGroup_.dSph()*Su;
285 
286  break;
287  }
288 
289  case sgParkRogak:
290  {
291  const fractal& sourceShape =
292  refCast<const fractal>
293  (
294  fu.shapeModelPtr()()
295  );
296 
297  volScalarField dp(6.0/sourceKappa);
298  const volScalarField a(sourceKappa*fu.x());
299  const dimensionedScalar dv(sizeGroup_.x() - fu.x());
300 
301  const volScalarField da1
302  (
303  (2.0/3.0)*dv
304  *(
305  sourceKappa
306  + sourceShape.Df_*(1.0/sourceShape.d() - 1.0/dp)
307  )
308  );
309 
310  dp += 6.0*(dv*a - fu.x()*da1)/sqr(a);
311 
312  const volScalarField np(6.0*sizeGroup_.x()/pi/pow3(dp));
313  const volScalarField dc(dp*pow(np/alphaC_, 1.0/Df_));
314 
315  const volScalarField da2
316  (
317  dv*(4.0/dp + 2.0*Df_/3.0*(1.0/dc - 1.0/dp))
318  );
319 
320  Su_ += (a + 0.5*da1 + 0.5*da2)/sizeGroup_.x()*Su;
321 
322  break;
323  }
324 
325  case sgConserved:
326  {
327  SecondaryPropertyModel::addDrift(Su, fu, model);
328 
329  break;
330  }
331 
332  default:
333  {
335  << "Unknown surface growth type. Valid types are:"
336  << sgTypeNames_ << nl << exit(FatalError);
337  }
338  }
339 }
340 
341 
342 // ************************************************************************* //
surfaceGrowthTypes
Surface growth type enumeration.
Definition: fractal.H:111
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< GeometricField< Type, fvPatchField, volMesh > > SuSp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:102
zeroField Su
Definition: alphaSuSp.H:1
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
fractal(const dictionary &dict, const sizeGroup &group)
Construct from dictionary and sizeGroup.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the collisional diameter.
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
static const scalar SMALL
Definition: scalar.H:115
virtual const volScalarField & fld() const
Return reference to secondary property field.
static const NamedEnum< surfaceGrowthTypes, 3 > sgTypeNames_
Surface growth type names.
Definition: fractal.H:119
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Calculate the first temporal derivative.
dynamicFvMesh & mesh
virtual const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
virtual volScalarField & src()
Access to secondary property source.
Calculate the matrix for the first temporal derivative.
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
Calculate the field for explicit evaluation of implicit and explicit sources.
Calculate the divergence of the given field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Calculate the matrix for the divergence of the given field and flux.
dimensionedScalar pow3(const dimensionedScalar &ds)
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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void addDrift(const volScalarField &Su, const sizeGroup &fu, const driftModel &model)
Add drift contribution to secondary property source.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.