FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

laminar.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "laminar.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <finiteVolume/fvcGrad.H>
00030 #include <finiteVolume/fvcDiv.H>
00031 #include <finiteVolume/fvmLaplacian.H>
00032 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00033 
00034 
00035 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00036 
00037 namespace Foam
00038 {
00039 namespace incompressible
00040 {
00041 
00042 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00043 
00044 defineTypeNameAndDebug(laminar, 0);
00045 addToRunTimeSelectionTable(turbulenceModel, laminar, turbulenceModel);
00046 
00047 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00048 
00049 laminar::laminar
00050 (
00051     const volVectorField& U,
00052     const surfaceScalarField& phi,
00053     transportModel& lamTransportModel
00054 )
00055 :
00056     turbulenceModel(U, phi, lamTransportModel)
00057 {}
00058 
00059 
00060 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
00061 
00062 autoPtr<laminar> laminar::New
00063 (
00064     const volVectorField& U,
00065     const surfaceScalarField& phi,
00066     transportModel& lamTransportModel
00067 )
00068 {
00069     return autoPtr<laminar>(new laminar(U, phi, lamTransportModel));
00070 }
00071 
00072 
00073 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00074 
00075 tmp<volScalarField> laminar::nut() const
00076 {
00077     return tmp<volScalarField>
00078     (
00079         new volScalarField
00080         (
00081             IOobject
00082             (
00083                 "nut",
00084                 runTime_.timeName(),
00085                 mesh_,
00086                 IOobject::NO_READ,
00087                 IOobject::NO_WRITE
00088             ),
00089             mesh_,
00090             dimensionedScalar("nut", nu().dimensions(), 0.0)
00091         )
00092     );
00093 }
00094 
00095 
00096 tmp<volScalarField> laminar::nuEff() const
00097 {
00098     return tmp<volScalarField>(new volScalarField("nuEff", nu()));
00099 }
00100 
00101 
00102 tmp<volScalarField> laminar::k() const
00103 {
00104     return tmp<volScalarField>
00105     (
00106         new volScalarField
00107         (
00108             IOobject
00109             (
00110                 "k",
00111                 runTime_.timeName(),
00112                 mesh_,
00113                 IOobject::NO_READ,
00114                 IOobject::NO_WRITE
00115             ),
00116             mesh_,
00117             dimensionedScalar("k", sqr(U_.dimensions()), 0.0)
00118         )
00119     );
00120 }
00121 
00122 
00123 tmp<volScalarField> laminar::epsilon() const
00124 {
00125     return tmp<volScalarField>
00126     (
00127         new volScalarField
00128         (
00129             IOobject
00130             (
00131                 "epsilon",
00132                 runTime_.timeName(),
00133                 mesh_,
00134                 IOobject::NO_READ,
00135                 IOobject::NO_WRITE
00136             ),
00137             mesh_,
00138             dimensionedScalar
00139             (
00140                 "epsilon", sqr(U_.dimensions())/dimTime, 0.0
00141             )
00142         )
00143     );
00144 }
00145 
00146 
00147 tmp<volSymmTensorField> laminar::R() const
00148 {
00149     return tmp<volSymmTensorField>
00150     (
00151         new volSymmTensorField
00152         (
00153             IOobject
00154             (
00155                 "R",
00156                 runTime_.timeName(),
00157                 mesh_,
00158                 IOobject::NO_READ,
00159                 IOobject::NO_WRITE
00160             ),
00161             mesh_,
00162             dimensionedSymmTensor
00163             (
00164                 "R", sqr(U_.dimensions()), symmTensor::zero
00165             )
00166         )
00167     );
00168 }
00169 
00170 
00171 tmp<volSymmTensorField> laminar::devReff() const
00172 {
00173     return tmp<volSymmTensorField>
00174     (
00175         new volSymmTensorField
00176         (
00177             IOobject
00178             (
00179                 "devRhoReff",
00180                 runTime_.timeName(),
00181                 mesh_,
00182                 IOobject::NO_READ,
00183                 IOobject::NO_WRITE
00184             ),
00185            -nu()*dev(twoSymm(fvc::grad(U_)))
00186         )
00187     );
00188 }
00189 
00190 
00191 tmp<fvVectorMatrix> laminar::divDevReff(volVectorField& U) const
00192 {
00193     return
00194     (
00195       - fvm::laplacian(nuEff(), U)
00196       - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00197     );
00198 }
00199 
00200 
00201 bool laminar::read()
00202 {
00203     return true;
00204 }
00205 
00206 
00207 void laminar::correct()
00208 {
00209     turbulenceModel::correct();
00210 }
00211 
00212 
00213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00214 
00215 } // End namespace incompressible
00216 } // End namespace Foam
00217 
00218 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines