/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2020-2024 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 "volFields.H" #include "surfaceFields.H" #include "fvcGrad.H" #include "coupledFvPatchFields.H" // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * // template class LimitFunc> void Foam::LimitedScheme::calcLimiter ( const GeometricField& phi, surfaceScalarField& limiterField ) const { typedef GeometricField VolFieldType; typedef GeometricField GradVolFieldType; const fvMesh& mesh = this->mesh(); tmp tlPhi = LimitFunc()(phi); const VolFieldType& lPhi = tlPhi(); tmp tgradc(fvc::grad(lPhi)); const GradVolFieldType& gradc = tgradc(); const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights(); const labelUList& owner = mesh.owner(); const labelUList& neighbour = mesh.neighbour(); const vectorField& C = mesh.C(); scalarField& pLim = limiterField.primitiveFieldRef(); forAll(pLim, face) { label own = owner[face]; label nei = neighbour[face]; pLim[face] = Limiter::limiter ( CDweights[face], this->faceFlux_[face], lPhi[own], lPhi[nei], gradc[own], gradc[nei], C[nei] - C[own] ); } surfaceScalarField::Boundary& bLim = limiterField.boundaryFieldRef(); forAll(bLim, patchi) { scalarField& pLim = bLim[patchi]; if (bLim[patchi].coupled()) { const scalarField& pCDweights = CDweights.boundaryField()[patchi]; const scalarField& pFaceFlux = this->faceFlux_.boundaryField()[patchi]; const Field plPhiP ( lPhi.boundaryField()[patchi].patchInternalField() ); const Field plPhiN ( lPhi.boundaryField()[patchi].patchNeighbourField() ); const Field pGradcP ( gradc.boundaryField()[patchi].patchInternalField() ); const Field pGradcN ( gradc.boundaryField()[patchi].patchNeighbourField() ); // Build the d-vectors vectorField pd(CDweights.boundaryField()[patchi].patch().delta()); forAll(pLim, face) { pLim[face] = Limiter::limiter ( pCDweights[face], pFaceFlux[face], plPhiP[face], plPhiN[face], pGradcP[face], pGradcN[face], pd[face] ); } } else { pLim = 1.0; } } limiterField.setOriented(); } // * * * * * * * * * * * * Public Member Functions * * * * * * * * * * * * // template class LimitFunc> Foam::tmp Foam::LimitedScheme::limiter ( const GeometricField& phi ) const { const fvMesh& mesh = this->mesh(); const word limiterFieldName(type() + "Limiter(" + phi.name() + ')'); if (this->mesh().cache("limiter")) { auto* fldptr = mesh.getObjectPtr(limiterFieldName); if (!fldptr) { fldptr = new surfaceScalarField ( IOobject ( limiterFieldName, mesh.time().timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE, IOobject::REGISTER ), mesh, dimless ); regIOobject::store(fldptr); } auto& limiterField = *fldptr; calcLimiter(phi, limiterField); return tmp::New ( limiterFieldName, limiterField ); } else { auto tlimiterField = surfaceScalarField::New ( limiterFieldName, mesh, dimless ); calcLimiter(phi, tlimiterField.ref()); return tlimiterField; } } // ************************************************************************* //