AnisothermalPhaseModel.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 | Copyright (C) 2015-2017 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::AnisothermalPhaseModel
26 
27 Description
28  Class which represents a phase for which the temperature (strictly energy)
29  varies. Returns the energy equation and corrects the thermodynamic model.
30 
31 SourceFiles
32  AnisothermalPhaseModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef AnisothermalPhaseModel_H
37 #define AnisothermalPhaseModel_H
38 
39 #include "phaseModel.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class phaseModel Declaration
48 \*---------------------------------------------------------------------------*/
49 
50 template<class BasePhaseModel>
52 :
53  public BasePhaseModel
54 {
55  // Private data
56 
57  //- Kinetic energy
58  volScalarField K_;
59 
60 
61  // Private member functions
62 
63  //- Optionally filter the pressure work term as the phase-fraction -> 0
64  tmp<volScalarField> filterPressureWork
65  (
66  const tmp<volScalarField>& pressureWork
67  ) const;
68 
69 
70 public:
71 
72  // Constructors
73 
75  (
76  const phaseSystem& fluid,
77  const word& phaseName,
78  const label index
79  );
80 
81 
82  //- Destructor
83  virtual ~AnisothermalPhaseModel();
84 
85 
86  // Member Functions
87 
88  //- Return true if the phase is compressible otherwise false
89  virtual bool compressible() const;
90 
91  //- Correct the kinematics
92  virtual void correctKinematics();
93 
94  //- Correct the thermodynamics
95  virtual void correctThermo();
96 
97  //- Return the enthalpy equation
98  virtual tmp<fvScalarMatrix> heEqn();
99 
100  //- Return the phase kinetic energy
101  virtual const volScalarField& K() const;
102 };
103 
104 
105 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
106 
107 } // End namespace Foam
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 #ifdef NoRepository
112  #include "AnisothermalPhaseModel.C"
113 #endif
114 
115 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116 
117 #endif
118 
119 // ************************************************************************* //
virtual void correctThermo()
Correct the thermodynamics.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
multiphaseSystem & fluid
Definition: createFields.H:11
Class which represents a phase for which the temperature (strictly energy) varies. Returns the energy equation and corrects the thermodynamic model.
virtual const volScalarField & K() const
Return the phase kinetic energy.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
A class for handling words, derived from string.
Definition: word.H:59
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
virtual void correctKinematics()
Correct the kinematics.
virtual bool compressible() const
Return true if the phase is compressible otherwise false.
virtual ~AnisothermalPhaseModel()
Destructor.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
Namespace for OpenFOAM.