All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
QGDCoeffs.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 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Class
29  Foam::qgd::QGDCoeffs
30 
31 Group
32  grpQGDCoeffs
33 
34 Description
35  Base class for all classes describing QGD model coefficients.
36  Provides interfaces for accessing QGD/QHD models coefficients \mu^{QGD},
37  \tau^{QGD}, \alpha^{QGD}, \alpha_u^{QGD}, h^{QGD}
38 
39 SourceFiles
40  QGDCoeffs.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #include "fvMesh.H"
45 #include "runTimeSelectionTables.H"
46 #include "regIOobject.H"
47 #include "volFields.H"
48 #include "surfaceFields.H"
49 
50 #ifndef QGDCoeffs_H
51 #define QGDCoeffs_H
52 
53 namespace Foam
54 {
55 
56 //forward declaration of psiQGDThermo
57 class QGDThermo;
58 
59 namespace qgd
60 {
61 
62 class QGDCoeffs : public regIOobject, public refCount
63 {
64 
65 protected:
66 
67  //-
68  const fvMesh& mesh_;
69 
70  //-
71  const Time& runTime_;
72 
73  //-
74  volScalarField muQGD_;
75 
76  //-
77  volScalarField alphauQGD_;
78 
79  //-
80  surfaceScalarField hQGDf_;
81 
82  //-
83  volScalarField hQGD_;
84 
85  //-
86  tmp<volScalarField> taQGD_;
87 
88  //-
89  volScalarField& aQGD_;
90 
91  //-
92  volScalarField tauQGD_;
93 
94  //-
95  surfaceScalarField tauQGDf_;
96 
97  //-
98  volScalarField PrQGD_;
99 
100  //-
101  volScalarField ScQGD_;
102 
103 protected:
104 
105  //-
106  virtual void updateQGDLength(const fvMesh&);
107 
108  //-
109  tmp<volScalarField> readOrCreateAlphaQGD(const fvMesh& );
111 public:
112 
113  //-
114  TypeName("QGDCoeffs");
116  //-
118  (
119  autoPtr,
121  dictionary,
122  (
123  const IOobject& io,
124  const fvMesh& mesh,
125  const dictionary& dict
126  ),
127  (io, mesh, dict)
128  );
129 
130  //- Construct from components
131  QGDCoeffs
132  (
133  const IOobject& io,
134  const fvMesh& mesh,
135  const dictionary& dict
136  );
137 
138  //- Return a reference to the selected fvscStencil model
139  static autoPtr<QGDCoeffs> New
140  (
141  const word& qgdCoeffsType,
142  const fvMesh& mesh,
143  const dictionary& dict
144  );
145 
146  //-
147  virtual ~QGDCoeffs();
148 
149  //-
150  virtual bool writeData(Ostream&) const
151  {
152  return true;
153  }
154 
155  //-
156  virtual void correct(const QGDThermo&);
157 
158  //-
159  const volScalarField& hQGD() const;
160 
161  //-
162  const volScalarField& tauQGD() const;
163 
164  //-
165  const surfaceScalarField& hQGDf() const;
166 
167  //-
168  const surfaceScalarField& tauQGDf() const;
169 
170  //-
171  const volScalarField& muQGD() const;
172 
173  //-
174  const volScalarField& alphauQGD() const;
175 
176 };
177 
178 }
179 
180 }
181 
182 #endif
183 
184 //
185 //END-OF-FILE
186 //
volScalarField PrQGD_
Definition: QGDCoeffs.H:115
virtual void updateQGDLength(const fvMesh &)
Definition: QGDCoeffs.C:291
surfaceScalarField hQGDf_
Definition: QGDCoeffs.H:85
volScalarField alphauQGD_
Definition: QGDCoeffs.H:80
QGDCoeffs(const IOobject &io, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: QGDCoeffs.C:156
const volScalarField & muQGD() const
Definition: QGDCoeffs.C:386
Base class for all classes describing QGD model coefficients. Provides interfaces for accessing QGD/Q...
Definition: QGDCoeffs.H:57
const surfaceScalarField & hQGDf() const
Definition: QGDCoeffs.C:371
declareRunTimeSelectionTable(autoPtr, QGDCoeffs, dictionary,(const IOobject &io, const fvMesh &mesh, const dictionary &dict),(io, mesh, dict))
surfaceScalarField tauQGDf_
Definition: QGDCoeffs.H:110
TypeName("QGDCoeffs")
virtual bool writeData(Ostream &) const
Definition: QGDCoeffs.H:185
volScalarField & aQGD_
Definition: QGDCoeffs.H:100
static autoPtr< QGDCoeffs > New(const word &qgdCoeffsType, const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected fvscStencil model.
Definition: QGDCoeffs.C:52
const volScalarField & tauQGD() const
Definition: QGDCoeffs.C:391
const surfaceScalarField & tauQGDf() const
Definition: QGDCoeffs.C:396
volScalarField tauQGD_
Definition: QGDCoeffs.H:105
const volScalarField & alphauQGD() const
Definition: QGDCoeffs.C:381
volScalarField hQGD_
Definition: QGDCoeffs.H:90
tmp< volScalarField > readOrCreateAlphaQGD(const fvMesh &)
Definition: QGDCoeffs.C:112
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
volScalarField muQGD_
Definition: QGDCoeffs.H:75
virtual void correct(const QGDThermo &)
Definition: QGDCoeffs.C:263
const fvMesh & mesh_
Definition: QGDCoeffs.H:65
virtual ~QGDCoeffs()
Definition: QGDCoeffs.C:259
const volScalarField & hQGD() const
Definition: QGDCoeffs.C:376
volScalarField ScQGD_
Definition: QGDCoeffs.H:120
tmp< volScalarField > taQGD_
Definition: QGDCoeffs.H:95
const Time & runTime_
Definition: QGDCoeffs.H:70