All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
H2bynuQHD.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 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 
29 #include "H2bynuQHD.H"
30 #include "QGDThermo.H"
31 #include "addToRunTimeSelectionTable.H"
32 #include "linear.H"
33 
34 namespace Foam
35 {
36 namespace qgd
37 {
38  defineTypeNameAndDebug(H2bynuQHD,0);
40  (
41  QGDCoeffs,
42  H2bynuQHD,
43  dictionary
44  );
45 }
46 }
47 
50 (
51  const IOobject& io,
52  const fvMesh& mesh,
53  const dictionary& dict
54 )
55 :
56  QGDCoeffs(io, mesh, dict)
57 {
58  scalar ScQGD = 0.0, PrQGD = 1.0;
59 //
60 
61  ScQGD_.primitiveFieldRef() = ScQGD;
62  PrQGD_.primitiveFieldRef() = PrQGD;
63  muQGD_.primitiveFieldRef() = 0.0;
64  alphauQGD_.primitiveFieldRef() = 0.0;
65 
66  ScQGD_.boundaryFieldRef() = ScQGD;
67  PrQGD_.boundaryFieldRef() = PrQGD;
68  muQGD_.boundaryFieldRef() = 0.0;
69  alphauQGD_.boundaryFieldRef() = 0.0;
70 }
71 
74 {
75 }
76 
77 void Foam::qgd::
78 H2bynuQHD::correct(const QGDThermo& qgdThermo)
79 {
80  const volScalarField nu = qgdThermo.mu()/qgdThermo.rho();
81  this->tauQGD_ = this->aQGD_ * sqr(this->hQGD_) / nu;
82  this->tauQGDf_= linearInterpolate(this->tauQGD_);
83 }
84 
85 //
86 //END-OF-FILE
87 //
volScalarField PrQGD_
Definition: QGDCoeffs.H:115
volScalarField alphauQGD_
Definition: QGDCoeffs.H:80
virtual tmp< volScalarField > rho() const =0
H2bynuQHD(const IOobject &io, const fvMesh &mesh, const dictionary &dict)
Definition: H2bynuQHD.C:43
virtual tmp< volScalarField > mu() const =0
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
volScalarField muQGD_
Definition: QGDCoeffs.H:75
void correct(const QGDThermo &)
Definition: H2bynuQHD.C:71
addToRunTimeSelectionTable(QGDCoeffs, constScPrModel1, dictionary)
defineTypeNameAndDebug(constScPrModel1, 0)
volScalarField ScQGD_
Definition: QGDCoeffs.H:120