BUG: Corrected reference to tmp in incompressible RAS models (mantis #360)
[OpenFOAM-2.0.x.git] / src / turbulenceModels / incompressible / RAS / derivedFvPatchFields / wallFunctions / nutWallFunctions / nutkRoughWallFunction / nutkRoughWallFunctionFvPatchScalarField.C
blob9dc82a03ad2e15f4d6adfdcd43d4f681b19d656a
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 "nutkRoughWallFunctionFvPatchScalarField.H"
27 #include "RASModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "addToRunTimeSelectionTable.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
36 namespace incompressible
38 namespace RASModels
41 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
43 scalar nutkRoughWallFunctionFvPatchScalarField::fnRough
45     const scalar KsPlus,
46     const scalar Cs
47 ) const
49     // Return fn based on non-dimensional roughness height
51     if (KsPlus < 90.0)
52     {
53         return pow
54         (
55             (KsPlus - 2.25)/87.75 + Cs*KsPlus,
56             sin(0.4258*(log(KsPlus) - 0.811))
57         );
58     }
59     else
60     {
61         return (1.0 + Cs*KsPlus);
62     }
66 tmp<scalarField> nutkRoughWallFunctionFvPatchScalarField::calcNut() const
68     const label patchI = patch().index();
70     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
71     const scalarField& y = rasModel.y()[patchI];
72     const tmp<volScalarField> tk = rasModel.k();
73     const volScalarField& k = tk();
74     const tmp<volScalarField> tnu = rasModel.nu();
75     const volScalarField& nu = tnu();
76     const scalarField& nuw = nu.boundaryField()[patchI];
78     const scalar Cmu25 = pow025(Cmu_);
80     tmp<scalarField> tnutw(new scalarField(*this));
81     scalarField& nutw = tnutw();
83     forAll(nutw, faceI)
84     {
85         label faceCellI = patch().faceCells()[faceI];
87         scalar uStar = Cmu25*sqrt(k[faceCellI]);
88         scalar yPlus = uStar*y[faceI]/nuw[faceI];
89         scalar KsPlus = uStar*Ks_[faceI]/nuw[faceI];
91         scalar Edash = E_;
92         if (KsPlus > 2.25)
93         {
94             Edash /= fnRough(KsPlus, Cs_[faceI]);
95         }
97         scalar limitingNutw = max(nutw[faceI], nuw[faceI]);
99         // To avoid oscillations limit the change in the wall viscosity
100         // which is particularly important if it temporarily becomes zero
101         nutw[faceI] =
102             max
103             (
104                 min
105                 (
106                     nuw[faceI]
107                    *(yPlus*kappa_/log(max(Edash*yPlus, 1+1e-4)) - 1),
108                     2*limitingNutw
109                 ), 0.5*limitingNutw
110             );
112         if (debug)
113         {
114             Info<< "yPlus = " << yPlus
115                 << ", KsPlus = " << KsPlus
116                 << ", Edash = " << Edash
117                 << ", nutw = " << nutw[faceI]
118                 << endl;
119         }
120     }
122     return tnutw;
126 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
128 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
130     const fvPatch& p,
131     const DimensionedField<scalar, volMesh>& iF
134     nutkWallFunctionFvPatchScalarField(p, iF),
135     Ks_(p.size(), 0.0),
136     Cs_(p.size(), 0.0)
140 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
142     const nutkRoughWallFunctionFvPatchScalarField& ptf,
143     const fvPatch& p,
144     const DimensionedField<scalar, volMesh>& iF,
145     const fvPatchFieldMapper& mapper
148     nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
149     Ks_(ptf.Ks_, mapper),
150     Cs_(ptf.Cs_, mapper)
154 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
156     const fvPatch& p,
157     const DimensionedField<scalar, volMesh>& iF,
158     const dictionary& dict
161     nutkWallFunctionFvPatchScalarField(p, iF, dict),
162     Ks_("Ks", dict, p.size()),
163     Cs_("Cs", dict, p.size())
167 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
169     const nutkRoughWallFunctionFvPatchScalarField& rwfpsf
172     nutkWallFunctionFvPatchScalarField(rwfpsf),
173     Ks_(rwfpsf.Ks_),
174     Cs_(rwfpsf.Cs_)
178 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
180     const nutkRoughWallFunctionFvPatchScalarField& rwfpsf,
181     const DimensionedField<scalar, volMesh>& iF
184     nutkWallFunctionFvPatchScalarField(rwfpsf, iF),
185     Ks_(rwfpsf.Ks_),
186     Cs_(rwfpsf.Cs_)
190 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
192 void nutkRoughWallFunctionFvPatchScalarField::autoMap
194     const fvPatchFieldMapper& m
197     nutkWallFunctionFvPatchScalarField::autoMap(m);
198     Ks_.autoMap(m);
199     Cs_.autoMap(m);
203 void nutkRoughWallFunctionFvPatchScalarField::rmap
205     const fvPatchScalarField& ptf,
206     const labelList& addr
209     nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
211     const nutkRoughWallFunctionFvPatchScalarField& nrwfpsf =
212         refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
214     Ks_.rmap(nrwfpsf.Ks_, addr);
215     Cs_.rmap(nrwfpsf.Cs_, addr);
219 void nutkRoughWallFunctionFvPatchScalarField::write(Ostream& os) const
221     fvPatchField<scalar>::write(os);
222     writeLocalEntries(os);
223     Cs_.writeEntry("Cs", os);
224     Ks_.writeEntry("Ks", os);
225     writeEntry("value", os);
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 makePatchTypeField
233     fvPatchScalarField,
234     nutkRoughWallFunctionFvPatchScalarField
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 } // End namespace RASModels
240 } // End namespace incompressible
241 } // End namespace Foam
243 // ************************************************************************* //