KochFriedlanderSintering.H
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 Class
25  Foam::fv::KochFriedlanderSintering
26 
27 Description
28  Sintering model of Koch and Friedlander (1990). The characteristic time for
29  sintering is given by
30 
31  \f[
32  \tau = c_s d_{p_i}^n T^m \exp(T_a/T \cdot [1 - d_{p,min}/d_{p_i}])\;.
33  \f]
34 
35  Note that the correction factor in the exponential function can be
36  eliminated by setting \f$d_{p,min}\f$ to zero which is done by default.
37 
38  Reference:
39  \verbatim
40  Koch, W., & Friedlander, S. K. (1990).
41  The effect of particle coalescence on the surface area of a coagulating
42  aerosol.
43  Journal of Colloid and Interface Science, 140(2), 419-427.
44  \endverbatim
45 
46 Usage
47  \table
48  Property | Description | Required | Default value
49  Cs | Sintering time coefficient | yes | none
50  n | Particle diameter exponent | yes | none
51  m | Temperature exponent | yes | none
52  Ta | Activation temperature | yes | none
53  dpMin | Minimum primary particle diameter | no | 0
54  \endtable
55 
56  Example usage:
57  \verbatim
58  sintering
59  {
60  type KochFriedlanderSintering;
61  libs ("libmultiphaseEulerFvModels.so");
62 
63  populationBalance aggregates;
64 
65  Cs 8.3e24;
66  n 4.0;
67  m 1.0;
68  Ta 3700.0;
69  }
70  \endverbatim
71 
72 SourceFiles
73  KochFriedlanderSintering.C
74 
75 \*---------------------------------------------------------------------------*/
76 
77 #ifndef KochFriedlanderSintering_H
78 #define KochFriedlanderSintering_H
79 
80 #include "populationBalanceModel.H"
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 namespace Foam
85 {
86 namespace fv
87 {
88 
89 /*---------------------------------------------------------------------------*\
90  Class KochFriedlanderSintering Declaration
91 \*---------------------------------------------------------------------------*/
92 
93 class KochFriedlanderSintering
94 :
95  public fvModel
96 {
97 private:
98 
99  // Private Data
100 
101  //- Population balance model
102  const diameterModels::populationBalanceModel& popBal_;
103 
104  //- Sintering time primary particle diameter exponent
105  scalar n_;
106 
107  //- Sintering time temperature exponent
108  scalar m_;
109 
110  //- Sintering time coefficient
111  dimensionedScalar Cs_;
112 
113  //- Activation temperature
114  dimensionedScalar Ta_;
115 
116  //- Minimum primary particle diameter, defaults to Zero
117  dimensionedScalar dpMin_;
118 
119  //- Map from the name of the surface-area-volume-ratio fields to the
120  // index of the size group
121  HashTable<label> kappaNameToSizeGroupIndices_;
122 
123 
124  // Private Member Functions
125 
126  //- Return the dimensions of the sintering time coefficient
127  dimensionSet CsDims() const;
128 
129  //- Non-virtual read
130  void readCoeffs(const dictionary& dict);
131 
132 
133 public:
134 
135  //- Runtime type information
136  TypeName("KochFriedlanderSintering");
137 
138 
139  // Constructors
140 
141  //- Construct from explicit source name and mesh
143  (
144  const word& name,
145  const word& modelType,
146  const fvMesh& mesh,
147  const dictionary& dict
148  );
149 
150 
151  // Member Functions
152 
153  // Checks
154 
155  //- Return true if the fvModel adds a source term to the given
156  // field's transport equation
157  virtual bool addsSupToField(const word& fieldName) const;
158 
159 
160  // Sources
161 
162  //- Return the characteristic time for sintering
164  (
166  ) const;
167 
168  //- Add a source term to the surface-area-volume-ratio equation
169  void addSup
170  (
171  const volScalarField& alphaFi,
172  const volScalarField& rho,
173  const volScalarField& kappa,
174  fvMatrix<scalar>& eqn
175  ) const;
176 
177 
178  // Mesh changes
179 
180  //- Update for mesh motion
181  virtual bool movePoints();
182 
183  //- Update topology using the given map
184  virtual void topoChange(const polyTopoChangeMap&);
185 
186  //- Update from another mesh using the given map
187  virtual void mapMesh(const polyMeshMap&);
188 
189  //- Redistribute or update using the given distribution map
190  virtual void distribute(const polyDistributionMap&);
191 
192 
193  // IO
194 
195  //- Read source dictionary
196  virtual bool read(const dictionary& dict);
197 };
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 } // End namespace fv
203 } // End namespace Foam
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 #endif
208 
209 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
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
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
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
static const word & fieldName()
Return the name of the field associated with a source term. Special.
Definition: fvModelI.H:49
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.
TypeName("KochFriedlanderSintering")
Runtime type information.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
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
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict