initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / compressible / turbulenceModel / laminar / laminar.C
blob4a391c3db5b1d2f4d509a898ab4308e5d2a6cc31
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "laminar.H"
28 #include "Time.H"
29 #include "volFields.H"
30 #include "fvcGrad.H"
31 #include "fvcDiv.H"
32 #include "fvmLaplacian.H"
33 #include "addToRunTimeSelectionTable.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 namespace Foam
40 namespace compressible
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 defineTypeNameAndDebug(laminar, 0);
46 addToRunTimeSelectionTable(turbulenceModel, laminar, turbulenceModel);
48 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
50 laminar::laminar
52     const volScalarField& rho,
53     const volVectorField& U,
54     const surfaceScalarField& phi,
55     const basicThermo& thermophysicalModel
58     turbulenceModel(rho, U, phi, thermophysicalModel)
62 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
64 autoPtr<laminar> laminar::New
66     const volScalarField& rho,
67     const volVectorField& U,
68     const surfaceScalarField& phi,
69     const basicThermo& thermophysicalModel
72     return autoPtr<laminar>(new laminar(rho, U, phi, thermophysicalModel));
76 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
78 tmp<volScalarField> laminar::mut() const
80     return tmp<volScalarField>
81     (
82         new volScalarField
83         (
84             IOobject
85             (
86                 "mut",
87                 runTime_.timeName(),
88                 mesh_,
89                 IOobject::NO_READ,
90                 IOobject::NO_WRITE
91             ),
92             mesh_,
93             dimensionedScalar("mut", mu().dimensions(), 0.0)
94         )
95     );
99 tmp<volScalarField> laminar::k() const
101     return tmp<volScalarField>
102     (
103         new volScalarField
104         (
105             IOobject
106             (
107                 "k",
108                 runTime_.timeName(),
109                 mesh_,
110                 IOobject::NO_READ,
111                 IOobject::NO_WRITE
112             ),
113             mesh_,
114             dimensionedScalar("k", sqr(U_.dimensions()), 0.0)
115         )
116     );
120 tmp<volScalarField> laminar::epsilon() const
122     return tmp<volScalarField>
123     (
124         new volScalarField
125         (
126             IOobject
127             (
128                 "epsilon",
129                 runTime_.timeName(),
130                 mesh_,
131                 IOobject::NO_READ,
132                 IOobject::NO_WRITE
133             ),
134             mesh_,
135             dimensionedScalar
136             (
137                 "epsilon", sqr(U_.dimensions())/dimTime, 0.0
138             )
139         )
140     );
144 tmp<volSymmTensorField> laminar::R() const
146     return tmp<volSymmTensorField>
147     (
148         new volSymmTensorField
149         (
150             IOobject
151             (
152                 "R",
153                 runTime_.timeName(),
154                 mesh_,
155                 IOobject::NO_READ,
156                 IOobject::NO_WRITE
157             ),
158             mesh_,
159             dimensionedSymmTensor
160             (
161                 "R", sqr(U_.dimensions()), symmTensor::zero
162             )
163         )
164     );
168 tmp<volSymmTensorField> laminar::devRhoReff() const
170     return tmp<volSymmTensorField>
171     (
172         new volSymmTensorField
173         (
174             IOobject
175             (
176                 "devRhoReff",
177                 runTime_.timeName(),
178                 mesh_,
179                 IOobject::NO_READ,
180                 IOobject::NO_WRITE
181             ),
182            -mu()*dev(twoSymm(fvc::grad(U_)))
183         )
184     );
188 tmp<fvVectorMatrix> laminar::divDevRhoReff(volVectorField& U) const
190     return
191     (
192       - fvm::laplacian(muEff(), U)
193       - fvc::div(muEff()*dev2(fvc::grad(U)().T()))
194     );
198 bool laminar::read()
200     return true;
204 void laminar::correct()
206     turbulenceModel::correct();
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 } // End namespace incompressible
213 } // End namespace Foam
215 // ************************************************************************* //