All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
QGDThermo.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  QGDThermo
31 
32 SourceFiles
33  QGDThermo.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "QGDThermo.H"
38 #include "volFields.H"
39 #include "surfaceFields.H"
40 #include "fvMesh.H"
41 #include "QGDCoeffs.H"
42 
43 namespace Foam
44 {
46 }
47 
48 Foam::QGDThermo::QGDThermo(const fvMesh& mesh, const dictionary& dict)
49 :
50 mesh_(mesh),
51 dict_(dict),
52 qgdCoeffsPtr_
53 (
54  Foam::qgd::QGDCoeffs::New
55  (
56  dict_.subDict("QGD").get<word>("QGDCoeffs"),
57  mesh,
58  dict_.subDict("QGD")
59  )
60 ),
61 implicitDiffusion_(true)
62 {
63 }
64 
66 {
67 
68 }
69 
71 {
72  if (dict_.subDict("QGD").found("implicitDiffusion"))
73  {
74  dict_.subDict("QGD").lookup("implicitDiffusion") >> implicitDiffusion_;
75  }
76  else
77  {
78  implicitDiffusion_ = true;
79  }
80 
81  return true;
82 }
83 
84 void Foam::QGDThermo::correctQGD(volScalarField& mu, volScalarField& alphau)
85 {
86  qgdCoeffsPtr_->correct(*this);
87 
88  const volScalarField& muQGD = this->muQGD();
89  const volScalarField& alphauQGD = this->alphauQGD();
90 
91  forAll(mu.primitiveField(), celli)
92  {
93  mu.primitiveFieldRef()[celli] +=
94  muQGD.primitiveField()[celli];
95 
96  alphau.primitiveFieldRef()[celli] +=
97  alphauQGD.primitiveField()[celli];
98  }
99 
100  forAll(mu.boundaryField(), patchi)
101  {
102  forAll(mu.boundaryField()[patchi], facei)
103  {
104  mu.boundaryFieldRef()[patchi][facei] +=
105  muQGD.boundaryField()[patchi][facei];
106 
107  alphau.boundaryFieldRef()[patchi][facei] +=
108  alphauQGD.boundaryField()[patchi][facei];
109  }
110  }
111 }
112 
113 const Foam::volScalarField& Foam::QGDThermo::hQGD() const
114 {
115  return qgdCoeffsPtr_->hQGD();
116 }
117 
118 const Foam::volScalarField& Foam::QGDThermo::tauQGD() const
119 {
120  return qgdCoeffsPtr_->tauQGD();
121 }
122 
123 const Foam::surfaceScalarField& Foam::QGDThermo::hQGDf() const
124 {
125  return qgdCoeffsPtr_->hQGDf();
126 }
127 
128 const Foam::surfaceScalarField& Foam::QGDThermo::tauQGDf() const
129 {
130  return qgdCoeffsPtr_->tauQGDf();
131 }
132 
133 const Foam::volScalarField& Foam::QGDThermo::muQGD() const
134 {
135  return qgdCoeffsPtr_->muQGD();
136 }
137 
138 const Foam::volScalarField& Foam::QGDThermo::alphauQGD() const
139 {
140  return qgdCoeffsPtr_->alphauQGD();
141 }
142 
144 {
145  return implicitDiffusion_;
146 }
147 
149 {
150  return qgdCoeffsPtr_();
151 }
152 
153 //
154 //END-OF-FILE
155 //
const volScalarField & tauQGD() const
Definition: QGDThermo.C:113
virtual ~QGDThermo()
Definition: QGDThermo.C:60
virtual bool read()
Definition: QGDThermo.C:65
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
const volScalarField & alphauQGD() const
Definition: QGDThermo.C:133
forAll(Y, i)
Definition: QGDYEqn.H:36
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
defineTypeNameAndDebug(psiQGDReactionThermo, 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