All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
rhoQGDThermo.C
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 Group
30  grpRhoQGDThermo
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "rhoQGDThermo.H"
35 #include "QGDCoeffs.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 Foam::rhoQGDThermo::rhoQGDThermo(const fvMesh& mesh, const word& phaseName)
49 :
50  rhoThermo(mesh, phaseName),
51  QGDThermo(mesh, *this),
52  c_
53  (
54  "thermo:c",
55  qgdCoeffs().hQGD()/mesh.time().deltaT()
56  ),
57  gamma_
58  (
59  "thermo:gamma",
60  c_ / c_
61  )
62 {
63  this->read();
64 }
65 
66 Foam::rhoQGDThermo::rhoQGDThermo(const fvMesh& mesh, const word& phaseName, const word& dictName)
67 :
68  rhoThermo(mesh, phaseName,dictName),
69  QGDThermo(mesh, *this),
70  c_
71  (
72  "thermo:c",
73  qgdCoeffs().hQGD()/mesh.time().deltaT()
74  ),
75  gamma_
76  (
77  "thermo:gamma",
78  c_ / c_
79  )
80 {
81  this->read();
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
86 
87 Foam::autoPtr<Foam::rhoQGDThermo> Foam::rhoQGDThermo::New
88 (
89  const fvMesh& mesh,
90  const word& phaseName
91 )
92 {
93  return basicThermo::New<rhoQGDThermo>(mesh, phaseName);
94 }
95 
96 Foam::autoPtr<Foam::rhoQGDThermo> Foam::rhoQGDThermo::New
97 (
98  const fvMesh& mesh,
99  const word& phaseName,
100  const word& dictName
101 )
102 {
103  return basicThermo::New<rhoQGDThermo>(mesh, phaseName);
104 }
105 
106 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
107 
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 
115 
117 {
118  if (!basicThermo::read())
119  {
120  return false;
121  }
122 
123  if (!QGDThermo::read())
124  {
125  return false;
126  }
127 
128  return true;
129 }
130 
131 const Foam::volScalarField& Foam::rhoQGDThermo::c() const
132 {
133  return this->c_;
134 }
135 
136 const Foam::volScalarField& Foam::rhoQGDThermo::p() const
137 {
138  return rhoThermo::p();
139 }
140 
141 Foam::volScalarField& Foam::rhoQGDThermo::p()
142 {
143  return rhoThermo::p();
144 }
145 
146 Foam::tmp<Foam::volScalarField> Foam::rhoQGDThermo::rho() const
147 {
148  return rhoThermo::rho();
149 }
150 
151 Foam::tmp<Foam::volScalarField> Foam::rhoQGDThermo::mu() const
152 {
153  return rhoThermo::mu();
154 }
155 
156 // ************************************************************************* //
virtual tmp< volScalarField > rho() const
Definition: rhoQGDThermo.C:139
virtual bool read()
Definition: QGDThermo.C:65
virtual ~rhoQGDThermo()
Destructor.
Definition: rhoQGDThermo.C:101
virtual const volScalarField & p() const
Definition: rhoQGDThermo.C:129
const volScalarField & c() const
Definition: rhoQGDThermo.C:124
Basic thermodynamic properties based on density.
Definition: rhoQGDThermo.H:55
defineRunTimeSelectionTable(psiQGDReactionThermo, fvMesh)
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
defineTypeNameAndDebug(psiQGDReactionThermo, 0)
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())
static autoPtr< rhoQGDThermo > New(const fvMesh &mesh, const word &phaseName=word::null)
Selector.
Definition: rhoQGDThermo.C:81
rhoQGDThermo(const rhoQGDThermo &)
Construct as copy (not implemented)
volScalarField & p
Definition: createFields.H:52
virtual tmp< volScalarField > mu() const
Definition: rhoQGDThermo.C:144