FreeFOAM The Cross-Platform CFD Toolkit
oneEqEddy.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 "oneEqEddy.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace incompressible
34 {
35 namespace LESModels
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 defineTypeNameAndDebug(oneEqEddy, 0);
41 addToRunTimeSelectionTable(LESModel, oneEqEddy, dictionary);
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void oneEqEddy::updateSubGridScaleFields()
47 {
48  nuSgs_ = ck_*sqrt(k_)*delta();
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 oneEqEddy::oneEqEddy
56 (
57  const volVectorField& U,
58  const surfaceScalarField& phi,
59  transportModel& transport
60 )
61 :
62  LESModel(typeName, U, phi, transport),
63  GenEddyVisc(U, phi, transport),
64 
65  k_
66  (
67  IOobject
68  (
69  "k",
70  runTime_.timeName(),
71  mesh_,
74  ),
75  mesh_
76  ),
77 
78  ck_
79  (
81  (
82  "ck",
83  coeffDict_,
84  0.094
85  )
86  )
87 {
88  updateSubGridScaleFields();
89 
90  printCoeffs();
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 {
98  GenEddyVisc::correct(gradU);
99 
100  volScalarField G = 2.0*nuSgs_*magSqr(symm(gradU));
101 
102  fvScalarMatrix kEqn
103  (
104  fvm::ddt(k_)
105  + fvm::div(phi(), k_)
106  - fvm::laplacian(DkEff(), k_)
107  ==
108  G
109  - fvm::Sp(ce_*sqrt(k_)/delta(), k_)
110  );
111 
112  kEqn.relax();
113  kEqn.solve();
114 
115  bound(k_, k0());
116 
117  updateSubGridScaleFields();
118 }
119 
120 
122 {
123  if (GenEddyVisc::read())
124  {
125  ck_.readIfPresent(coeffDict());
126 
127  return true;
128  }
129  else
130  {
131  return false;
132  }
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace LESModels
139 } // End namespace incompressible
140 } // End namespace Foam
141 
142 // ************************ vim: set sw=4 sts=4 et: ************************ //