Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
turbulenceModels
incompressible
RAS
laminar
laminar.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 "
laminar.H
"
27
#include <
OpenFOAM/addToRunTimeSelectionTable.H
>
28
29
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30
31
namespace
Foam
32
{
33
namespace
incompressible
34
{
35
namespace
RASModels
36
{
37
38
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40
defineTypeNameAndDebug
(laminar, 0);
41
addToRunTimeSelectionTable
(
RASModel
, laminar, dictionary);
42
43
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45
laminar::laminar
46
(
47
const
volVectorField
&
U
,
48
const
surfaceScalarField
&
phi
,
49
transportModel
& lamTransportModel
50
)
51
:
52
RASModel
(typeName, U, phi, lamTransportModel)
53
{}
54
55
56
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57
58
tmp<volScalarField>
laminar::nut
()
const
59
{
60
return
tmp<volScalarField>
61
(
62
new
volScalarField
63
(
64
IOobject
65
(
66
"nut"
,
67
runTime_
.
timeName
(),
68
mesh_
,
69
IOobject::NO_READ
,
70
IOobject::NO_WRITE
71
),
72
mesh_
,
73
dimensionedScalar
(
"nut"
,
nu
().
dimensions
(), 0.0)
74
)
75
);
76
}
77
78
79
tmp<volScalarField>
laminar::k
()
const
80
{
81
return
tmp<volScalarField>
82
(
83
new
volScalarField
84
(
85
IOobject
86
(
87
"k"
,
88
runTime_
.
timeName
(),
89
mesh_
,
90
IOobject::NO_READ
,
91
IOobject::NO_WRITE
92
),
93
mesh_
,
94
dimensionedScalar
(
"k"
,
sqr
(
U_
.
dimensions
()), 0.0)
95
)
96
);
97
}
98
99
100
tmp<volScalarField>
laminar::epsilon
()
const
101
{
102
return
tmp<volScalarField>
103
(
104
new
volScalarField
105
(
106
IOobject
107
(
108
"epsilon"
,
109
runTime_
.
timeName
(),
110
mesh_
,
111
IOobject::NO_READ
,
112
IOobject::NO_WRITE
113
),
114
mesh_
,
115
dimensionedScalar
116
(
117
"epsilon"
,
sqr
(
U_
.
dimensions
())/
dimTime
, 0.0
118
)
119
)
120
);
121
}
122
123
124
tmp<volSymmTensorField>
laminar::R
()
const
125
{
126
return
tmp<volSymmTensorField>
127
(
128
new
volSymmTensorField
129
(
130
IOobject
131
(
132
"R"
,
133
runTime_
.
timeName
(),
134
mesh_
,
135
IOobject::NO_READ
,
136
IOobject::NO_WRITE
137
),
138
mesh_
,
139
dimensionedSymmTensor
140
(
141
"R"
,
sqr
(
U_
.
dimensions
()),
symmTensor::zero
142
)
143
)
144
);
145
}
146
147
148
tmp<volSymmTensorField>
laminar::devReff
()
const
149
{
150
return
tmp<volSymmTensorField>
151
(
152
new
volSymmTensorField
153
(
154
IOobject
155
(
156
"devRhoReff"
,
157
runTime_
.
timeName
(),
158
mesh_
,
159
IOobject::NO_READ
,
160
IOobject::NO_WRITE
161
),
162
-
nu
()*
dev
(
twoSymm
(
fvc::grad
(
U_
)))
163
)
164
);
165
}
166
167
168
tmp<fvVectorMatrix>
laminar::divDevReff
(
volVectorField
& U)
const
169
{
170
return
171
(
172
-
fvm::laplacian
(
nuEff
(), U)
173
-
fvc::div
(
nuEff
()*
dev
(
fvc::grad
(U)().
T
()))
174
);
175
}
176
177
178
bool
laminar::read
()
179
{
180
return
RASModel::read
();
181
}
182
183
184
void
laminar::correct
()
185
{
186
turbulenceModel::correct
();
187
}
188
189
190
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191
192
}
// End namespace RASModels
193
}
// End namespace incompressible
194
}
// End namespace Foam
195
196
// ************************ vim: set sw=4 sts=4 et: ************************ //