EddyDiffusivity.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 | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "EddyDiffusivity.H"
27 
28 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
29 
30 template<class BasicTurbulenceModel>
32 {
33  // Read Prt if provided
35  (
36  "Prt",
37  this->coeffDict(),
38  1.0
39  );
40 
41  alphat_ = this->rho_*this->nut()/Prt_;
42  alphat_.correctBoundaryConditions();
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<class BasicTurbulenceModel>
50 (
51  const word& type,
52  const alphaField& alpha,
53  const volScalarField& rho,
54  const volVectorField& U,
55  const surfaceScalarField& alphaRhoPhi,
56  const surfaceScalarField& phi,
57  const transportModel& transport,
58  const word& propertiesName
59 )
60 :
61  BasicTurbulenceModel
62  (
63  type,
64  alpha,
65  rho,
66  U,
67  alphaRhoPhi,
68  phi,
69  transport,
70  propertiesName
71  ),
72 
73  // Cannot read Prt yet
74  Prt_("Prt", dimless, 1.0),
75 
76  alphat_
77  (
78  IOobject
79  (
80  IOobject::groupName("alphat", U.group()),
81  this->runTime_.timeName(),
82  this->mesh_,
83  IOobject::MUST_READ,
84  IOobject::AUTO_WRITE
85  ),
86  this->mesh_
87  )
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 template<class BasicTurbulenceModel>
95 {
97  {
98  Prt_.readIfPresent(this->coeffDict());
99 
100  return true;
101  }
102  else
103  {
104  return false;
105  }
106 }
107 
108 
109 template<class BasicTurbulenceModel>
111 {
113 }
114 
115 
116 // ************************************************************************* //
virtual void correctNut()
EddyDiffusivity(const word &type, const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &trasport, const word &propertiesName)
Construct.
word group() const
Return group (extension part of name)
Definition: IOobject.C:349
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
BasicTurbulenceModel::alphaField alphaField
Templated abstract base class for single-phase compressible turbulence models.
A class for handling words, derived from string.
Definition: word.H:59
virtual void correctEnergyTransport()
Correct the turbulence thermal diffusivity for energy transport.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
scalar nut