phaseModel.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) 2011-2019 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::phaseModel
26 
27 SourceFiles
28  phaseModel.C
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef phaseModel_H
33 #define phaseModel_H
34 
35 #include "dictionary.H"
36 #include "dictionaryEntry.H"
37 #include "dimensionedScalar.H"
38 #include "volFields.H"
39 #include "surfaceFields.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // Forward declarations
47 class diameterModel;
48 
49 /*---------------------------------------------------------------------------*\
50  Class phaseModel Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 class phaseModel
54 :
55  public volScalarField
56 {
57  // Private Data
58 
59  //- Name of phase
60  word name_;
61 
62  dictionary phaseDict_;
63 
64  //- Kinematic viscosity
66 
67  //- Thermal conductivity
68  dimensionedScalar kappa_;
69 
70  //- Heat capacity
72 
73  //- Density
74  dimensionedScalar rho_;
75 
76  //- Velocity
77  volVectorField U_;
78 
79  //- Substantive derivative of the velocity
80  volVectorField DDtU_;
81 
82  //- Volumetric flux of the phase
83  surfaceScalarField alphaPhi_;
84 
85  //- Volumetric flux for the phase
86  autoPtr<surfaceScalarField> phiPtr_;
87 
88  //- Diameter model
89  autoPtr<diameterModel> dPtr_;
90 
91 
92 public:
93 
94  // Constructors
95 
97  (
98  const word& phaseName,
99  const dictionary& phaseDict,
100  const fvMesh& mesh
101  );
102 
103  //- Return clone
104  autoPtr<phaseModel> clone() const;
105 
106  //- Return a pointer to a new phase created on freestore
107  // from Istream
108  class iNew
109  {
110  const fvMesh& mesh_;
111 
112  public:
113 
115  (
116  const fvMesh& mesh
117  )
118  :
119  mesh_(mesh)
120  {}
123  {
125  return autoPtr<phaseModel>
126  (
127  new phaseModel(ent.keyword(), ent, mesh_)
128  );
129  }
130  };
131 
132 
133  //- Destructor
134  virtual ~phaseModel();
135 
136 
137  // Member Functions
139  const word& name() const
140  {
141  return name_;
142  }
144  const word& keyword() const
145  {
146  return name();
147  }
148 
149  tmp<volScalarField> d() const;
151  const dimensionedScalar& nu() const
152  {
153  return nu_;
154  }
156  const dimensionedScalar& kappa() const
157  {
158  return kappa_;
159  }
161  const dimensionedScalar& Cp() const
162  {
163  return Cp_;
164  }
166  const dimensionedScalar& rho() const
167  {
168  return rho_;
169  }
171  const volVectorField& U() const
172  {
173  return U_;
174  }
176  volVectorField& U()
177  {
178  return U_;
179  }
181  const volVectorField& DDtU() const
182  {
183  return DDtU_;
184  }
187  {
188  return DDtU_;
189  }
191  const surfaceScalarField& phi() const
192  {
193  return phiPtr_();
194  }
197  {
198  return phiPtr_();
199  }
201  const surfaceScalarField& alphaPhi() const
202  {
203  return alphaPhi_;
204  }
207  {
208  return alphaPhi_;
209  }
210 
211  //- Ensure that the flux at inflow/outflow BCs is preserved
213 
214  //- Correct the phase properties
215  void correct();
216 
217  //-Inherit read from volScalarField
218  using volScalarField::read;
219 
220  //- Read base transportProperties dictionary
221  bool read(const dictionary& phaseDict);
222 };
223 
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 } // End namespace Foam
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 #endif
232 
233 // ************************************************************************* //
Foam::surfaceFields.
const dimensionedScalar & Cp() const
Definition: phaseModel.H:160
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
const keyType & keyword() const
Return keyword.
Definition: entry.H:123
virtual bool read()
Read object.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static const dictionary null
Null dictionary.
Definition: dictionary.H:223
autoPtr< phaseModel > operator()(Istream &is) const
Definition: phaseModel.H:100
autoPtr< phaseModel > clone() const
Return clone.
const word & name() const
Definition: phaseModel.H:109
tmp< volScalarField > d() const
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const surfaceScalarField & phi() const
Definition: phaseModel.H:190
A keyword and a list of tokens is a &#39;dictionaryEntry&#39;.
iNew(const volScalarField &p, const volScalarField &T)
Definition: phaseModel.H:91
const volVectorField & DDtU() const
Definition: phaseModel.H:180
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:200
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const volVectorField & U() const
Definition: phaseModel.H:170
virtual ~phaseModel()
Destructor.
virtual bool read()
Read phase properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
const word & keyword() const
Definition: phaseModel.H:114
void correct()
Correct the laminar viscosity.
const dimensionedScalar & kappa() const
Definition: phaseModel.H:155
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
const Mesh & mesh() const
Return mesh.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.