IATE.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) 2013-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 Class
25  Foam::diameterModels::IATE
26 
27 Description
28  IATE (Interfacial Area Transport Equation) bubble diameter model.
29 
30  Solves for the interfacial curvature per unit volume of the phase rather
31  than interfacial area per unit volume to avoid stability issues relating to
32  the consistency requirements between the phase fraction and interfacial area
33  per unit volume. In every other respect this model is as presented in the
34  paper:
35 
36  Reference:
37  \verbatim
38  Ishii, M., Kim, S., & Kelly, J. (2005).
39  Development of interfacial area transport equation.
40  Nuclear Engineering and Technology, 37(6), 525-536.
41  \endverbatim
42 
43 SourceFiles
44  IATE.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef IATE_H
49 #define IATE_H
50 
51 #include "diameterModel.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 namespace diameterModels
58 {
59 
60 // Forward declaration of classes
61 class IATEsource;
62 
63 /*---------------------------------------------------------------------------*\
64  Class IATE Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class IATE
68 :
69  public diameterModel
70 {
71  // Private Data
72 
73  //- Interfacial curvature (alpha*interfacial area)
74  volScalarField kappai_;
75 
76  //- Maximum diameter used for stabilisation in the limit kappai->0
77  dimensionedScalar dMax_;
78 
79  //- Minimum diameter used for stabilisation in the limit kappai->inf
80  dimensionedScalar dMin_;
81 
82  //- Residual phase fraction
83  dimensionedScalar residualAlpha_;
84 
85  //- The Sauter-mean diameter of the phase
86  volScalarField d_;
87 
88  //- IATE sources
89  PtrList<IATEsource> sources_;
90 
91 
92  // Private Member Functions
93 
94  tmp<volScalarField> dsm() const;
95 
96 
97 public:
98 
99  friend class IATEsource;
100 
101  //- Runtime type information
102  TypeName("IATE");
103 
104 
105  // Constructors
106 
107  //- Construct from dictionary and phase
108  IATE
109  (
111  const phaseModel& phase
112  );
113 
114 
115  //- Destructor
116  virtual ~IATE();
117 
118 
119  // Member Functions
120 
121  //- Return the interfacial curvature
122  inline const volScalarField& kappai() const
123  {
124  return kappai_;
125  }
126 
127  //- Get the diameter field
128  virtual tmp<volScalarField> d() const;
129 
130  //- Get the surface area per unit volume field
131  virtual tmp<volScalarField> Av() const;
132 
133  //- Correct the model
134  virtual void correct();
135 
136  //- Read phaseProperties dictionary
137  virtual bool read(const dictionary& phaseProperties);
138 };
139 
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 } // End namespace diameterModels
144 } // End namespace Foam
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 #endif
149 
150 // ************************************************************************* //
Generic GeometricField class.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
const phaseModel & phase() const
Return the phase.
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:69
const volScalarField & kappai() const
Return the interfacial curvature.
Definition: IATE.H:121
virtual ~IATE()
Destructor.
Definition: IATE.C:91
virtual void correct()
Correct the model.
Definition: IATE.C:109
TypeName("IATE")
Runtime type information.
virtual tmp< volScalarField > d() const
Get the diameter field.
Definition: IATE.C:97
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: IATE.C:173
virtual tmp< volScalarField > Av() const
Get the surface area per unit volume field.
Definition: IATE.C:103
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
Definition: IATE.C:63
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:54
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Helper class to manage multi-specie phase properties.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.