All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2022 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 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace diameterModels
39 {
40 namespace shapeModels
41 {
42  defineTypeNameAndDebug(fractal, 0);
43  addToRunTimeSelectionTable(shapeModel, fractal, dictionary);
44 }
45 }
46 }
47 
48 template<>
49 const char*
51 <
53  3
54 >::names[] = {"hardSphere", "ParkRogak", "conserved"};
55 
56 const Foam::NamedEnum
57 <
59  3
61 
62 
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const dictionary& dict,
71  const sizeGroup& group
72 )
73 :
74  SecondaryPropertyModel<shapeModel>(dict, group),
75  kappa_
76  (
77  IOobject
78  (
79  "kappa" + group.name().substr(1),
80  group.mesh().time().timeName(),
81  group.mesh(),
82  IOobject::READ_IF_PRESENT,
83  IOobject::AUTO_WRITE
84  ),
85  group.mesh(),
87  (
88  "kappa",
89  inv(dimLength),
90  group.dict()
91  ),
92  group.VelocityGroup().f().boundaryField().types()
93  ),
94  Df_("Df", dimless, group.dict()),
95  alphaC_("alphaC", dimless, group.dict()),
96  dColl_
97  (
98  IOobject
99  (
100  "dColl" + group.name().substr(1),
101  group.mesh().time().timeName(),
102  group.mesh()
103  ),
104  this->dColl()
105  ),
106  Su_
107  (
108  IOobject
109  (
110  IOobject::groupName("Su", kappa_.name()),
111  group.mesh().time().timeName(),
112  group.mesh()
113  ),
114  group.mesh(),
115  dimensionedScalar(kappa_.dimensions()/dimTime, Zero)
116  )
117 {
118  // Adjust refValue at mixedFvPatchField boundaries
119  forAll(kappa_.boundaryField(), patchi)
120  {
121  typedef mixedFvPatchField<scalar> mixedFvPatchScalarField;
122 
123  if
124  (
125  isA<const mixedFvPatchScalarField>(kappa_.boundaryField()[patchi])
126  )
127  {
128  mixedFvPatchScalarField& kappa =
129  refCast<mixedFvPatchScalarField>
130  (
131  kappa_.boundaryFieldRef()[patchi]
132  );
133 
134  kappa.refValue() = sizeGroup_.dict().lookup<scalar>("kappa");
135  }
136  }
137 
138  sinteringModel_ =
139  sinteringModel::New(dict.subDict(type() + "Coeffs"), *this);
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
153 {
154  return kappa_;
155 }
156 
157 
160 {
161  return Su_;
162 }
163 
164 
166 Foam::diameterModels::shapeModels::fractal::dColl() const
167 {
168  tmp<volScalarField> tDColl
169  (
171  (
172  "dColl",
173  sizeGroup_.mesh(),
175  )
176  );
177 
178  volScalarField& dColl = tDColl.ref();
179 
180  dColl =
181  6/kappa_
182  *pow(sizeGroup_.x()*pow3(kappa_)/(36*pi*alphaC_), 1/Df_);
183 
184  return tDColl;
185 }
186 
187 
189 {
190  const sizeGroup& fi = sizeGroup_;
191  const phaseModel& phase = fi.phase();
192  const volScalarField& alpha = phase;
193 
194  const populationBalanceModel& popBal =
195  sizeGroup_.mesh().lookupObject<populationBalanceModel>
196  (
197  sizeGroup_.VelocityGroup().popBalName()
198  );
199 
200  surfaceScalarField fAlphaPhi
201  (
202  "fAlphaPhi",
203  max(fvc::interpolate(fi, "fi"), small)*phase.alphaPhi()
204  );
205 
206  fvScalarMatrix kappaEqn
207  (
208  fvm::ddt(alpha, fi, kappa_)
209  + fvm::div(fAlphaPhi, kappa_)
210  ==
211  - sinteringModel_->R()
212  + Su_
213  - fvm::Sp(popBal.Sp(fi.i())*fi, kappa_)
214  - correction
215  (
216  fvm::Sp
217  (
218  max(phase.residualAlpha() - alpha*fi, scalar(0))
219  /sizeGroup_.mesh().time().deltaT(),
220  kappa_
221  )
222  )
223  );
224 
225  kappaEqn.relax();
226 
227  kappaEqn.solve();
228 
229  // Bounding of kappa assuming first sizeGroup to represent one primary
230  // particle
231  kappa_ =
232  min
233  (
234  max(kappa_, 6/sizeGroup_.dSph()),
235  6/popBal.sizeGroups().first().dSph()
236  );
237 
238  kappa_.correctBoundaryConditions();
239 
240  // Update the collisional diameter
241  dColl_ = dColl();
242 }
243 
244 
247 {
248  return kappa_*sizeGroup_.x();
249 }
250 
251 
253 (
254  const volScalarField &Su,
255  const sizeGroup &fu,
256  const driftModel &model
257 )
258 {
259  surfaceGrowthTypes sgType
260  (
261  sgTypeNames_.read
262  (
263  model.dict().lookup("surfaceGrowthType")
264  )
265  );
266 
267  const volScalarField& sourceKappa =
268  SecondaryPropertyModelTable()[SecondaryPropertyName(fu)]->fld();
269 
270  switch (sgType)
271  {
272  case sgHardSphere:
273  {
274  Su_ += sourceKappa*fu.dSph()/sizeGroup_.dSph()*Su;
275 
276  break;
277  }
278 
279  case sgParkRogak:
280  {
281  const fractal& sourceShape =
282  refCast<const fractal>
283  (
284  fu.shapeModelPtr()()
285  );
286 
287  volScalarField dp(6/sourceKappa);
288  const volScalarField a(sourceKappa*fu.x());
289  const dimensionedScalar dv(sizeGroup_.x() - fu.x());
290 
291  const volScalarField da1
292  (
293  (2.0/3.0)*dv
294  *(
295  sourceKappa
296  + sourceShape.Df_*(1/sourceShape.d() - 1/dp)
297  )
298  );
299 
300  dp += 6*(dv*a - fu.x()*da1)/sqr(a);
301 
302  const volScalarField np(6*sizeGroup_.x()/pi/pow3(dp));
303  const volScalarField dc(dp*pow(np/alphaC_, 1/Df_));
304 
305  const volScalarField da2
306  (
307  dv*(4/dp + 2*Df_/3*(1/dc - 1/dp))
308  );
309 
310  Su_ += (a + 0.5*da1 + 0.5*da2)/sizeGroup_.x()*Su;
311 
312  break;
313  }
314 
315  case sgConserved:
316  {
317  SecondaryPropertyModel::addDrift(Su, fu, model);
318 
319  break;
320  }
321 
322  default:
323  {
325  << "Unknown surface growth type. Valid types are:"
326  << sgTypeNames_ << nl << exit(FatalError);
327  }
328  }
329 }
330 
331 
332 // ************************************************************************* //
surfaceGrowthTypes
Surface growth type enumeration.
Definition: fractal.H:112
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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)
virtual const volScalarField & fld() const
Return reference to secondary property field.
const dimensionSet dimless
static const NamedEnum< surfaceGrowthTypes, 3 > sgTypeNames_
Surface growth type names.
Definition: fractal.H:120
fvMesh & mesh
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
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:58
const dimensionSet dimTime
virtual const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
virtual volScalarField & src()
Access to secondary property source.
Calculate the matrix for the first temporal derivative.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
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.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
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
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 tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
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.