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

maxhxhyhzDelta.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) 2010-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 "maxhxhyhzDelta.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 
00034 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00035 
00036 defineTypeNameAndDebug(maxhxhyhzDelta, 0);
00037 addToRunTimeSelectionTable(LESdelta, maxhxhyhzDelta, dictionary);
00038 
00039 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00040 
00041 void maxhxhyhzDelta::calcDelta()
00042 {
00043     label nD = mesh().nGeometricD();
00044 
00045     tmp<volScalarField> hmax
00046     (
00047         new volScalarField
00048         (
00049             IOobject
00050             (
00051                 "hmax",
00052                 mesh().time().timeName(),
00053                 mesh(),
00054                 IOobject::NO_READ,
00055                 IOobject::NO_WRITE
00056             ),
00057             mesh(),
00058             dimensionedScalar("zrero", dimLength, 0.0)
00059         )
00060     );
00061 
00062     const cellList& cells = mesh().cells();
00063 
00064     forAll(cells,cellI)
00065     {
00066         scalar deltaMaxTmp = 0.0;
00067         const labelList& cFaces = mesh().cells()[cellI];
00068         const point& centrevector = mesh().cellCentres()[cellI];
00069 
00070         forAll(cFaces, cFaceI)
00071         {
00072             label faceI = cFaces[cFaceI];
00073             const point& facevector = mesh().faceCentres()[faceI];
00074             scalar tmp = mag(facevector - centrevector);
00075             if (tmp > deltaMaxTmp)
00076             {
00077                 deltaMaxTmp = tmp;
00078             }
00079         }
00080         hmax()[cellI] = deltaCoeff_*deltaMaxTmp;
00081     }
00082 
00083     if (nD == 3)
00084     {
00085         delta_.internalField() = hmax();
00086     }
00087     else if (nD == 2)
00088     {
00089         WarningIn("maxhxhyhzDelta::calcDelta()")
00090             << "Case is 2D, LES is not strictly applicable\n"
00091             << endl;
00092 
00093         delta_.internalField() = hmax();
00094     }
00095     else
00096     {
00097         FatalErrorIn("maxhxhyhzDelta::calcDelta()")
00098             << "Case is not 3D or 2D, LES is not applicable"
00099             << exit(FatalError);
00100     }
00101 }
00102 
00103 
00104 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00105 
00106 maxhxhyhzDelta::maxhxhyhzDelta
00107 (
00108     const word& name,
00109     const fvMesh& mesh,
00110     const dictionary& dd
00111 )
00112 :
00113     LESdelta(name, mesh),
00114     deltaCoeff_(readScalar(dd.subDict(type() + "Coeffs").lookup("deltaCoeff")))
00115 {
00116     calcDelta();
00117 }
00118 
00119 
00120 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00121 
00122 void maxhxhyhzDelta::read(const dictionary& dd)
00123 {
00124     dd.subDict(type() + "Coeffs").lookup("deltaCoeff") >> deltaCoeff_;
00125     calcDelta();
00126 }
00127 
00128 
00129 void maxhxhyhzDelta::correct()
00130 {
00131     if (mesh_.changing())
00132     {
00133         calcDelta();
00134     }
00135 }
00136 
00137 
00138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00139 
00140 } // End namespace Foam
00141 
00142 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines