1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2008-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "nutRoughWallFunctionFvPatchScalarField.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace incompressible
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 scalar nutRoughWallFunctionFvPatchScalarField::fnRough
50 // Return fn based on non-dimensional roughness height
56 (KsPlus - 2.25)/87.75 + Cs*KsPlus,
57 sin(0.4258*(log(KsPlus) - 0.811))
62 return (1.0 + Cs*KsPlus);
67 tmp<scalarField> nutRoughWallFunctionFvPatchScalarField::calcNut() const
69 const label patchI = patch().index();
71 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
72 const scalarField& y = rasModel.y()[patchI];
73 const tmp<volScalarField> tk = rasModel.k();
74 const volScalarField& k = tk();
75 const scalarField& nuw = rasModel.nu().boundaryField()[patchI];
77 const scalar Cmu25 = pow(Cmu_, 0.25);
79 tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
80 scalarField& nutw = tnutw();
84 label faceCellI = patch().faceCells()[faceI];
86 scalar uStar = Cmu25*sqrt(k[faceCellI]);
87 scalar yPlus = uStar*y[faceI]/nuw[faceI];
88 scalar KsPlus = uStar*Ks_[faceI]/nuw[faceI];
94 Edash /= fnRough(KsPlus, Cs_[faceI]);
97 if (yPlus > yPlusLam_)
99 scalar limitingNutw = max(nutw[faceI], nuw[faceI]);
101 // To avoid oscillations limit the change in the wall viscosity
102 // which is particularly important if it temporarily becomes zero
109 *(yPlus*kappa_/log(max(Edash*yPlus, 1+1e-4)) - 1),
117 Info<< "yPlus = " << yPlus
118 << ", KsPlus = " << KsPlus
119 << ", Edash = " << Edash
120 << ", nutw = " << nutw[faceI]
129 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
131 nutRoughWallFunctionFvPatchScalarField::nutRoughWallFunctionFvPatchScalarField
134 const DimensionedField<scalar, volMesh>& iF
137 nutWallFunctionFvPatchScalarField(p, iF),
143 nutRoughWallFunctionFvPatchScalarField::nutRoughWallFunctionFvPatchScalarField
145 const nutRoughWallFunctionFvPatchScalarField& ptf,
147 const DimensionedField<scalar, volMesh>& iF,
148 const fvPatchFieldMapper& mapper
151 nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
152 Ks_(ptf.Ks_, mapper),
157 nutRoughWallFunctionFvPatchScalarField::nutRoughWallFunctionFvPatchScalarField
160 const DimensionedField<scalar, volMesh>& iF,
161 const dictionary& dict
164 nutWallFunctionFvPatchScalarField(p, iF, dict),
165 Ks_("Ks", dict, p.size()),
166 Cs_("Cs", dict, p.size())
170 nutRoughWallFunctionFvPatchScalarField::nutRoughWallFunctionFvPatchScalarField
172 const nutRoughWallFunctionFvPatchScalarField& rwfpsf
175 nutWallFunctionFvPatchScalarField(rwfpsf),
181 nutRoughWallFunctionFvPatchScalarField::nutRoughWallFunctionFvPatchScalarField
183 const nutRoughWallFunctionFvPatchScalarField& rwfpsf,
184 const DimensionedField<scalar, volMesh>& iF
187 nutWallFunctionFvPatchScalarField(rwfpsf, iF),
193 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 void nutRoughWallFunctionFvPatchScalarField::autoMap
197 const fvPatchFieldMapper& m
200 nutWallFunctionFvPatchScalarField::autoMap(m);
206 void nutRoughWallFunctionFvPatchScalarField::rmap
208 const fvPatchScalarField& ptf,
209 const labelList& addr
212 nutWallFunctionFvPatchScalarField::rmap(ptf, addr);
214 const nutRoughWallFunctionFvPatchScalarField& nrwfpsf =
215 refCast<const nutRoughWallFunctionFvPatchScalarField>(ptf);
217 Ks_.rmap(nrwfpsf.Ks_, addr);
218 Cs_.rmap(nrwfpsf.Cs_, addr);
222 void nutRoughWallFunctionFvPatchScalarField::write(Ostream& os) const
224 fvPatchField<scalar>::write(os);
225 writeLocalEntries(os);
226 Cs_.writeEntry("Cs", os);
227 Ks_.writeEntry("Ks", os);
228 writeEntry("value", os);
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 makePatchTypeField(fvPatchScalarField, nutRoughWallFunctionFvPatchScalarField);
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 } // End namespace RASModels
239 } // End namespace incompressible
240 } // End namespace Foam
242 // ************************************************************************* //