BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / turbulenceModel / laminar / laminar.C
blobc102b52cb1e10e2ae7a5525b590be169360b9a1c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "laminar.H"
27 #include "Time.H"
28 #include "volFields.H"
29 #include "fvcGrad.H"
30 #include "fvcDiv.H"
31 #include "fvmLaplacian.H"
32 #include "addToRunTimeSelectionTable.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace Foam
39 namespace incompressible
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 defineTypeNameAndDebug(laminar, 0);
45 addToRunTimeSelectionTable(turbulenceModel, laminar, turbulenceModel);
47 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
49 laminar::laminar
51     const volVectorField& U,
52     const surfaceScalarField& phi,
53     transportModel& transport,
54     const word& turbulenceModelName
57     turbulenceModel(U, phi, transport, turbulenceModelName)
61 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
63 autoPtr<laminar> laminar::New
65     const volVectorField& U,
66     const surfaceScalarField& phi,
67     transportModel& transport,
68     const word& turbulenceModelName
71     return autoPtr<laminar>
72     (
73         new laminar(U, phi, transport, turbulenceModelName)
74     );
78 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
80 tmp<volScalarField> laminar::nut() const
82     return tmp<volScalarField>
83     (
84         new volScalarField
85         (
86             IOobject
87             (
88                 "nut",
89                 runTime_.timeName(),
90                 mesh_,
91                 IOobject::NO_READ,
92                 IOobject::NO_WRITE
93             ),
94             mesh_,
95             dimensionedScalar("nut", nu()().dimensions(), 0.0)
96         )
97     );
101 tmp<volScalarField> laminar::nuEff() const
103     return tmp<volScalarField>(new volScalarField("nuEff", nu()));
107 tmp<volScalarField> laminar::k() const
109     return tmp<volScalarField>
110     (
111         new volScalarField
112         (
113             IOobject
114             (
115                 "k",
116                 runTime_.timeName(),
117                 mesh_,
118                 IOobject::NO_READ,
119                 IOobject::NO_WRITE
120             ),
121             mesh_,
122             dimensionedScalar("k", sqr(U_.dimensions()), 0.0)
123         )
124     );
128 tmp<volScalarField> laminar::epsilon() const
130     return tmp<volScalarField>
131     (
132         new volScalarField
133         (
134             IOobject
135             (
136                 "epsilon",
137                 runTime_.timeName(),
138                 mesh_,
139                 IOobject::NO_READ,
140                 IOobject::NO_WRITE
141             ),
142             mesh_,
143             dimensionedScalar
144             (
145                 "epsilon", sqr(U_.dimensions())/dimTime, 0.0
146             )
147         )
148     );
152 tmp<volSymmTensorField> laminar::R() const
154     return tmp<volSymmTensorField>
155     (
156         new volSymmTensorField
157         (
158             IOobject
159             (
160                 "R",
161                 runTime_.timeName(),
162                 mesh_,
163                 IOobject::NO_READ,
164                 IOobject::NO_WRITE
165             ),
166             mesh_,
167             dimensionedSymmTensor
168             (
169                 "R", sqr(U_.dimensions()), symmTensor::zero
170             )
171         )
172     );
176 tmp<volSymmTensorField> laminar::devReff() const
178     return tmp<volSymmTensorField>
179     (
180         new volSymmTensorField
181         (
182             IOobject
183             (
184                 "devRhoReff",
185                 runTime_.timeName(),
186                 mesh_,
187                 IOobject::NO_READ,
188                 IOobject::NO_WRITE
189             ),
190            -nu()*dev(twoSymm(fvc::grad(U_)))
191         )
192     );
196 tmp<fvVectorMatrix> laminar::divDevReff(volVectorField& U) const
198     return
199     (
200       - fvm::laplacian(nuEff(), U)
201       - fvc::div(nuEff()*dev(T(fvc::grad(U))))
202     );
206 void laminar::correct()
208     turbulenceModel::correct();
212 bool laminar::read()
214     return true;
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 } // End namespace incompressible
221 } // End namespace Foam
223 // ************************************************************************* //