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-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::diameterModels::IATE
28
29Description
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 Reference:
39 \verbatim
40 Ishii, M., Kim, S., & Kelly, J. (2005).
41 Development of interfacial area transport equation.
42 Nuclear Engineering and Technology, 37(6), 525-536.
43 \endverbatim
44
45SourceFiles
46 IATE.C
47
48\*---------------------------------------------------------------------------*/
49
50#ifndef IATE_H
51#define IATE_H
52
53#include "diameterModel.H"
54
55// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57namespace Foam
58{
59namespace diameterModels
60{
61
62// Forward declaration of classes
63class IATEsource;
64
65/*---------------------------------------------------------------------------*\
66 Class IATE Declaration
67\*---------------------------------------------------------------------------*/
69class IATE
70:
71 public diameterModel
72{
73 // Private data
74
75 //- Interfacial curvature (alpha*interfacial area)
76 volScalarField kappai_;
77
78 //- Maximum diameter used for stabilisation in the limit kappai->0
80
81 //- Minimum diameter used for stabilisation in the limit kappai->inf
83
84 //- Residual phase fraction
85 dimensionedScalar residualAlpha_;
86
87 //- The Sauter-mean diameter of the phase
89
90 //- IATE sources
91 PtrList<IATEsource> sources_;
92
93
94 // Private member functions
95
96 tmp<volScalarField> dsm() const;
97
98
99public:
101 friend class IATEsource;
102
103 //- Runtime type information
104 TypeName("IATE");
105
106
107 // Constructors
108
109 //- Construct from components
110 IATE
111 (
113 const phaseModel& phase
114 );
115
116
117 //- Destructor
118 virtual ~IATE();
119
120
121 // Member Functions
122
123 //- Return the interfacial curvature
124 const volScalarField& kappai() const
125 {
126 return kappai_;
127 }
128
129 //- Return the interfacial area
130 tmp<volScalarField> a() const
131 {
132 return phase_*kappai_;
133 }
134
135 //- Return the Sauter-mean diameter
136 virtual tmp<volScalarField> d() const
137 {
138 return d_;
139 }
140
141 //- Correct the diameter field
142 virtual void correct();
143
144 //- Read phaseProperties dictionary
145 virtual bool read(const dictionary& phaseProperties);
146};
147
148
149// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150
151} // End namespace diameterModels
152} // End namespace Foam
153
154// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155
156#endif
157
158// ************************************************************************* //
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:54
const phaseModel & phase_
Definition: diameterModel.H:61
const dictionary & diameterProperties() const
Return the phase diameter properties dictionary.
IATE (Interfacial Area Transport Equation) bubble diameter model.
Definition: IATE.H:71
const volScalarField & kappai() const
Return the interfacial curvature.
Definition: IATE.H:123
virtual ~IATE()
Destructor.
Definition: IATE.C:106
virtual void correct()
Correct the diameter field.
Definition: IATE.C:117
TypeName("IATE")
Runtime type information.
tmp< volScalarField > a() const
Return the interfacial area.
Definition: IATE.H:129
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: IATE.C:171
virtual tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: IATE.H:135
IATE (Interfacial Area Transport Equation) bubble diameter model run-time selectable sources.
Definition: IATEsource.H:57
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Helper class to manage multi-specie phase properties.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73