1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
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
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"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "addToRunTimeSelectionTable.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace incompressible
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 scalar nutkRoughWallFunctionFvPatchScalarField::fnRough
49 // Return fn based on non-dimensional roughness height
55 (KsPlus - 2.25)/87.75 + Cs*KsPlus,
56 sin(0.4258*(log(KsPlus) - 0.811))
61 return (1.0 + Cs*KsPlus);
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 scalarField& nuw = rasModel.nu()().boundaryField()[patchI];
76 const scalar Cmu25 = pow025(Cmu_);
78 tmp<scalarField> tnutw(new scalarField(*this));
79 scalarField& nutw = tnutw();
83 label faceCellI = patch().faceCells()[faceI];
85 scalar uStar = Cmu25*sqrt(k[faceCellI]);
86 scalar yPlus = uStar*y[faceI]/nuw[faceI];
87 scalar KsPlus = uStar*Ks_[faceI]/nuw[faceI];
93 Edash /= fnRough(KsPlus, Cs_[faceI]);
96 if (yPlus > yPlusLam_)
98 scalar limitingNutw = max(nutw[faceI], nuw[faceI]);
100 // To avoid oscillations limit the change in the wall viscosity
101 // which is particularly important if it temporarily becomes zero
108 *(yPlus*kappa_/log(max(Edash*yPlus, 1+1e-4)) - 1),
116 Info<< "yPlus = " << yPlus
117 << ", KsPlus = " << KsPlus
118 << ", Edash = " << Edash
119 << ", nutw = " << nutw[faceI]
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
130 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
133 const DimensionedField<scalar, volMesh>& iF
136 nutkWallFunctionFvPatchScalarField(p, iF),
142 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
144 const nutkRoughWallFunctionFvPatchScalarField& ptf,
146 const DimensionedField<scalar, volMesh>& iF,
147 const fvPatchFieldMapper& mapper
150 nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
151 Ks_(ptf.Ks_, mapper),
156 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
159 const DimensionedField<scalar, volMesh>& iF,
160 const dictionary& dict
163 nutkWallFunctionFvPatchScalarField(p, iF, dict),
164 Ks_("Ks", dict, p.size()),
165 Cs_("Cs", dict, p.size())
169 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
171 const nutkRoughWallFunctionFvPatchScalarField& rwfpsf
174 nutkWallFunctionFvPatchScalarField(rwfpsf),
180 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
182 const nutkRoughWallFunctionFvPatchScalarField& rwfpsf,
183 const DimensionedField<scalar, volMesh>& iF
186 nutkWallFunctionFvPatchScalarField(rwfpsf, iF),
192 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194 void nutkRoughWallFunctionFvPatchScalarField::autoMap
196 const fvPatchFieldMapper& m
199 nutkWallFunctionFvPatchScalarField::autoMap(m);
205 void nutkRoughWallFunctionFvPatchScalarField::rmap
207 const fvPatchScalarField& ptf,
208 const labelList& addr
211 nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
213 const nutkRoughWallFunctionFvPatchScalarField& nrwfpsf =
214 refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
216 Ks_.rmap(nrwfpsf.Ks_, addr);
217 Cs_.rmap(nrwfpsf.Cs_, addr);
221 void nutkRoughWallFunctionFvPatchScalarField::write(Ostream& os) const
223 fvPatchField<scalar>::write(os);
224 writeLocalEntries(os);
225 Cs_.writeEntry("Cs", os);
226 Ks_.writeEntry("Ks", os);
227 writeEntry("value", os);
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 nutkRoughWallFunctionFvPatchScalarField
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 } // End namespace RASModels
242 } // End namespace incompressible
243 } // End namespace Foam
245 // ************************************************************************* //