Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / derivedFvPatchFields / wallFunctions / nutWallFunctions / nutkWallFunction / nutkWallFunctionFvPatchScalarField.C
blob7a095ce55eebd7c432e6598e95b6974147ee73d9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2010 OpenCFD Ltd.
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 scalarField& nuw = rasModel.nu()().boundaryField()[patchI];
85     const scalar Cmu25 = pow025(Cmu_);
87     tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
88     scalarField& nutw = tnutw();
90     forAll(nutw, faceI)
91     {
92         label faceCellI = patch().faceCells()[faceI];
94         scalar yPlus = Cmu25*y[faceI]*sqrt(k[faceCellI])/nuw[faceI];
96         if (yPlus > yPlusLam_)
97         {
98             nutw[faceI] = nuw[faceI]*(yPlus*kappa_/log(E_*yPlus) - 1.0);
99         }
100     }
102     return tnutw;
106 void nutkWallFunctionFvPatchScalarField::writeLocalEntries(Ostream& os) const
108     os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
109     os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
110     os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
114 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
116 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
118     const fvPatch& p,
119     const DimensionedField<scalar, volMesh>& iF
122     fixedValueFvPatchScalarField(p, iF),
123     Cmu_(0.09),
124     kappa_(0.41),
125     E_(9.8),
126     yPlusLam_(calcYPlusLam(kappa_, E_))
128     checkType();
132 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
134     const nutkWallFunctionFvPatchScalarField& ptf,
135     const fvPatch& p,
136     const DimensionedField<scalar, volMesh>& iF,
137     const fvPatchFieldMapper& mapper
140     fixedValueFvPatchScalarField(ptf, p, iF, mapper),
141     Cmu_(ptf.Cmu_),
142     kappa_(ptf.kappa_),
143     E_(ptf.E_),
144     yPlusLam_(ptf.yPlusLam_)
146     checkType();
150 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
152     const fvPatch& p,
153     const DimensionedField<scalar, volMesh>& iF,
154     const dictionary& dict
157     fixedValueFvPatchScalarField(p, iF, dict),
158     Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
159     kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
160     E_(dict.lookupOrDefault<scalar>("E", 9.8)),
161     yPlusLam_(calcYPlusLam(kappa_, E_))
163     checkType();
167 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
169     const nutkWallFunctionFvPatchScalarField& wfpsf
172     fixedValueFvPatchScalarField(wfpsf),
173     Cmu_(wfpsf.Cmu_),
174     kappa_(wfpsf.kappa_),
175     E_(wfpsf.E_),
176     yPlusLam_(wfpsf.yPlusLam_)
178     checkType();
182 nutkWallFunctionFvPatchScalarField::nutkWallFunctionFvPatchScalarField
184     const nutkWallFunctionFvPatchScalarField& wfpsf,
185     const DimensionedField<scalar, volMesh>& iF
188     fixedValueFvPatchScalarField(wfpsf, iF),
189     Cmu_(wfpsf.Cmu_),
190     kappa_(wfpsf.kappa_),
191     E_(wfpsf.E_),
192     yPlusLam_(wfpsf.yPlusLam_)
194     checkType();
198 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
200 void nutkWallFunctionFvPatchScalarField::updateCoeffs()
202     if (updated())
203     {
204         return;
205     }
207     operator==(calcNut());
209     fixedValueFvPatchScalarField::updateCoeffs();
213 tmp<scalarField> nutkWallFunctionFvPatchScalarField::yPlus() const
215     const label patchI = patch().index();
217     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
218     const scalarField& y = rasModel.y()[patchI];
220     const tmp<volScalarField> tk = rasModel.k();
221     const volScalarField& k = tk();
222     tmp<scalarField> kwc = k.boundaryField()[patchI].patchInternalField();
223     const scalarField& nuw = rasModel.nu()().boundaryField()[patchI];
225     return pow025(Cmu_)*y*sqrt(kwc)/nuw;
229 void nutkWallFunctionFvPatchScalarField::write(Ostream& os) const
231     fvPatchField<scalar>::write(os);
232     writeLocalEntries(os);
233     writeEntry("value", os);
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 makePatchTypeField
241     fvPatchScalarField,
242     nutkWallFunctionFvPatchScalarField
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 } // End namespace RASModels
248 } // End namespace incompressible
249 } // End namespace Foam
251 // ************************************************************************* //