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 Group
29  QGDThermo
30 Class
31  Foam::QGDThermo
32 
33 Description
34  Abstract base class for classes implementing thermophysical properties
35  of gases and fluids governed by regularized equations (QGD & QHD)
36 
37 SourceFiles
38  QGDThermo.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef QGDTHERMO_H
43 #define QGDTHERMO_H
44 
45 #include "volFieldsFwd.H"
46 #include "surfaceFieldsFwd.H"
47 #include "Switch.H"
48 
49 namespace Foam
50 {
51 class fvMesh;
52 class dictionary;
53 namespace qgd
54 {
55  class QGDCoeffs;
56 }
57 
58 class QGDThermo
59 {
60 
61  //-
62  const fvMesh& mesh_;
63 
64  //-
65  const dictionary& dict_;
66 
67  //-
68  autoPtr<Foam::qgd::QGDCoeffs> qgdCoeffsPtr_;
69 
70  //-
71  Switch implicitDiffusion_;
72 
73 private:
74 
75  //-
76  QGDThermo();
77 
78  //-
79  QGDThermo(const QGDThermo& );
80 
81 protected:
82 
83  //-
85 
86  //-
87  void correctQGD(volScalarField& mu, volScalarField& alphau);
88 
89 public:
90 
91  TypeName ("QGDThermo");
92 
93  //- Construct from mesh and dictionary
94  QGDThermo(const fvMesh& mesh, const dictionary& dict);
95 
96  virtual ~QGDThermo();
97 
98  //-
99  virtual const volScalarField& c() const = 0;
100 
101  //-
102  virtual const volScalarField& p() const = 0;
103 
104  //-
105  virtual tmp<volScalarField> rho() const = 0;
106 
107  //-
108  virtual tmp<volScalarField> mu() const = 0;
109 
110  //-
111  const volScalarField& hQGD() const;
112 
113  //-
114  const volScalarField& tauQGD() const;
115 
116  //-
117  const surfaceScalarField& hQGDf() const;
118 
119  //-
120  const surfaceScalarField& tauQGDf() const;
121 
122  //-
123  const volScalarField& muQGD() const;
124 
125  //-
126  const volScalarField& alphauQGD() const;
127 
128  //-
129  const volScalarField& aQGD() const;
130 
131  //-
132  const surfaceScalarField& aQGDf() const;
133 
134  //-
135  Switch implicitDiffusion() const;
136 
137  //-
138  virtual bool read();
139 
140 };
141 }
142 
143 #endif
144 
145 //
146 //END-OF-FILE
147 //
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:53
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