.. _program_listing_file_wallModels_ExplicitWallModelFvPatchScalarField.C: Program Listing for File ExplicitWallModelFvPatchScalarField.C ============================================================== |exhale_lsh| :ref:`Return to documentation for file ` (``wallModels/ExplicitWallModelFvPatchScalarField.C``) .. |exhale_lsh| unicode:: U+021B0 .. UPWARDS ARROW WITH TIP LEFTWARDS .. code-block:: cpp /*---------------------------------------------------------------------------* \ License This file is part of libWallModelledLES. libWallModelledLES 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. libWallModelledLES 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 libWallModelledLES. If not, see . \*---------------------------------------------------------------------------*/ #include "ExplicitWallModelFvPatchScalarField.H" #include "fvPatchFieldMapper.H" #include "addToRunTimeSelectionTable.H" #include "codeRules.H" #include "scalarListIOList.H" #include "helpers.H" #include "ExplicitLawOfTheWall.H" #include "SingleCellSampler.H" // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * // void Foam::ExplicitWallModelFvPatchScalarField::writeLocalEntries(Ostream& os) const { wallModelFvPatchScalarField::writeLocalEntries(os); law_->write(os); sampler_->write(os); } Foam::tmp Foam::ExplicitWallModelFvPatchScalarField::calcNut() const { if (debug) { Info<< "Updating nut for patch " << patch().name() << nl; } const label patchi = patch().index(); tmp nuw = this->nu(patchi); const scalarListIOList & wallGradU = sampler_->db().lookupObject("wallGradU"); scalarField magGradU(Helpers::mag(wallGradU)); return max ( scalar(0), sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw ); } Foam::tmp Foam::ExplicitWallModelFvPatchScalarField:: calcUTau(const scalarField & magGradU) const { const label patchi = patch().index(); const label patchSize = patch().size(); tmp tnuw = this->nu(patchi); const scalarField& nuw = tnuw(); // Computed uTau tmp tuTau(new scalarField(patchSize, 0.0)); scalarField & uTau = tuTau.ref(); // Grab global uTau field volScalarField & uTauField = const_cast ( db().lookupObject("uTauPredicted") ); // Compute uTau for each face forAll(uTau, faceI) { scalar nuwI = nuw[faceI]; // Compute uTau uTau[faceI] = max(0.0, law_().uTau(sampler_(), faceI, nuwI)); } // Assign computed uTau to the boundary field of the global field uTauField.boundaryFieldRef()[patchi] == uTau; return tuTau; } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::ExplicitWallModelFvPatchScalarField:: ExplicitWallModelFvPatchScalarField ( const fvPatch & p, const DimensionedField & iF ) : wallModelFvPatchScalarField(p, iF), law_(), sampler_(nullptr) { if (debug) { Info<< "Constructing ExplicitWallModelFvPatchScalarField (lotw1) " << "from fvPatch and DimensionedField for patch " << patch().name() << nl; } } Foam::ExplicitWallModelFvPatchScalarField:: ExplicitWallModelFvPatchScalarField ( const ExplicitWallModelFvPatchScalarField & orig, const fvPatch & p, const DimensionedField & iF, const fvPatchFieldMapper & mapper ) : wallModelFvPatchScalarField(orig, p, iF, mapper), #ifdef FOAM_AUTOPTR_HAS_CLONE_METHOD law_(orig.law_.clone()), #else law_(orig.law_, false), #endif sampler_(new SingleCellSampler(orig.sampler())) { if (debug) { Info<< "Constructing ExplicitWallModelFvPatchScalarField (lotw2) " << "from copy, fvPatch, DimensionedField, and fvPatchFieldMapper" << " for patch " << patch().name() << nl; } law_->addFieldsToSampler(sampler()); } Foam::ExplicitWallModelFvPatchScalarField:: ExplicitWallModelFvPatchScalarField ( const fvPatch & p, const DimensionedField & iF, const dictionary & dict ) : wallModelFvPatchScalarField(p, iF, dict), law_(ExplicitLawOfTheWall::New(dict.subDict("Law"))), sampler_ ( new SingleCellSampler ( p, averagingTime(), dict.lookupOrDefault("interpolationType", "cell"), dict.lookupOrDefault("sampler", "Tree"), dict.lookupOrDefault("lengthScale", "CubeRootVol"), dict.lookupOrDefault("hIsIndex", false) ) ) { if (debug) { Info<< "Constructing ExplicitWallModelFvPatchScalarField (lotw3) " << "from fvPatch, DimensionedField, and dictionary for patch " << patch().name() << nl; } law_->addFieldsToSampler(sampler()); } #ifdef FOAM_FVPATCHFIELD_NO_COPY #else Foam::ExplicitWallModelFvPatchScalarField:: ExplicitWallModelFvPatchScalarField ( const ExplicitWallModelFvPatchScalarField & orig ) : wallModelFvPatchScalarField(orig), #ifdef FOAM_AUTOPTR_HAS_CLONE_METHOD law_(orig.law_.clone()), #else law_(orig.law_, false), #endif sampler_(new SingleCellSampler(orig.sampler_())) { if (debug) { Info<< "Constructing ExplicitWallModelFvPatchScalarField (lotw4)" << "from copy for patch " << patch().name() << nl; } law_->addFieldsToSampler(sampler()); } #endif Foam::ExplicitWallModelFvPatchScalarField:: ExplicitWallModelFvPatchScalarField ( const ExplicitWallModelFvPatchScalarField & orig, const DimensionedField & iF ) : wallModelFvPatchScalarField(orig, iF), #ifdef FOAM_AUTOPTR_HAS_CLONE_METHOD law_(orig.law_.clone()), #else law_(orig.law_, false), #endif sampler_(new SingleCellSampler(orig.sampler_())) { if (debug) { Info<< "Constructing ExplicitLOTWModelFvPatchScalarField (lotw5) " << "from copy and DimensionedField for patch " << patch().name() << nl; } law_->addFieldsToSampler(sampler()); } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::ExplicitWallModelFvPatchScalarField::write(Ostream& os) const { wallModelFvPatchScalarField::write(os); } void Foam::ExplicitWallModelFvPatchScalarField::updateCoeffs() { if (updated()) { return; } sampler().recomputeFields(); sampler().sample(); wallModelFvPatchScalarField::updateCoeffs(); } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // #if !defined(DOXYGEN_SHOULD_SKIP_THIS) namespace Foam { makePatchTypeField ( fvPatchScalarField, ExplicitWallModelFvPatchScalarField ); } #endif // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //