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

SCOPEXiEq.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 "SCOPEXiEq.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace XiEqModels
00034 {
00035     defineTypeNameAndDebug(SCOPEXiEq, 0);
00036     addToRunTimeSelectionTable(XiEqModel, SCOPEXiEq, dictionary);
00037 };
00038 };
00039 
00040 
00041 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00042 
00043 Foam::XiEqModels::SCOPEXiEq::SCOPEXiEq
00044 (
00045     const dictionary& XiEqProperties,
00046     const hhuCombustionThermo& thermo,
00047     const compressible::RASModel& turbulence,
00048     const volScalarField& Su
00049 )
00050 :
00051     XiEqModel(XiEqProperties, thermo, turbulence, Su),
00052     XiEqCoef(readScalar(XiEqModelCoeffs_.lookup("XiEqCoef"))),
00053     XiEqExp(readScalar(XiEqModelCoeffs_.lookup("XiEqExp"))),
00054     lCoef(readScalar(XiEqModelCoeffs_.lookup("lCoef"))),
00055     SuMin(0.01*Su.average()),
00056     MaModel
00057     (
00058         IOdictionary
00059         (
00060             IOobject
00061             (
00062                 "combustionProperties",
00063                 Su.mesh().time().constant(),
00064                 Su.mesh(),
00065                 IOobject::MUST_READ
00066             )
00067         ),
00068         thermo
00069     )
00070 {}
00071 
00072 
00073 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
00074 
00075 Foam::XiEqModels::SCOPEXiEq::~SCOPEXiEq()
00076 {}
00077 
00078 
00079 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
00080 
00081 Foam::tmp<Foam::volScalarField> Foam::XiEqModels::SCOPEXiEq::XiEq() const
00082 {
00083     const volScalarField& k = turbulence_.k();
00084     const volScalarField& epsilon = turbulence_.epsilon();
00085 
00086     volScalarField up = sqrt((2.0/3.0)*k);
00087     volScalarField l = (lCoef*sqrt(3.0/2.0))*up*k/epsilon;
00088     volScalarField Rl = up*l*thermo_.rhou()/thermo_.muu();
00089 
00090     volScalarField upBySu = up/(Su_ + SuMin);
00091     volScalarField K = 0.157*upBySu/sqrt(Rl);
00092     volScalarField Ma = MaModel.Ma();
00093 
00094     tmp<volScalarField> tXiEq
00095     (
00096         new volScalarField
00097         (
00098             IOobject
00099             (
00100                 "XiEq",
00101                 epsilon.time().timeName(),
00102                 epsilon.db(),
00103                 IOobject::NO_READ,
00104                 IOobject::NO_WRITE
00105             ),
00106             epsilon.mesh(),
00107             dimensionedScalar("XiEq", dimless, 0.0)
00108         )
00109     );
00110     volScalarField& xieq = tXiEq();
00111 
00112     forAll(xieq, celli)
00113     {
00114         if (Ma[celli] > 0.01)
00115         {
00116             xieq[celli] =
00117                 XiEqCoef*pow(K[celli]*Ma[celli], -XiEqExp)*upBySu[celli];
00118         }
00119     }
00120 
00121     forAll(xieq.boundaryField(), patchi)
00122     {
00123         scalarField& xieqp = xieq.boundaryField()[patchi];
00124         const scalarField& Kp = K.boundaryField()[patchi];
00125         const scalarField& Map = Ma.boundaryField()[patchi];
00126         const scalarField& upBySup = upBySu.boundaryField()[patchi];
00127 
00128         forAll(xieqp, facei)
00129         {
00130             if (Ma[facei] > 0.01)
00131             {
00132                 xieqp[facei] =
00133                     XiEqCoef*pow(Kp[facei]*Map[facei], -XiEqExp)*upBySup[facei];
00134             }
00135         }
00136     }
00137 
00138     return tXiEq;
00139 }
00140 
00141 
00142 bool Foam::XiEqModels::SCOPEXiEq::read(const dictionary& XiEqProperties)
00143 {
00144     XiEqModel::read(XiEqProperties);
00145 
00146     XiEqModelCoeffs_.lookup("XiEqCoef") >> XiEqCoef;
00147     XiEqModelCoeffs_.lookup("XiEqExp") >> XiEqExp;
00148     XiEqModelCoeffs_.lookup("lCoef") >> lCoef;
00149 
00150     return true;
00151 }
00152 
00153 
00154 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines