All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
QGDThermo.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 
13 License
14  This file is part of QGDsolver, based on OpenFOAM library.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Class
30  Foam::QGDThermo
31 
32 Description
33  Abstract base class for classes implementing thermophysical properties
34  of gases and fluids governed by regularized equations (QGD & QHD)
35 
36 SourceFiles
37  QGDThermo.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef QGDTHERMO_H
42 #define QGDTHERMO_H
43 
44 #include "volFieldsFwd.H"
45 #include "surfaceFieldsFwd.H"
46 #include "Switch.H"
47 
48 namespace Foam
49 {
50 class fvMesh;
51 class dictionary;
52 namespace qgd
53 {
54  class QGDCoeffs;
55 }
56 
57 class QGDThermo
58 {
59 
60  //-
61  const fvMesh& mesh_;
62 
63  //-
64  const dictionary& dict_;
65 
66  //-
67  autoPtr<Foam::qgd::QGDCoeffs> qgdCoeffsPtr_;
68 
69  //-
70  Switch implicitDiffusion_;
71 
72 private:
73 
74  //-
75  QGDThermo();
76 
77  //-
78  QGDThermo(const QGDThermo& );
79 
80 protected:
81 
82  //-
84 
85  //-
86  void correctQGD(volScalarField& mu, volScalarField& alphau);
87 
88 public:
89 
90  TypeName ("QGDThermo");
91 
92  //- Construct from mesh and dictionary
93  QGDThermo(const fvMesh& mesh, const dictionary& dict);
94 
95  virtual ~QGDThermo();
96 
97  //-
98  virtual const volScalarField& c() const = 0;
99 
100  //-
101  virtual const volScalarField& p() const = 0;
102 
103  //-
104  virtual tmp<volScalarField> rho() const = 0;
105 
106  //-
107  virtual tmp<volScalarField> mu() const = 0;
108 
109  //-
110  const volScalarField& hQGD() const;
111 
112  //-
113  const volScalarField& tauQGD() const;
114 
115  //-
116  const surfaceScalarField& hQGDf() const;
117 
118  //-
119  const surfaceScalarField& tauQGDf() const;
120 
121  //-
122  const volScalarField& muQGD() const;
123 
124  //-
125  const volScalarField& alphauQGD() const;
126 
127  //-
128  const volScalarField& aQGD() const;
129 
130  //-
131  const surfaceScalarField& aQGDf() const;
132 
133  //-
134  Switch implicitDiffusion() const;
135 
136  //-
137  virtual bool read();
138 
139 };
140 }
141 
142 #endif
143 
144 //
145 //END-OF-FILE
146 //
const volScalarField & tauQGD() const
Definition: QGDThermo.C:113
const volScalarField & aQGD() const
const surfaceScalarField & aQGDf() const
virtual ~QGDThermo()
Definition: QGDThermo.C:60
virtual bool read()
Definition: QGDThermo.C:65
virtual const volScalarField & c() const =0
const surfaceScalarField & hQGDf() const
Definition: QGDThermo.C:118
Base class for all classes describing QGD model coefficients. Provides interfaces for accessing QGD/Q...
Definition: QGDCoeffs.H:57
Switch implicitDiffusion() const
Definition: QGDThermo.C:138
const volScalarField & muQGD() const
Definition: QGDThermo.C:128
virtual tmp< volScalarField > rho() const =0
const volScalarField & alphauQGD() const
Definition: QGDThermo.C:133
TypeName("QGDThermo")
virtual tmp< volScalarField > mu() const =0
qgd::QGDCoeffs & qgdCoeffs()
Definition: QGDThermo.C:143
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:52
virtual const volScalarField & p() const =0
const volScalarField & hQGD() const
Definition: QGDThermo.C:108
void correctQGD(volScalarField &mu, volScalarField &alphau)
Definition: QGDThermo.C:79
const surfaceScalarField & tauQGDf() const
Definition: QGDThermo.C:123