KochFriedlanderSintering.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) 2024-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 "fvmSup.H"
28 #include "fractal.H"
30 #include "sizeGroup.H"
31 #include "volFieldsFwd.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
41 }
42 }
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 Foam::dimensionSet Foam::fv::KochFriedlanderSintering::CsDims() const
47 {
48  return pow(dimLength, -n_)*pow(dimTemperature, -m_)*dimTime;
49 }
50 
51 
52 void Foam::fv::KochFriedlanderSintering::readCoeffs(const dictionary& dict)
53 {
54  if (dict.lookup<word>("populationBalance") != popBal_.name())
55  {
57  << "Cannot change the population balance of a " << type()
58  << " model at run time" << exit(FatalIOError);
59  }
60 
61  n_ = dict.lookup<scalar>("n");
62  m_ = dict.lookup<scalar>("m");
63  Cs_.read(dict, CsDims());
64  Ta_.read(dict);
65  dpMin_.readIfPresent(dict);
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const word& name,
74  const word& modelType,
75  const fvMesh& mesh,
76  const dictionary& dict
77 )
78 :
79  fvModel(name, modelType, mesh, dict),
80  popBal_
81  (
82  mesh().lookupObject<diameterModels::populationBalanceModel>
83  (
84  coeffs(dict).lookup("populationBalance")
85  )
86  ),
87  n_(coeffs(dict).lookup<scalar>("n")),
88  m_(coeffs(dict).lookup<scalar>("m")),
89  Cs_("Cs", CsDims(), coeffs(dict)),
90  Ta_("Ta", dimTemperature, coeffs(dict)),
91  dpMin_("dpMin", dimLength, coeffs(dict), scalar(0))
92 {
93  readCoeffs(coeffs(dict));
94 
95  forAll(popBal_.sizeGroups(), sizeGroupi)
96  {
97  const diameterModels::shapeModel& shapeModel =
98  popBal_.sizeGroups()[sizeGroupi].shape();
99 
100  if (!isA<diameterModels::shapeModels::fractal>(shapeModel)) continue;
101 
102  const diameterModels::shapeModels::fractal& fractal =
103  refCast<const diameterModels::shapeModels::fractal>(shapeModel);
104 
105  kappaNameToSizeGroupIndices_.insert(fractal.fld().name(), sizeGroupi);
106  }
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 (
114  const word& fieldName
115 ) const
116 {
117  return kappaNameToSizeGroupIndices_.found(fieldName);
118 }
119 
120 
123 (
125 ) const
126 {
127  const diameterModels::sizeGroup& fi =
128  popBal_.sizeGroups()[kappaNameToSizeGroupIndices_[kappa.name()]];
129 
130  const volScalarField::Internal& T = fi.phase().thermo().T();
131 
132  const volScalarField::Internal dp(6/max(6/fi.dSph(), kappa));
133 
134  return Cs_*pow(dp, n_)*pow(T, m_)*exp(Ta_/T*(1 - dpMin_/dp));
135 }
136 
137 
139 (
140  const volScalarField& alphaFi,
141  const volScalarField& rho,
142  const volScalarField& kappa,
143  fvMatrix<scalar>& eqn
144 ) const
145 {
146  const diameterModels::sizeGroup& fi =
147  popBal_.sizeGroups()[kappaNameToSizeGroupIndices_[kappa.name()]];
148 
149  const volScalarField::Internal R(alphaFi()*rho()/tau(kappa));
150 
151  eqn += R*6/fi.dSph() - fvm::Sp(R, kappa);
152 }
153 
154 
156 {
157  return true;
158 }
159 
160 
162 {}
163 
164 
166 {}
167 
168 
170 {}
171 
172 
174 {
175  if (fvModel::read(dict))
176  {
177  readCoeffs(coeffs(dict));
178  return true;
179  }
180  else
181  {
182  return false;
183  }
184 }
185 
186 
187 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
const word & name() const
Return name.
Definition: IOobject.H:307
virtual const volScalarField & T() const =0
Temperature [K].
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
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:94
virtual const volScalarField & fld() const
Return reference to secondary property field.
Definition: fractal.C:139
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 phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
Sintering model of Koch and Friedlander (1990). The characteristic time for sintering is given by.
virtual bool movePoints()
Update for mesh motion.
virtual tmp< volScalarField::Internal > tau(const volScalarField::Internal &kappa) const
Return the characteristic time for sintering.
void addSup(const volScalarField &alphaFi, const volScalarField &rho, const volScalarField &kappa, fvMatrix< scalar > &eqn) const
Add a source term to the surface-area-volume-ratio equation.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
KochFriedlanderSintering(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
const vector tau
Calculate the matrix for implicit and explicit sources.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar exp(const dimensionedScalar &ds)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimLength
const dimensionSet dimTemperature
const dimensionSet dimTime
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList fv(nPoints)
dictionary dict