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-2018 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  \verbatim
37  "Development of Interfacial Area Transport Equation"
38  Ishii, M., Kim, S. and Kelly, J.,
39  Nuclear Engineering and Technology, Vol.37 No.6 December 2005
40  \endverbatim
41 
42 SourceFiles
43  IATE.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef IATE_H
48 #define IATE_H
49 
50 #include "diameterModel.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 namespace diameterModels
57 {
58 
59 // Forward declaration of classes
60 class IATEsource;
61 
62 /*---------------------------------------------------------------------------*\
63  Class IATE Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class IATE
67 :
68  public diameterModel
69 {
70  // Private data
71 
72  //- Interfacial curvature (alpha*interfacial area)
73  volScalarField kappai_;
74 
75  //- Maximum diameter used for stabilisation in the limit kappai->0
76  dimensionedScalar dMax_;
77 
78  //- Minimum diameter used for stabilisation in the limit kappai->inf
79  dimensionedScalar dMin_;
80 
81  //- Residual phase fraction
82  dimensionedScalar residualAlpha_;
83 
84  //- The Sauter-mean diameter of the phase
85  volScalarField d_;
86 
87  //- IATE sources
88  PtrList<IATEsource> sources_;
89 
90 
91  // Private member functions
92 
93  tmp<volScalarField> dsm() const;
94 
95 
96 public:
97 
98  friend class IATEsource;
99 
100  //- Runtime type information
101  TypeName("IATE");
102 
103 
104  // Constructors
105 
106  //- Construct from components
107  IATE
108  (
109  const dictionary& diameterProperties,
110  const phaseModel& phase
111  );
112 
113 
114  //- Destructor
115  virtual ~IATE();
116 
117 
118  // Member Functions
119 
120  //- Return the interfacial curvature
121  const volScalarField& kappai() const
122  {
123  return kappai_;
124  }
125 
126  //- Return the interfacial area
127  tmp<volScalarField> a() const
128  {
129  return phase_*kappai_;
130  }
131 
132  //- Return the Sauter-mean diameter
133  virtual tmp<volScalarField> d() const
134  {
135  return d_;
136  }
137 
138  //- Correct the diameter field
139  virtual void correct();
140 
141  //- Read phaseProperties dictionary
142  virtual bool read(const dictionary& phaseProperties);
143 };
144 
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 } // End namespace diameterModels
149 } // End namespace Foam
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 #endif
154 
155 // ************************************************************************* //
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
TypeName("IATE")
Runtime type information.
Helper class to manage multi-specie phase properties.
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: IATE.H:132
const volScalarField & kappai() const
Return the interfacial curvature.
Definition: IATE.H:120
friend class IATEsource
Definition: IATE.H:97
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
const phaseModel & phase() const
Return the phase.
virtual ~IATE()
Destructor.
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< volScalarField > a() const
Return the interfacial area.
Definition: IATE.H:126
virtual void correct()
Correct the diameter field.
A class for managing temporary objects.
Definition: PtrList.H:53
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Namespace for OpenFOAM.
const phaseModel & phase_
Definition: diameterModel.H:58