FreeFOAM The Cross-Platform CFD Toolkit
lowReOneEqEddy.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) 1991-2010 OpenCFD Ltd.
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 "lowReOneEqEddy.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace compressible
34 {
35 namespace LESModels
36 {
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(lowReOneEqEddy, 0);
41 addToRunTimeSelectionTable(LESModel, lowReOneEqEddy, dictionary);
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void lowReOneEqEddy::updateSubGridScaleFields()
46 {
47  // High Re eddy viscosity
48  muSgs_ = ck_*rho()*sqrt(k_)*delta();
49 
50  // low Re no corrected eddy viscosity
51  muSgs_ -= (mu()/beta_)*(scalar(1) - exp(-beta_*muSgs_/mu()));
52  muSgs_.correctBoundaryConditions();
53 
54  alphaSgs_ = muSgs_/Prt_;
55  alphaSgs_.correctBoundaryConditions();
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 lowReOneEqEddy::lowReOneEqEddy
62 (
63  const volScalarField& rho,
64  const volVectorField& U,
65  const surfaceScalarField& phi,
66  const basicThermo& thermoPhysicalModel
67 )
68 :
69  LESModel(typeName, rho, U, phi, thermoPhysicalModel),
70  GenEddyVisc(rho, U, phi, thermoPhysicalModel),
71 
72  ck_
73  (
75  (
76  "ck",
77  coeffDict_,
78  0.07
79  )
80  ),
81  beta_
82  (
84  (
85  "beta",
86  coeffDict_,
87  0.01
88  )
89  )
90 {
91  updateSubGridScaleFields();
92 
93  printCoeffs();
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 {
101  const volTensorField& gradU = tgradU();
102 
103  GenEddyVisc::correct(gradU);
104 
106  volScalarField G = 2*muSgs_*(gradU && dev(symm(gradU)));
107 
108  solve
109  (
110  fvm::ddt(rho(), k_)
111  + fvm::div(phi(), k_)
112  - fvm::laplacian(DkEff(), k_)
113  ==
114  G
115  - fvm::SuSp(2.0/3.0*rho()*divU, k_)
116  - fvm::Sp(ce_*rho()*sqrt(k_)/delta(), k_)
117  );
118 
119  bound(k_, k0());
120 
121  updateSubGridScaleFields();
122 }
123 
124 
125 bool lowReOneEqEddy::read()
126 {
127  if (GenEddyVisc::read())
128  {
129  ck_.readIfPresent(coeffDict());
130  beta_.readIfPresent(coeffDict());
131 
132  return true;
133  }
134  else
135  {
136  return false;
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 } // End namespace LESModels
144 } // End namespace compressible
145 } // End namespace Foam
146 
147 // ************************ vim: set sw=4 sts=4 et: ************************ //