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 compressible
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 volScalarField& rho,
00052     const volVectorField& U,
00053     const surfaceScalarField& phi,
00054     const basicThermo& thermophysicalModel
00055 )
00056 :
00057     turbulenceModel(rho, U, phi, thermophysicalModel)
00058 {}
00059 
00060 
00061 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
00062 
00063 autoPtr<laminar> laminar::New
00064 (
00065     const volScalarField& rho,
00066     const volVectorField& U,
00067     const surfaceScalarField& phi,
00068     const basicThermo& thermophysicalModel
00069 )
00070 {
00071     return autoPtr<laminar>(new laminar(rho, U, phi, thermophysicalModel));
00072 }
00073 
00074 
00075 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00076 
00077 tmp<volScalarField> laminar::mut() const
00078 {
00079     return tmp<volScalarField>
00080     (
00081         new volScalarField
00082         (
00083             IOobject
00084             (
00085                 "mut",
00086                 runTime_.timeName(),
00087                 mesh_,
00088                 IOobject::NO_READ,
00089                 IOobject::NO_WRITE
00090             ),
00091             mesh_,
00092             dimensionedScalar("mut", mu().dimensions(), 0.0)
00093         )
00094     );
00095 }
00096 
00097 
00098 tmp<volScalarField> laminar::k() const
00099 {
00100     return tmp<volScalarField>
00101     (
00102         new volScalarField
00103         (
00104             IOobject
00105             (
00106                 "k",
00107                 runTime_.timeName(),
00108                 mesh_,
00109                 IOobject::NO_READ,
00110                 IOobject::NO_WRITE
00111             ),
00112             mesh_,
00113             dimensionedScalar("k", sqr(U_.dimensions()), 0.0)
00114         )
00115     );
00116 }
00117 
00118 
00119 tmp<volScalarField> laminar::epsilon() const
00120 {
00121     return tmp<volScalarField>
00122     (
00123         new volScalarField
00124         (
00125             IOobject
00126             (
00127                 "epsilon",
00128                 runTime_.timeName(),
00129                 mesh_,
00130                 IOobject::NO_READ,
00131                 IOobject::NO_WRITE
00132             ),
00133             mesh_,
00134             dimensionedScalar
00135             (
00136                 "epsilon", sqr(U_.dimensions())/dimTime, 0.0
00137             )
00138         )
00139     );
00140 }
00141 
00142 
00143 tmp<volSymmTensorField> laminar::R() const
00144 {
00145     return tmp<volSymmTensorField>
00146     (
00147         new volSymmTensorField
00148         (
00149             IOobject
00150             (
00151                 "R",
00152                 runTime_.timeName(),
00153                 mesh_,
00154                 IOobject::NO_READ,
00155                 IOobject::NO_WRITE
00156             ),
00157             mesh_,
00158             dimensionedSymmTensor
00159             (
00160                 "R", sqr(U_.dimensions()), symmTensor::zero
00161             )
00162         )
00163     );
00164 }
00165 
00166 
00167 tmp<volSymmTensorField> laminar::devRhoReff() const
00168 {
00169     return tmp<volSymmTensorField>
00170     (
00171         new volSymmTensorField
00172         (
00173             IOobject
00174             (
00175                 "devRhoReff",
00176                 runTime_.timeName(),
00177                 mesh_,
00178                 IOobject::NO_READ,
00179                 IOobject::NO_WRITE
00180             ),
00181            -mu()*dev(twoSymm(fvc::grad(U_)))
00182         )
00183     );
00184 }
00185 
00186 
00187 tmp<fvVectorMatrix> laminar::divDevRhoReff(volVectorField& U) const
00188 {
00189     return
00190     (
00191       - fvm::laplacian(muEff(), U)
00192       - fvc::div(muEff()*dev2(fvc::grad(U)().T()))
00193     );
00194 }
00195 
00196 
00197 bool laminar::read()
00198 {
00199     return true;
00200 }
00201 
00202 
00203 void laminar::correct()
00204 {
00205     turbulenceModel::correct();
00206 }
00207 
00208 
00209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00210 
00211 } // End namespace incompressible
00212 } // End namespace Foam
00213 
00214 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines