IATE.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2014 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::diameterModels::IATE
28 
29 Description
30  IATE (Interfacial Area Transport Equation) bubble diameter model.
31 
32  Solves for the interfacial curvature per unit volume of the phase rather
33  than interfacial area per unit volume to avoid stability issues relating to
34  the consistency requirements between the phase fraction and interfacial area
35  per unit volume. In every other respect this model is as presented in the
36  paper:
37 
38  \verbatim
39  "Development of Interfacial Area Transport Equation"
40  Ishii, M., Kim, S. and Kelly, J.,
41  Nuclear Engineering and Technology, Vol.37 No.6 December 2005
42  \endverbatim
43 
44 SourceFiles
45  IATE.C
46 
47 \*---------------------------------------------------------------------------*/
48 
49 #ifndef IATE_H
50 #define IATE_H
51 
52 #include "diameterModel.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 namespace diameterModels
59 {
60 
61 // Forward declaration of classes
62 class IATEsource;
63 
64 /*---------------------------------------------------------------------------*\
65  Class IATE Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class IATE
69 :
70  public diameterModel
71 {
72  // Private data
73 
74  //- Interfacial curvature (alpha*interfacial area)
75  volScalarField kappai_;
76 
77  //- Maximum diameter used for stabilisation in the limit kappai->0
78  dimensionedScalar dMax_;
79 
80  //- Minimum diameter used for stabilisation in the limit kappai->inf
81  dimensionedScalar dMin_;
82 
83  //- Residual phase fraction
84  dimensionedScalar residualAlpha_;
85 
86  //- The Sauter-mean diameter of the phase
87  volScalarField d_;
88 
89  //- IATE sources
90  PtrList<IATEsource> sources_;
91 
92 
93  // Private member functions
94 
95  tmp<volScalarField> dsm() const;
96 
97 
98 public:
99 
100  friend class IATEsource;
101 
102  //- Runtime type information
103  TypeName("IATE");
104 
105 
106  // Constructors
107 
108  //- Construct from components
109  IATE
110  (
111  const dictionary& diameterProperties,
112  const phaseModel& phase
113  );
114 
115 
116  //- Destructor
117  virtual ~IATE();
118 
119 
120  // Member Functions
121 
122  //- Return the interfacial curvature
123  const volScalarField& kappai() const
124  {
125  return kappai_;
126  }
127 
128  //- Return the interfacial area
129  tmp<volScalarField> a() const
130  {
131  return phase_*kappai_;
132  }
133 
134  //- Return the Sauter-mean diameter
135  virtual tmp<volScalarField> d() const
136  {
137  return d_;
138  }
139 
140  //- Correct the diameter field
141  virtual void correct();
142 
143  //- Read phaseProperties dictionary
144  virtual bool read(const dictionary& phaseProperties);
145 };
146 
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 } // End namespace diameterModels
151 } // End namespace Foam
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 #endif
156 
157 // ************************************************************************* //
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::diameterModel::diameterProperties
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
Definition: diameterModel.H:110
Foam::phaseProperties
Helper class to manage multi-specie phase properties.
Definition: phaseProperties.H:60
Foam::diameterModels::IATE::TypeName
TypeName("IATE")
Runtime type information.
Foam::diameterModels::IATE::d
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: IATE.H:134
Foam::diameterModels::IATE::kappai
const volScalarField & kappai() const
Return the interfacial curvature.
Definition: IATE.H:122
Foam::diameterModels::IATE::correct
virtual void correct()
Correct the diameter field.
Definition: IATE.C:117
Foam::diameterModels::IATE::read
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: IATE.C:171
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::diameterModels::IATE::a
tmp< volScalarField > a() const
Return the interfacial area.
Definition: IATE.H:128
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::diameterModel::phase_
const phaseModel & phase_
Definition: diameterModel.H:61
Foam::diameterModels::IATE::~IATE
virtual ~IATE()
Destructor.
Definition: IATE.C:106
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::IATE::IATEsource
friend class IATEsource
Definition: IATE.H:100
Foam::diameterModel::phase
const phaseModel & phase() const
Return the phase.
Definition: diameterModel.H:116
Foam::diameterModels::IATE::IATE
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
Definition: IATE.C:63