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 -------------------------------------------------------------------------------
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::phaseModel
28 
29 SourceFiles
30  phaseModel.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef phaseModel_H
35 #define phaseModel_H
36 
37 #include "dictionary.H"
38 #include "dictionaryEntry.H"
39 #include "dimensionedScalar.H"
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 // Forward declarations
49 class diameterModel;
50 
51 /*---------------------------------------------------------------------------*\
52  Class phaseModel Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class phaseModel
56 :
57  public volScalarField
58 {
59  // Private data
60 
61  //- Name of phase
62  word name_;
63 
64  dictionary phaseDict_;
65 
66  //- Kinematic viscosity
68 
69  //- Thermal conductivity
70  dimensionedScalar kappa_;
71 
72  //- Heat capacity
74 
75  //- Density
76  dimensionedScalar rho_;
77 
78  //- Velocity
79  volVectorField U_;
80 
81  //- Substantive derivative of the velocity
82  volVectorField DDtU_;
83 
84  //- Volumetric flux of the phase
85  surfaceScalarField alphaPhi_;
86 
87  //- Volumetric flux for the phase
89 
90  //- Diameter model
92 
93 
94 public:
95 
96  // Constructors
97 
99  (
100  const word& phaseName,
101  const dictionary& phaseDict,
102  const fvMesh& mesh
103  );
104 
105  //- Return clone
106  autoPtr<phaseModel> clone() const;
107 
108  //- Return a pointer to a new phase created on freestore
109  // from Istream
110  class iNew
111  {
112  const fvMesh& mesh_;
113 
114  public:
115 
117  (
118  const fvMesh& mesh
119  )
120  :
121  mesh_(mesh)
122  {}
123 
125  {
127  return autoPtr<phaseModel>
128  (
129  new phaseModel(ent.keyword(), ent, mesh_)
130  );
131  }
132  };
133 
134 
135  //- Destructor
136  virtual ~phaseModel();
137 
138 
139  // Member Functions
140 
141  const word& name() const
142  {
143  return name_;
144  }
145 
146  const word& keyword() const
147  {
148  return name();
149  }
150 
151  tmp<volScalarField> d() const;
152 
153  const dimensionedScalar& nu() const
154  {
155  return nu_;
156  }
157 
158  const dimensionedScalar& kappa() const
159  {
160  return kappa_;
161  }
162 
163  const dimensionedScalar& Cp() const
164  {
165  return Cp_;
166  }
167 
168  const dimensionedScalar& rho() const
169  {
170  return rho_;
171  }
172 
173  const volVectorField& U() const
174  {
175  return U_;
176  }
177 
178  volVectorField& U()
179  {
180  return U_;
181  }
182 
183  const volVectorField& DDtU() const
184  {
185  return DDtU_;
186  }
187 
189  {
190  return DDtU_;
191  }
192 
193  const surfaceScalarField& phi() const
194  {
195  return *phiPtr_;
196  }
197 
199  {
200  return *phiPtr_;
201  }
202 
203  const surfaceScalarField& alphaPhi() const
204  {
205  return alphaPhi_;
206  }
207 
209  {
210  return alphaPhi_;
211  }
212 
213  //- Ensure that the flux at inflow/outflow BCs is preserved
215 
216  //- Correct the phase properties
217  void correct();
218 
219  //-Inherit read from volScalarField
220  using volScalarField::read;
221 
222  //- Read base transportProperties dictionary
223  bool read(const dictionary& phaseDict);
224 };
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 } // End namespace Foam
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 #endif
234 
235 // ************************************************************************* //
volFields.H
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::dictionaryEntry
A keyword and a list of tokens is a 'dictionaryEntry'.
Definition: dictionaryEntry.H:65
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::phaseModel::U
volVectorField & U()
Definition: phaseModel.H:177
Foam::phaseModel::kappa
const dimensionedScalar & kappa() const
Definition: phaseModel.H:157
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::phaseModel::d
tmp< volScalarField > d() const
Definition: phaseModel.C:257
Foam::phaseModel::alphaPhi
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:202
Foam::phaseModel::keyword
const word & keyword() const
Definition: phaseModel.H:145
Foam::phaseModel::correct
void correct()
Correct the phase properties.
Definition: phaseModel.C:215
Foam::entry::keyword
const keyType & keyword() const noexcept
Return keyword.
Definition: entry.H:195
surfaceFields.H
Foam::surfaceFields.
Foam::phaseModel::Cp
const dimensionedScalar & Cp() const
Definition: phaseModel.H:162
Foam::phaseModel::~phaseModel
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:201
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:392
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::phaseModel::phaseModel
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition: phaseModel.C:39
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::phaseModel::iNew::operator()
autoPtr< phaseModel > operator()(Istream &is) const
Definition: phaseModel.H:123
Foam::phaseModel::DDtU
const volVectorField & DDtU() const
Definition: phaseModel.H:182
Foam::phaseModel::phi
surfaceScalarField & phi()
Definition: phaseModel.H:197
Foam::phaseModel::nu
const dimensionedScalar & nu() const
Definition: phaseModel.H:152
Foam::phaseModel::iNew
Return a pointer to a new phase created on freestore.
Definition: phaseModel.H:109
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
dictionaryEntry.H
Foam::phaseModel::phi
const surfaceScalarField & phi() const
Definition: phaseModel.H:192
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
dimensionedScalar.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::phaseModel::U
const volVectorField & U() const
Definition: phaseModel.H:172
Foam::phaseModel::correctInflowOutflow
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:239
Foam::phaseModel::read
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:321
dictionary.H
Foam::phaseModel::DDtU
volVectorField & DDtU()
Definition: phaseModel.H:187
Foam::phaseModel::rho
const dimensionedScalar & rho() const
Definition: phaseModel.H:167
Foam::phaseModel::alphaPhi
surfaceScalarField & alphaPhi()
Definition: phaseModel.H:207
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::phaseModel::clone
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:207
Foam::phaseModel::iNew::iNew
iNew(const fvMesh &mesh)
Definition: phaseModel.H:116