FreeFOAM The Cross-Platform CFD Toolkit
turbulenceModel.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 "turbulenceModel.H"
27 #include <finiteVolume/volFields.H>
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 defineTypeNameAndDebug(turbulenceModel, 0);
40 defineRunTimeSelectionTable(turbulenceModel, turbulenceModel);
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 turbulenceModel::turbulenceModel
45 (
46  const volScalarField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& phi,
49  const basicThermo& thermophysicalModel
50 )
51 :
52  runTime_(U.time()),
53  mesh_(U.mesh()),
54 
55  rho_(rho),
56  U_(U),
57  phi_(phi),
58  thermophysicalModel_(thermophysicalModel)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
63 
64 autoPtr<turbulenceModel> turbulenceModel::New
65 (
66  const volScalarField& rho,
67  const volVectorField& U,
68  const surfaceScalarField& phi,
69  const basicThermo& thermophysicalModel
70 )
71 {
72  word modelName;
73 
74  // Enclose the creation of the dictionary to ensure it is deleted
75  // before the turbulenceModel is created otherwise the dictionary is
76  // entered in the database twice
77  {
78  IOdictionary dict
79  (
80  IOobject
81  (
82  "turbulenceProperties",
83  U.time().constant(),
84  U.db(),
87  )
88  );
89 
90  dict.lookup("simulationType") >> modelName;
91  }
92 
93  Info<< "Selecting turbulence model type " << modelName << endl;
94 
95  turbulenceModelConstructorTable::iterator cstrIter =
96  turbulenceModelConstructorTablePtr_->find(modelName);
97 
98  if (cstrIter == turbulenceModelConstructorTablePtr_->end())
99  {
101  (
102  "turbulenceModel::New(const volScalarField&, "
103  "const volVectorField&, const surfaceScalarField&, "
104  "basicThermo&)"
105  ) << "Unknown turbulenceModel type " << modelName
106  << endl << endl
107  << "Valid turbulenceModel types are :" << endl
108  << turbulenceModelConstructorTablePtr_->sortedToc()
109  << exit(FatalError);
110  }
111 
113  (
114  cstrIter()(rho, U, phi, thermophysicalModel)
115  );
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 
127 } // End namespace compressible
128 } // End namespace Foam
129 
130 // ************************ vim: set sw=4 sts=4 et: ************************ //