All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
constScPrModel1.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 
30 #include "QGDThermo.H"
31 #include "addToRunTimeSelectionTable.H"
32 #include "linear.H"
33 
34 namespace Foam
35 {
36 namespace qgd
37 {
40  (
41  QGDCoeffs,
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 PrQGD = 1.0;
59  if (dict.found("PrQGD"))
60  {
61  dict.lookup("PrQGD") >> PrQGD;
62  }
63  PrQGD_.primitiveFieldRef() = PrQGD;
64  PrQGD_.boundaryFieldRef() = PrQGD;
65 
66  IOobject ScHeader
67  (
68  "ScQGD",
69  mesh.time().timeName(),
70  mesh,
71  IOobject::READ_IF_PRESENT,
72  IOobject::NO_WRITE
73  );
74 
75  if (ScHeader.typeHeaderOk<volScalarField>())
76  {
77  //do nothing, ScQGD field is present
78  ScQGD_.writeOpt() = IOobject::AUTO_WRITE;
79  }
80  else
81  {
82  scalar ScQGD = 1.0;
83  if (dict.found("ScQGD"))
84  {
85  dict.lookup("ScQGD") >> ScQGD;
86  }
87  ScQGD_.primitiveFieldRef() = ScQGD;
88  ScQGD_.boundaryFieldRef() = ScQGD;
89  }
90 }
91 
94 {
95 }
96 
97 void Foam::qgd::
99 {
100  const volScalarField& cSound = qgdThermo.c();
101  const volScalarField& p = qgdThermo.p();
102 
103  this->tauQGDf_= linearInterpolate(this->aQGD_ / cSound) * hQGDf_;
104  this->tauQGD_ = this->aQGD_ * this->hQGD_ / cSound;
105 
106  forAll(p.primitiveField(), celli)
107  {
108  muQGD_.primitiveFieldRef()[celli] =
109  p.primitiveField()[celli] *
110  ScQGD_.primitiveField()[celli] *
111  tauQGD_.primitiveField()[celli];
112 
113  alphauQGD_.primitiveFieldRef()[celli] = muQGD_.primitiveField()[celli] /
114  PrQGD_.primitiveField()[celli];
115  }
116 
117  forAll(p.boundaryField(), patchi)
118  {
119  forAll(p.boundaryField()[patchi], facei)
120  {
121  muQGD_.boundaryFieldRef()[patchi][facei] =
122  p.boundaryField()[patchi][facei] *
123  ScQGD_.boundaryField()[patchi][facei] *
124  tauQGD_.boundaryField()[patchi][facei];
125 
126  alphauQGD_.boundaryFieldRef()[patchi][facei] =
127  muQGD_.boundaryFieldRef()[patchi][facei] /
128  PrQGD_.boundaryField()[patchi][facei];
129  }
130  }
131 }
132 
133 //
134 //END-OF-FILE
135 //
volScalarField PrQGD_
Definition: QGDCoeffs.H:115
virtual const volScalarField & c() const =0
Base class for all classes describing QGD model coefficients. Provides interfaces for accessing QGD/Q...
Definition: QGDCoeffs.H:57
void correct(const QGDThermo &)
constScPrModel1(const IOobject &io, const fvMesh &mesh, const dictionary &dict)
forAll(Y, i)
Definition: QGDYEqn.H:36
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
Class describing QGD model coefficients using constant values for Sc^{QGD}, Pr^{QGD} and specified f...
virtual const volScalarField & p() const =0
addToRunTimeSelectionTable(QGDCoeffs, constScPrModel1, dictionary)
defineTypeNameAndDebug(constScPrModel1, 0)
volScalarField & p
Definition: createFields.H:52
volScalarField ScQGD_
Definition: QGDCoeffs.H:120