/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 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 "volFields.H" #include "surfaceFields.H" #include "fvcGrad.H" #include "coupledFvPatchFields.H" #include "surfaceInterpolate.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // template Foam::tmp Foam::PhiScheme::limiter ( const GeometricField& phi ) const { const fvMesh& mesh = this->mesh(); tmp tLimiter ( new surfaceScalarField ( IOobject ( "PhiLimiter", mesh.time().timeName(), mesh ), mesh, dimless ) ); surfaceScalarField& Limiter = tLimiter.ref(); const surfaceScalarField& CDweights = mesh.surfaceInterpolation::weights(); const surfaceVectorField& Sf = mesh.Sf(); const surfaceScalarField& magSf = mesh.magSf(); const labelUList& owner = mesh.owner(); const labelUList& neighbour = mesh.neighbour(); tmp tUflux = this->faceFlux_; if (this->faceFlux_.dimensions() == dimMass/dimTime) { const volScalarField& rho = phi.db().objectRegistry::template lookupObject ("rho"); tUflux = this->faceFlux_/fvc::interpolate(rho); } else if (this->faceFlux_.dimensions() != dimVolume/dimTime) { FatalErrorInFunction << "dimensions of faceFlux are not correct" << exit(FatalError); } const surfaceScalarField& Uflux = tUflux(); scalarField& pLimiter = Limiter.primitiveFieldRef(); forAll(pLimiter, face) { pLimiter[face] = PhiLimiter::limiter ( CDweights[face], Uflux[face], phi[owner[face]], phi[neighbour[face]], Sf[face], magSf[face] ); } surfaceScalarField::Boundary& bLimiter = Limiter.boundaryFieldRef(); forAll(bLimiter, patchi) { scalarField& pLimiter = bLimiter[patchi]; if (bLimiter[patchi].coupled()) { const scalarField& pCDweights = CDweights.boundaryField()[patchi]; const vectorField& pSf = Sf.boundaryField()[patchi]; const scalarField& pmagSf = magSf.boundaryField()[patchi]; const scalarField& pFaceFlux = Uflux.boundaryField()[patchi]; const Field pphiP ( phi.boundaryField()[patchi].patchInternalField() ); const Field pphiN ( phi.boundaryField()[patchi].patchNeighbourField() ); forAll(pLimiter, face) { pLimiter[face] = PhiLimiter::limiter ( pCDweights[face], pFaceFlux[face], pphiP[face], pphiN[face], pSf[face], pmagSf[face] ); } } else { pLimiter = 1.0; } } return tLimiter; } // ************************************************************************* //