/*---------------------------------------------------------------------------* \ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | \\/ M anipulation | ------------------------------------------------------------------------------- | Copyright (C) 2011-2016, 2019 OpenFOAM Foundation ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . \*---------------------------------------------------------------------------*/ #include "nutWallFunctionFvPatchScalarField.H" #include "fvPatchFieldMapper.H" #include "volFields.H" #include "wallFvPatch.H" #include "turbulenceModel.H" #include "addToRunTimeSelectionTable.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(nutWallFunctionFvPatchScalarField, 0); } // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * // void Foam::nutWallFunctionFvPatchScalarField::checkType() { if (!isA(patch())) { FatalErrorInFunction << "Invalid wall function specification" << nl << " Patch type for patch " << patch().name() << " must be wall" << nl << " Current patch type is " << patch().type() << nl << endl << abort(FatalError); } } const Foam::volVectorField& Foam::nutWallFunctionFvPatchScalarField::U ( const turbulenceModel& turb ) const { if (UName_ == word::null) { return turb.U(); } else { return db().lookupObject(UName_); } } void Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries ( Ostream& os ) const { if (UName_ != word::null) { os.writeEntry("U", UName_); } os.writeEntry("Cmu", Cmu_); os.writeEntry("kappa", kappa_); os.writeEntry("E", E_); } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField ( const fvPatch& p, const DimensionedField& iF ) : fixedValueFvPatchScalarField(p, iF), UName_(word::null), Cmu_(0.09), kappa_(0.41), E_(9.8), yPlusLam_(yPlusLam(kappa_, E_)) { checkType(); } Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField ( const nutWallFunctionFvPatchScalarField& ptf, const fvPatch& p, const DimensionedField& iF, const fvPatchFieldMapper& mapper ) : fixedValueFvPatchScalarField(ptf, p, iF, mapper), UName_(ptf.UName_), Cmu_(ptf.Cmu_), kappa_(ptf.kappa_), E_(ptf.E_), yPlusLam_(ptf.yPlusLam_) { checkType(); } Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField ( const fvPatch& p, const DimensionedField& iF, const dictionary& dict ) : fixedValueFvPatchScalarField(p, iF, dict), UName_(dict.lookupOrDefault("U", word::null)), Cmu_(dict.lookupOrDefault("Cmu", 0.09)), kappa_(dict.lookupOrDefault("kappa", 0.41)), E_(dict.lookupOrDefault("E", 9.8)), yPlusLam_(yPlusLam(kappa_, E_)) { checkType(); } Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField ( const nutWallFunctionFvPatchScalarField& wfpsf ) : fixedValueFvPatchScalarField(wfpsf), UName_(wfpsf.UName_), Cmu_(wfpsf.Cmu_), kappa_(wfpsf.kappa_), E_(wfpsf.E_), yPlusLam_(wfpsf.yPlusLam_) { checkType(); } Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField ( const nutWallFunctionFvPatchScalarField& wfpsf, const DimensionedField& iF ) : fixedValueFvPatchScalarField(wfpsf, iF), UName_(wfpsf.UName_), Cmu_(wfpsf.Cmu_), kappa_(wfpsf.kappa_), E_(wfpsf.E_), yPlusLam_(wfpsf.yPlusLam_) { checkType(); } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // const Foam::nutWallFunctionFvPatchScalarField& Foam::nutWallFunctionFvPatchScalarField::nutw ( const turbulenceModel& turbModel, const label patchi ) { return refCast ( turbModel.nut()().boundaryField()[patchi] ); } Foam::scalar Foam::nutWallFunctionFvPatchScalarField::yPlusLam ( const scalar kappa, const scalar E ) { scalar ypl = 11.0; for (int i=0; i<10; i++) { ypl = log(max(E*ypl, 1))/kappa; } return ypl; } Foam::scalar Foam::nutWallFunctionFvPatchScalarField::yPlusLam() const { return yPlusLam_; } void Foam::nutWallFunctionFvPatchScalarField::updateCoeffs() { if (updated()) { return; } operator==(calcNut()); fixedValueFvPatchScalarField::updateCoeffs(); } void Foam::nutWallFunctionFvPatchScalarField::write(Ostream& os) const { fvPatchField::write(os); writeLocalEntries(os); writeEntry("value", os); } // ************************************************************************* //