Loading...
Searching...
No Matches
phaseModel.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::phaseModel
29
30SourceFiles
31 phaseModel.C
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef Foam_multiphaseEuler_phaseModel_H
36#define Foam_multiphaseEuler_phaseModel_H
37
38#include "dictionary.H"
39#include "dictionaryEntry.H"
40#include "dimensionedScalar.H"
41#include "volFields.H"
42#include "surfaceFields.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
49namespace multiphaseEuler
50{
51// Forward declarations
52class diameterModel;
54
55/*---------------------------------------------------------------------------*\
56 Class phaseModel Declaration
57\*---------------------------------------------------------------------------*/
58
59class phaseModel
60:
61 public volScalarField
62{
63 // Private data
64
65 //- Name of phase
66 word name_;
67
68 dictionary phaseDict_;
69
70 //- Kinematic viscosity
72
73 //- Thermal conductivity
74 dimensionedScalar kappa_;
75
76 //- Heat capacity
78
79 //- Density
81
82 //- Velocity
84
85 //- Substantive derivative of the velocity
86 volVectorField DDtU_;
87
88 //- Volumetric flux of the phase
89 surfaceScalarField alphaPhi_;
90
91 //- Volumetric flux for the phase
93
94 //- Diameter model
96
97
98public:
99
100 // Constructors
101
103 (
104 const word& phaseName,
105 const dictionary& phaseDict,
106 const fvMesh& mesh
107 );
108
109 //- Return clone
111
112 //- Return a pointer to a new phase created on freestore
113 // from Istream
114 class iNew
115 {
116 const fvMesh& mesh_;
117
118 public:
119
120 iNew
121 (
122 const fvMesh& mesh
123 )
124 :
125 mesh_(mesh)
126 {}
127
129 {
132 (
133 new phaseModel(ent.keyword(), ent, mesh_)
134 );
135 }
136 };
137
138
139 //- Destructor
140 virtual ~phaseModel();
141
142
143 // Member Functions
144
145 const word& name() const
146 {
147 return name_;
148 }
149
150 const word& keyword() const
151 {
152 return name();
153 }
154
155 tmp<volScalarField> d() const;
156
157 const dimensionedScalar& nu() const
158 {
159 return nu_;
160 }
161
162 const dimensionedScalar& kappa() const
163 {
164 return kappa_;
165 }
167 const dimensionedScalar& Cp() const
168 {
169 return Cp_;
170 }
172 const dimensionedScalar& rho() const
173 {
174 return rho_;
175 }
176
177 const volVectorField& U() const
179 return U_;
180 }
181
184 return U_;
185 }
186
187 const volVectorField& DDtU() const
189 return DDtU_;
190 }
191
194 return DDtU_;
195 }
196
197 const surfaceScalarField& phi() const
199 return *phiPtr_;
200 }
201
204 return *phiPtr_;
205 }
206
207 const surfaceScalarField& alphaPhi() const
209 return alphaPhi_;
210 }
211
214 return alphaPhi_;
215 }
216
217 //- Ensure that the flux at inflow/outflow BCs is preserved
219
220 //- Correct the phase properties
221 void correct();
222
223 //-Inherit read from volScalarField
225
226 //- Read base transportProperties dictionary
227 bool read(const dictionary& phaseDict);
229
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233} // End namespace Foam
234
235// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236
237#endif
238
239// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A keyword and a list of tokens is a 'dictionaryEntry'.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
const keyType & keyword() const noexcept
Return keyword.
Definition entry.H:231
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base-class for dispersed-phase particle diameter models.
iNew(const fvMesh &mesh)
Definition phaseModel.H:140
autoPtr< phaseModel > operator()(Istream &is) const
Definition phaseModel.H:147
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const volVectorField & U() const
Definition phaseModel.H:198
volVectorField & U()
Definition phaseModel.H:203
void correct()
Correct the phase properties.
Definition phaseModel.C:190
virtual ~phaseModel()
Destructor.
Definition phaseModel.C:176
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition phaseModel.C:33
const surfaceScalarField & phi() const
Definition phaseModel.H:218
surfaceScalarField & phi()
Definition phaseModel.H:223
volVectorField & DDtU()
Definition phaseModel.H:213
const dimensionedScalar & rho() const
Definition phaseModel.H:193
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition phaseModel.H:178
tmp< volScalarField > d() const
Definition phaseModel.C:233
const surfaceScalarField & alphaPhi() const
Definition phaseModel.H:228
const dimensionedScalar & Cp() const
Definition phaseModel.H:188
const word & keyword() const
Definition phaseModel.H:171
surfaceScalarField & alphaPhi()
Definition phaseModel.H:233
const volVectorField & DDtU() const
Definition phaseModel.H:208
autoPtr< phaseModel > clone() const
Return clone.
Definition phaseModel.C:182
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition phaseModel.C:214
const word & name() const
Definition phaseModel.H:166
virtual bool read()
Read phase properties dictionary.
Definition phaseModel.C:193
const dimensionedScalar & kappa() const
Definition phaseModel.H:183
virtual bool read()
Read object.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Foam::surfaceFields.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))