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

turbulentInletFvPatchField.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 <finiteVolume/turbulentInletFvPatchField.H>
00027 
00028 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00029 
00030 namespace Foam
00031 {
00032 
00033 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00034 
00035 template<class Type>
00036 turbulentInletFvPatchField<Type>::turbulentInletFvPatchField
00037 (
00038     const fvPatch& p,
00039     const DimensionedField<Type, volMesh>& iF
00040 )
00041 :
00042     fixedValueFvPatchField<Type>(p, iF),
00043     ranGen_(label(0)),
00044     fluctuationScale_(pTraits<Type>::zero),
00045     referenceField_(p.size()),
00046     alpha_(0.1),
00047     curTimeIndex_(-1)
00048 {}
00049 
00050 
00051 template<class Type>
00052 turbulentInletFvPatchField<Type>::turbulentInletFvPatchField
00053 (
00054     const turbulentInletFvPatchField<Type>& ptf,
00055     const fvPatch& p,
00056     const DimensionedField<Type, volMesh>& iF,
00057     const fvPatchFieldMapper& mapper
00058 )
00059 :
00060     fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
00061     ranGen_(label(0)),
00062     fluctuationScale_(ptf.fluctuationScale_),
00063     referenceField_(ptf.referenceField_, mapper),
00064     alpha_(ptf.alpha_),
00065     curTimeIndex_(-1)
00066 {}
00067 
00068 
00069 template<class Type>
00070 turbulentInletFvPatchField<Type>::turbulentInletFvPatchField
00071 (
00072     const fvPatch& p,
00073     const DimensionedField<Type, volMesh>& iF,
00074     const dictionary& dict
00075 )
00076 :
00077     fixedValueFvPatchField<Type>(p, iF),
00078     ranGen_(label(0)),
00079     fluctuationScale_(pTraits<Type>(dict.lookup("fluctuationScale"))),
00080     referenceField_("referenceField", dict, p.size()),
00081     alpha_(dict.lookupOrDefault<scalar>("alpha", 0.1)),
00082     curTimeIndex_(-1)
00083 {
00084     if (dict.found("value"))
00085     {
00086         fixedValueFvPatchField<Type>::operator==
00087         (
00088             Field<Type>("value", dict, p.size())
00089         );
00090     }
00091     else
00092     {
00093         fixedValueFvPatchField<Type>::operator==(referenceField_);
00094     }
00095 }
00096 
00097 
00098 template<class Type>
00099 turbulentInletFvPatchField<Type>::turbulentInletFvPatchField
00100 (
00101     const turbulentInletFvPatchField<Type>& ptf
00102 )
00103 :
00104     fixedValueFvPatchField<Type>(ptf),
00105     ranGen_(ptf.ranGen_),
00106     fluctuationScale_(ptf.fluctuationScale_),
00107     referenceField_(ptf.referenceField_),
00108     alpha_(ptf.alpha_),
00109     curTimeIndex_(-1)
00110 {}
00111 
00112 
00113 template<class Type>
00114 turbulentInletFvPatchField<Type>::turbulentInletFvPatchField
00115 (
00116     const turbulentInletFvPatchField<Type>& ptf,
00117     const DimensionedField<Type, volMesh>& iF
00118 )
00119 :
00120     fixedValueFvPatchField<Type>(ptf, iF),
00121     ranGen_(ptf.ranGen_),
00122     fluctuationScale_(ptf.fluctuationScale_),
00123     referenceField_(ptf.referenceField_),
00124     alpha_(ptf.alpha_),
00125     curTimeIndex_(-1)
00126 {}
00127 
00128 
00129 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00130 
00131 template<class Type>
00132 void turbulentInletFvPatchField<Type>::autoMap
00133 (
00134     const fvPatchFieldMapper& m
00135 )
00136 {
00137     fixedValueFvPatchField<Type>::autoMap(m);
00138     referenceField_.autoMap(m);
00139 }
00140 
00141 
00142 template<class Type>
00143 void turbulentInletFvPatchField<Type>::rmap
00144 (
00145     const fvPatchField<Type>& ptf,
00146     const labelList& addr
00147 )
00148 {
00149     fixedValueFvPatchField<Type>::rmap(ptf, addr);
00150 
00151     const turbulentInletFvPatchField<Type>& tiptf =
00152         refCast<const turbulentInletFvPatchField<Type> >(ptf);
00153 
00154     referenceField_.rmap(tiptf.referenceField_, addr);
00155 }
00156 
00157 
00158 template<class Type>
00159 void turbulentInletFvPatchField<Type>::updateCoeffs()
00160 {
00161     if (this->updated())
00162     {
00163         return;
00164     }
00165 
00166     if (curTimeIndex_ != this->db().time().timeIndex())
00167     {
00168         Field<Type>& patchField = *this;
00169 
00170         Field<Type> randomField(this->size());
00171 
00172         forAll(patchField, facei)
00173         {
00174             ranGen_.randomise(randomField[facei]);
00175         }
00176 
00177         // Correction-factor to compensate for the loss of RMS fluctuation
00178         // due to the temporal correlation introduced by the alpha parameter.
00179         scalar rmsCorr = sqrt(12*(2*alpha_ - sqr(alpha_)))/alpha_;
00180 
00181         patchField =
00182             (1 - alpha_)*patchField
00183           + alpha_*
00184             (
00185                 referenceField_
00186               + rmsCorr*cmptMultiply
00187                 (
00188                     randomField - 0.5*pTraits<Type>::one,
00189                     fluctuationScale_
00190                 )*mag(referenceField_)
00191             );
00192 
00193         curTimeIndex_ = this->db().time().timeIndex();
00194     }
00195 
00196     fixedValueFvPatchField<Type>::updateCoeffs();
00197 }
00198 
00199 
00200 template<class Type>
00201 void turbulentInletFvPatchField<Type>::write(Ostream& os) const
00202 {
00203     fvPatchField<Type>::write(os);
00204     os.writeKeyword("fluctuationScale")
00205         << fluctuationScale_ << token::END_STATEMENT << nl;
00206     referenceField_.writeEntry("referenceField", os);
00207     os.writeKeyword("alpha") << alpha_ << token::END_STATEMENT << nl;
00208     this->writeEntry("value", os);
00209 }
00210 
00211 
00212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00213 
00214 } // End namespace Foam
00215 
00216 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines