ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / regionModels / surfaceFilmModels / derivedFvPatchFields / wallFunctions / mutkFilmWallFunction / mutkFilmWallFunctionFvPatchScalarField.C
blob348821df9bc853a7ff8a580bba1cc26b3b6c54e3
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 "mutkFilmWallFunctionFvPatchScalarField.H"
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "RASModel.H"
30 #include "addToRunTimeSelectionTable.H"
31 #include "surfaceFilmModel.H"
32 #include "directMappedWallPolyPatch.H"
33 #include "mapDistribute.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace Foam
39 namespace compressible
41 namespace RASModels
44 // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
46 tmp<scalarField> mutkFilmWallFunctionFvPatchScalarField::calcUTau
48     const scalarField& magGradU
49 ) const
51     tmp<scalarField> tuTau(new scalarField(patch().size(), 0.0));
52     scalarField& uTau = tuTau();
54     typedef regionModels::surfaceFilmModels::surfaceFilmModel modelType;
56     bool ok =
57         db().objectRegistry::foundObject<modelType>("surfaceFilmProperties");
59     if (!ok)
60     {
61         // do nothing on construction - film model doesn't exist yet
62         return tuTau;
63     }
65     const label patchI = patch().index();
67     // Retrieve phase change mass from surface film model
68     const modelType& filmModel =
69         db().objectRegistry::lookupObject<modelType>("surfaceFilmProperties");
71     const label filmPatchI = filmModel.regionPatchID(patchI);
73     const mapDistribute& distMap = filmModel.mappedPatches()[filmPatchI].map();
75     tmp<volScalarField> mDotFilm(filmModel.primaryMassTrans());
76     scalarField mDotFilmp = mDotFilm().boundaryField()[filmPatchI];
77     distMap.distribute(mDotFilmp);
80     // Retrieve RAS turbulence model
81     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
82     const scalarField& y = rasModel.y()[patchI];
83     const fvPatchScalarField& rhow = rasModel.rho().boundaryField()[patchI];
84     const fvPatchScalarField& muw = rasModel.mu().boundaryField()[patchI];
85     const tmp<volScalarField> tk = rasModel.k();
86     const volScalarField& k = tk();
88     const scalar Cmu25 = pow(Cmu_, 0.25);
90     forAll(uTau, faceI)
91     {
92         label faceCellI = patch().faceCells()[faceI];
94         scalar ut = Cmu25*sqrt(k[faceCellI]);
96         scalar yPlus = y[faceI]*ut/(muw[faceI]/rhow[faceI]);
98         scalar mStar = mDotFilmp[faceI]/(y[faceI]*ut);
100         scalar factor = 0.0;
101         if (yPlus > yPlusCrit_)
102         {
103             scalar expTerm = exp(min(50.0, B_*mStar));
104             scalar powTerm = pow(yPlus, mStar/kappa_);
105             factor = mStar/(expTerm*powTerm - 1.0 + ROOTVSMALL);
106         }
107         else
108         {
109             scalar expTerm = exp(min(50.0, mStar));
110             factor = mStar/(expTerm*yPlus - 1.0 + ROOTVSMALL);
111         }
113         uTau[faceI] = sqrt(max(0, magGradU[faceI]*ut*factor));
114     }
116     return tuTau;
120 tmp<scalarField> mutkFilmWallFunctionFvPatchScalarField::calcMut() const
122     const label patchI = patch().index();
124     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
125     const fvPatchVectorField& Uw = rasModel.U().boundaryField()[patchI];
126     const scalarField magGradU(mag(Uw.snGrad()));
127     const scalarField& rhow = rasModel.rho().boundaryField()[patchI];
128     const scalarField& muw = rasModel.mu().boundaryField()[patchI];
130     return max
131     (
132         scalar(0),
133         rhow*sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - muw
134     );
138 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
140 mutkFilmWallFunctionFvPatchScalarField::mutkFilmWallFunctionFvPatchScalarField
142     const fvPatch& p,
143     const DimensionedField<scalar, volMesh>& iF
146     mutkWallFunctionFvPatchScalarField(p, iF),
147     B_(5.5),
148     yPlusCrit_(11.05)
152 mutkFilmWallFunctionFvPatchScalarField::mutkFilmWallFunctionFvPatchScalarField
154     const mutkFilmWallFunctionFvPatchScalarField& ptf,
155     const fvPatch& p,
156     const DimensionedField<scalar, volMesh>& iF,
157     const fvPatchFieldMapper& mapper
160     mutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
161     B_(5.5),
162     yPlusCrit_(11.05)
166 mutkFilmWallFunctionFvPatchScalarField::mutkFilmWallFunctionFvPatchScalarField
168     const fvPatch& p,
169     const DimensionedField<scalar, volMesh>& iF,
170     const dictionary& dict
173     mutkWallFunctionFvPatchScalarField(p, iF, dict),
174     B_(dict.lookupOrDefault("B", 5.5)),
175     yPlusCrit_(dict.lookupOrDefault("yPlusCrit", 11.05))
179 mutkFilmWallFunctionFvPatchScalarField::mutkFilmWallFunctionFvPatchScalarField
181     const mutkFilmWallFunctionFvPatchScalarField& wfpsf
184     mutkWallFunctionFvPatchScalarField(wfpsf),
185     B_(wfpsf.B_),
186     yPlusCrit_(wfpsf.yPlusCrit_)
190 mutkFilmWallFunctionFvPatchScalarField::mutkFilmWallFunctionFvPatchScalarField
192     const mutkFilmWallFunctionFvPatchScalarField& wfpsf,
193     const DimensionedField<scalar, volMesh>& iF
196     mutkWallFunctionFvPatchScalarField(wfpsf, iF),
197     B_(wfpsf.B_),
198     yPlusCrit_(wfpsf.yPlusCrit_)
202 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
204 tmp<scalarField> mutkFilmWallFunctionFvPatchScalarField::yPlus() const
206     const label patchI = patch().index();
208     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
209     const scalarField& y = rasModel.y()[patchI];
210     const fvPatchVectorField& Uw = rasModel.U().boundaryField()[patchI];
211     const scalarField& rhow = rasModel.rho().boundaryField()[patchI];
212     const scalarField& muw = rasModel.mu().boundaryField()[patchI];
214     return y*calcUTau(mag(Uw.snGrad()))/(muw/rhow);
218 void mutkFilmWallFunctionFvPatchScalarField::write(Ostream& os) const
220     fvPatchField<scalar>::write(os);
221     writeLocalEntries(os);
222     os.writeKeyword("B") << B_ << token::END_STATEMENT << nl;
223     os.writeKeyword("yPlusCrit") << yPlusCrit_ << token::END_STATEMENT << nl;
224     writeEntry("value", os);
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 makePatchTypeField(fvPatchScalarField, mutkFilmWallFunctionFvPatchScalarField);
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 } // End namespace RASModels
235 } // End namespace compressible
236 } // End namespace Foam
238 // ************************************************************************* //