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 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();
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];
94 Edash /= fnRough(KsPlus, Cs_[faceI]);
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
107 *(yPlus*kappa_/log(max(Edash*yPlus, 1+1e-4)) - 1),
114 Info<< "yPlus = " << yPlus
115 << ", KsPlus = " << KsPlus
116 << ", Edash = " << Edash
117 << ", nutw = " << nutw[faceI]
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
128 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
131 const DimensionedField<scalar, volMesh>& iF
134 nutkWallFunctionFvPatchScalarField(p, iF),
140 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
142 const nutkRoughWallFunctionFvPatchScalarField& ptf,
144 const DimensionedField<scalar, volMesh>& iF,
145 const fvPatchFieldMapper& mapper
148 nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
149 Ks_(ptf.Ks_, mapper),
154 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
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),
178 nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
180 const nutkRoughWallFunctionFvPatchScalarField& rwfpsf,
181 const DimensionedField<scalar, volMesh>& iF
184 nutkWallFunctionFvPatchScalarField(rwfpsf, iF),
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 void nutkRoughWallFunctionFvPatchScalarField::autoMap
194 const fvPatchFieldMapper& m
197 nutkWallFunctionFvPatchScalarField::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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 nutkRoughWallFunctionFvPatchScalarField
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 } // End namespace RASModels
240 } // End namespace incompressible
241 } // End namespace Foam
243 // ************************************************************************* //