BUG: Corrected reference to tmp in incompressible RAS models (mantis #360)
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / derivedFvPatchFields / wallFunctions / nutWallFunctions / nutkWallFunction / nutkWallFunctionFvPatchScalarField.C
blob157b88731a567d58f1d505278ddc8ce074271c5c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "nutkWallFunctionFvPatchScalarField.H"
27 #include "RASModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "wallFvPatch.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
37 namespace incompressible
39 namespace RASModels
42 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
44 void nutkWallFunctionFvPatchScalarField::checkType()
46     if (!isA<wallFvPatch>(patch()))
47     {
48         FatalErrorIn("nutkWallFunctionFvPatchScalarField::checkType()")
49             << "Invalid wall function specification" << nl
50             << "    Patch type for patch " << patch().name()
51             << " must be wall" << nl
52             << "    Current patch type is " << patch().type() << nl << endl
53             << abort(FatalError);
54     }
58 scalar nutkWallFunctionFvPatchScalarField::calcYPlusLam
60     const scalar kappa,
61     const scalar E
62 ) const
64     scalar ypl = 11.0;
66     for (int i=0; i<10; i++)
67     {
68         ypl = log(E*ypl)/kappa;
69     }
71     return ypl;
75 tmp<scalarField> nutkWallFunctionFvPatchScalarField::calcNut() const
77     const label patchI = patch().index();
79     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
80     const scalarField& y = rasModel.y()[patchI];
81     const tmp<volScalarField> tk = rasModel.k();
82     const volScalarField& k = tk();
83     const tmp<volScalarField> tnu = rasModel.nu();
84     const volScalarField& nu = tnu();
85     const scalarField& nuw = nu.boundaryField()[patchI];
87     const scalar Cmu25 = pow025(Cmu_);
89     tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
90     scalarField& nutw = tnutw();
92     forAll(nutw, faceI)
93     {
94         label faceCellI = patch().faceCells()[faceI];
96         scalar yPlus = Cmu25*y[faceI]*sqrt(k[faceCellI])/nuw[faceI];
98         if (yPlus > yPlusLam_)
99         {
100             nutw[faceI] = nuw[faceI]*(yPlus*kappa_/log(E_*yPlus) - 1.0);
101         }
102     }
104     return tnutw;
108 void nutkWallFunctionFvPatchScalarField::writeLocalEntries(Ostream& os) const
110     os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
111     os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
112     os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
116 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
118 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
120     const fvPatch& p,
121     const DimensionedField<scalar, volMesh>& iF
124     fixedValueFvPatchScalarField(p, iF),
125     Cmu_(0.09),
126     kappa_(0.41),
127     E_(9.8),
128     yPlusLam_(calcYPlusLam(kappa_, E_))
130     checkType();
134 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
136     const nutkWallFunctionFvPatchScalarField& ptf,
137     const fvPatch& p,
138     const DimensionedField<scalar, volMesh>& iF,
139     const fvPatchFieldMapper& mapper
142     fixedValueFvPatchScalarField(ptf, p, iF, mapper),
143     Cmu_(ptf.Cmu_),
144     kappa_(ptf.kappa_),
145     E_(ptf.E_),
146     yPlusLam_(ptf.yPlusLam_)
148     checkType();
152 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
154     const fvPatch& p,
155     const DimensionedField<scalar, volMesh>& iF,
156     const dictionary& dict
159     fixedValueFvPatchScalarField(p, iF, dict),
160     Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
161     kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
162     E_(dict.lookupOrDefault<scalar>("E", 9.8)),
163     yPlusLam_(calcYPlusLam(kappa_, E_))
165     checkType();
169 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
171     const nutkWallFunctionFvPatchScalarField& wfpsf
174     fixedValueFvPatchScalarField(wfpsf),
175     Cmu_(wfpsf.Cmu_),
176     kappa_(wfpsf.kappa_),
177     E_(wfpsf.E_),
178     yPlusLam_(wfpsf.yPlusLam_)
180     checkType();
184 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
186     const nutkWallFunctionFvPatchScalarField& wfpsf,
187     const DimensionedField<scalar, volMesh>& iF
190     fixedValueFvPatchScalarField(wfpsf, iF),
191     Cmu_(wfpsf.Cmu_),
192     kappa_(wfpsf.kappa_),
193     E_(wfpsf.E_),
194     yPlusLam_(wfpsf.yPlusLam_)
196     checkType();
200 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
202 void nutkWallFunctionFvPatchScalarField::updateCoeffs()
204     if (updated())
205     {
206         return;
207     }
209     operator==(calcNut());
211     fixedValueFvPatchScalarField::updateCoeffs();
215 tmp<scalarField> nutkWallFunctionFvPatchScalarField::yPlus() const
217     const label patchI = patch().index();
219     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
220     const scalarField& y = rasModel.y()[patchI];
222     const tmp<volScalarField> tk = rasModel.k();
223     const volScalarField& k = tk();
224     tmp<scalarField> kwc = k.boundaryField()[patchI].patchInternalField();
225     const tmp<volScalarField> tnu = rasModel.nu();
226     const volScalarField& nu = tnu();
227     const scalarField& nuw = nu.boundaryField()[patchI];
229     return pow025(Cmu_)*y*sqrt(kwc)/nuw;
233 void nutkWallFunctionFvPatchScalarField::write(Ostream& os) const
235     fvPatchField<scalar>::write(os);
236     writeLocalEntries(os);
237     writeEntry("value", os);
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 makePatchTypeField
245     fvPatchScalarField,
246     nutkWallFunctionFvPatchScalarField
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 } // End namespace RASModels
252 } // End namespace incompressible
253 } // End namespace Foam
255 // ************************************************************************* //