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 "mutUSpaldingWallFunctionFvPatchScalarField.H"
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
30 #include "addToRunTimeSelectionTable.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace compressible
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 tmp<scalarField> mutUSpaldingWallFunctionFvPatchScalarField::calcUTau
45 const scalarField& magGradU
48 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
49 const scalarField& y = rasModel.y()[patch().index()];
51 const fvPatchVectorField& Uw =
52 rasModel.U().boundaryField()[patch().index()];
54 scalarField magUp(mag(Uw.patchInternalField() - Uw));
56 const fvPatchScalarField& rhow =
57 rasModel.rho().boundaryField()[patch().index()];
59 const fvPatchScalarField& muw =
60 rasModel.mu().boundaryField()[patch().index()];
61 const scalarField& mutw = *this;
63 tmp<scalarField> tuTau(new scalarField(patch().size(), 0.0));
64 scalarField& uTau = tuTau();
69 sqrt((mutw[faceI] + muw[faceI])*magGradU[faceI]/rhow[faceI]);
78 scalar kUu = min(kappa_*magUp[faceI]/ut, 50);
79 scalar fkUu = exp(kUu) - 1 - kUu*(1 + 0.5*kUu);
82 - ut*y[faceI]/(muw[faceI]/rhow[faceI])
84 + 1/E_*(fkUu - 1.0/6.0*kUu*sqr(kUu));
87 y[faceI]/(muw[faceI]/rhow[faceI])
88 + magUp[faceI]/sqr(ut)
91 scalar uTauNew = ut + f/df;
92 err = mag((ut - uTauNew)/ut);
95 } while (ut > ROOTVSMALL && err > 0.01 && ++iter < 10);
97 uTau[faceI] = max(0.0, ut);
105 tmp<scalarField> mutUSpaldingWallFunctionFvPatchScalarField::calcMut() const
107 const label patchI = patch().index();
109 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
110 const fvPatchVectorField& Uw = rasModel.U().boundaryField()[patchI];
111 const scalarField magGradU(mag(Uw.snGrad()));
112 const scalarField& rhow = rasModel.rho().boundaryField()[patchI];
113 const scalarField& muw = rasModel.mu().boundaryField()[patchI];
118 rhow*sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - muw
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125 mutUSpaldingWallFunctionFvPatchScalarField::
126 mutUSpaldingWallFunctionFvPatchScalarField
129 const DimensionedField<scalar, volMesh>& iF
132 mutkWallFunctionFvPatchScalarField(p, iF)
136 mutUSpaldingWallFunctionFvPatchScalarField::
137 mutUSpaldingWallFunctionFvPatchScalarField
139 const mutUSpaldingWallFunctionFvPatchScalarField& ptf,
141 const DimensionedField<scalar, volMesh>& iF,
142 const fvPatchFieldMapper& mapper
145 mutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
149 mutUSpaldingWallFunctionFvPatchScalarField::
150 mutUSpaldingWallFunctionFvPatchScalarField
153 const DimensionedField<scalar, volMesh>& iF,
154 const dictionary& dict
157 mutkWallFunctionFvPatchScalarField(p, iF, dict)
161 mutUSpaldingWallFunctionFvPatchScalarField::
162 mutUSpaldingWallFunctionFvPatchScalarField
164 const mutUSpaldingWallFunctionFvPatchScalarField& wfpsf
167 mutkWallFunctionFvPatchScalarField(wfpsf)
171 mutUSpaldingWallFunctionFvPatchScalarField::
172 mutUSpaldingWallFunctionFvPatchScalarField
174 const mutUSpaldingWallFunctionFvPatchScalarField& wfpsf,
175 const DimensionedField<scalar, volMesh>& iF
178 mutkWallFunctionFvPatchScalarField(wfpsf, iF)
182 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
184 tmp<scalarField> mutUSpaldingWallFunctionFvPatchScalarField::yPlus() const
186 const label patchI = patch().index();
188 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
189 const scalarField& y = rasModel.y()[patchI];
190 const fvPatchVectorField& Uw = rasModel.U().boundaryField()[patchI];
191 const scalarField& rhow = rasModel.rho().boundaryField()[patchI];
192 const scalarField& muw = rasModel.mu().boundaryField()[patchI];
194 return y*calcUTau(mag(Uw.snGrad()))/(muw/rhow);
198 void mutUSpaldingWallFunctionFvPatchScalarField::write(Ostream& os) const
200 fvPatchField<scalar>::write(os);
201 writeLocalEntries(os);
202 writeEntry("value", os);
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 mutUSpaldingWallFunctionFvPatchScalarField
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 } // End namespace RASModels
217 } // End namespace compressible
218 } // End namespace Foam
220 // ************************************************************************* //