FreeFOAM The Cross-Platform CFD Toolkit
basicXiSubXiEq.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 "basicXiSubXiEq.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiEqModels
34 {
35  defineTypeNameAndDebug(basicSubGrid, 0);
36  addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary);
37 };
38 };
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiEqModels::basicSubGrid::basicSubGrid
44 (
45  const dictionary& XiEqProperties,
46  const hhuCombustionThermo& thermo,
48  const volScalarField& Su
49 )
50 :
51  XiEqModel(XiEqProperties, thermo, turbulence, Su),
52 
53  N_
54  (
55  IOobject
56  (
57  "N",
58  Su.mesh().time().findInstance(polyMesh::meshSubDir, "N"),
59  polyMesh::meshSubDir,
60  Su.mesh(),
61  IOobject::MUST_READ,
62  IOobject::NO_WRITE
63  ),
64  Su.mesh()
65  ),
66 
67  ns_
68  (
69  IOobject
70  (
71  "ns",
72  Su.mesh().time().findInstance(polyMesh::meshSubDir, "ns"),
73  polyMesh::meshSubDir,
74  Su.mesh(),
75  IOobject::MUST_READ,
76  IOobject::NO_WRITE
77  ),
78  Su.mesh()
79  ),
80 
81  B_
82  (
83  IOobject
84  (
85  "B",
86  Su.mesh().time().findInstance(polyMesh::meshSubDir, "B"),
87  polyMesh::meshSubDir,
88  Su.mesh(),
89  IOobject::MUST_READ,
90  IOobject::NO_WRITE
91  ),
92  Su.mesh()
93  ),
94 
95  Lobs_
96  (
97  IOobject
98  (
99  "Lobs",
100  Su.mesh().time().findInstance(polyMesh::meshSubDir, "Lobs"),
101  polyMesh::meshSubDir,
102  Su.mesh(),
103  IOobject::MUST_READ,
104  IOobject::NO_WRITE
105  ),
106  Su.mesh()
107  ),
108 
109  XiEqModel_(XiEqModel::New(XiEqModelCoeffs_, thermo, turbulence, Su))
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
114 
116 {}
117 
118 
119 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
120 
122 {
123  const objectRegistry& db = Su_.db();
124  const volVectorField& U = db.lookupObject<volVectorField>("U");
125 
126  volScalarField magU = mag(U);
127  volVectorField Uhat =
128  U/(mag(U) + dimensionedScalar("Usmall", U.dimensions(), 1e-4));
129 
130  volScalarField n = max(N_ - (Uhat & ns_ & Uhat), scalar(1e-4));
131 
132  volScalarField b = (Uhat & B_ & Uhat)/n;
133 
134  volScalarField up = sqrt((2.0/3.0)*turbulence_.k());
135 
136  volScalarField XiSubEq =
137  scalar(1)
138  + max(2.2*sqrt(b), min(0.34*magU/up, scalar(1.6)))
139  *min(0.25*n, scalar(1));
140 
141  return XiSubEq*XiEqModel_->XiEq();
142 }
143 
144 
145 bool Foam::XiEqModels::basicSubGrid::read(const dictionary& XiEqProperties)
146 {
147  XiEqModel::read(XiEqProperties);
148 
149  return XiEqModel_->read(XiEqModelCoeffs_);
150 }
151 
152 
153 // ************************ vim: set sw=4 sts=4 et: ************************ //