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

turbulentMixingLengthDissipationRateInletFvPatchScalarField.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) 2006-2011 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 "turbulentMixingLengthDissipationRateInletFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/surfaceFields.H>
00030 #include <finiteVolume/volFields.H>
00031 #include <compressibleRASModels/RASModel.H>
00032 
00033 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00034 
00035 namespace Foam
00036 {
00037 namespace compressible
00038 {
00039 
00040 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00041 
00042 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
00043 turbulentMixingLengthDissipationRateInletFvPatchScalarField
00044 (
00045     const fvPatch& p,
00046     const DimensionedField<scalar, volMesh>& iF
00047 )
00048 :
00049     inletOutletFvPatchScalarField(p, iF),
00050     mixingLength_(0.0),
00051     phiName_("undefined-phi"),
00052     kName_("undefined-k")
00053 {
00054     this->refValue() = 0.0;
00055     this->refGrad() = 0.0;
00056     this->valueFraction() = 0.0;
00057 }
00058 
00059 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
00060 turbulentMixingLengthDissipationRateInletFvPatchScalarField
00061 (
00062     const turbulentMixingLengthDissipationRateInletFvPatchScalarField& ptf,
00063     const fvPatch& p,
00064     const DimensionedField<scalar, volMesh>& iF,
00065     const fvPatchFieldMapper& mapper
00066 )
00067 :
00068     inletOutletFvPatchScalarField(ptf, p, iF, mapper),
00069     mixingLength_(ptf.mixingLength_),
00070     phiName_(ptf.phiName_),
00071     kName_(ptf.kName_)
00072 {}
00073 
00074 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
00075 turbulentMixingLengthDissipationRateInletFvPatchScalarField
00076 (
00077     const fvPatch& p,
00078     const DimensionedField<scalar, volMesh>& iF,
00079     const dictionary& dict
00080 )
00081 :
00082     inletOutletFvPatchScalarField(p, iF),
00083     mixingLength_(readScalar(dict.lookup("mixingLength"))),
00084     phiName_(dict.lookupOrDefault<word>("phi", "phi")),
00085     kName_(dict.lookupOrDefault<word>("k", "k"))
00086 {
00087     fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
00088 
00089     this->refValue() = 0.0;
00090     this->refGrad() = 0.0;
00091     this->valueFraction() = 0.0;
00092 }
00093 
00094 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
00095 turbulentMixingLengthDissipationRateInletFvPatchScalarField
00096 (
00097     const turbulentMixingLengthDissipationRateInletFvPatchScalarField& ptf
00098 )
00099 :
00100     inletOutletFvPatchScalarField(ptf),
00101     mixingLength_(ptf.mixingLength_),
00102     phiName_(ptf.phiName_),
00103     kName_(ptf.kName_)
00104 {}
00105 
00106 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
00107 turbulentMixingLengthDissipationRateInletFvPatchScalarField
00108 (
00109     const turbulentMixingLengthDissipationRateInletFvPatchScalarField& ptf,
00110     const DimensionedField<scalar, volMesh>& iF
00111 )
00112 :
00113     inletOutletFvPatchScalarField(ptf, iF),
00114     mixingLength_(ptf.mixingLength_),
00115     phiName_(ptf.phiName_),
00116     kName_(ptf.kName_)
00117 {}
00118 
00119 
00120 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00121 
00122 void turbulentMixingLengthDissipationRateInletFvPatchScalarField::updateCoeffs()
00123 {
00124     if (updated())
00125     {
00126         return;
00127     }
00128 
00129     // Lookup Cmu corresponding to the turbulence model selected
00130     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
00131 
00132     const scalar Cmu =
00133         rasModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
00134 
00135     const scalar Cmu75 = pow(Cmu, 0.75);
00136 
00137     const fvPatchScalarField& kp =
00138         patch().lookupPatchField<volScalarField, scalar>(kName_);
00139 
00140     const fvsPatchScalarField& phip =
00141         patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
00142 
00143     this->refValue() = Cmu75*kp*sqrt(kp)/mixingLength_;
00144     this->valueFraction() = 1.0 - pos(phip);
00145 
00146     inletOutletFvPatchScalarField::updateCoeffs();
00147 }
00148 
00149 
00150 void turbulentMixingLengthDissipationRateInletFvPatchScalarField::write
00151 (
00152     Ostream& os
00153 ) const
00154 {
00155     fvPatchScalarField::write(os);
00156     os.writeKeyword("mixingLength")
00157         << mixingLength_ << token::END_STATEMENT << nl;
00158     os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
00159     os.writeKeyword("k") << kName_ << token::END_STATEMENT << nl;
00160     writeEntry("value", os);
00161 }
00162 
00163 
00164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00165 
00166 makePatchTypeField
00167 (
00168     fvPatchScalarField,
00169     turbulentMixingLengthDissipationRateInletFvPatchScalarField
00170 );
00171 
00172 
00173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00174 
00175 } // End namespace compressible
00176 } // End namespace Foam
00177 
00178 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines