openfoam/src/finiteVolume/interpolation/surfaceInterpolation/limitedSchemes/LimitedScheme/LimitedScheme.C
Mark Olesen 092db087c9 ENH: use tmp field factory methods [2] (#2723)
- src/finiteVolume, src/finiteArea
2024-02-21 14:31:39 +01:00

202 lines
5.8 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / 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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "volFields.H"
#include "surfaceFields.H"
#include "fvcGrad.H"
#include "coupledFvPatchFields.H"
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
template<class Type, class Limiter, template<class> class LimitFunc>
void Foam::LimitedScheme<Type, Limiter, LimitFunc>::calcLimiter
(
const GeometricField<Type, fvPatchField, volMesh>& phi,
surfaceScalarField& limiterField
) const
{
typedef GeometricField<typename Limiter::phiType, fvPatchField, volMesh>
VolFieldType;
typedef GeometricField<typename Limiter::gradPhiType, fvPatchField, volMesh>
GradVolFieldType;
const fvMesh& mesh = this->mesh();
tmp<VolFieldType> tlPhi = LimitFunc<Type>()(phi);
const VolFieldType& lPhi = tlPhi();
tmp<GradVolFieldType> 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<typename Limiter::phiType> plPhiP
(
lPhi.boundaryField()[patchi].patchInternalField()
);
const Field<typename Limiter::phiType> plPhiN
(
lPhi.boundaryField()[patchi].patchNeighbourField()
);
const Field<typename Limiter::gradPhiType> pGradcP
(
gradc.boundaryField()[patchi].patchInternalField()
);
const Field<typename Limiter::gradPhiType> 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 Type, class Limiter, template<class> class LimitFunc>
Foam::tmp<Foam::surfaceScalarField>
Foam::LimitedScheme<Type, Limiter, LimitFunc>::limiter
(
const GeometricField<Type, fvPatchField, volMesh>& phi
) const
{
const fvMesh& mesh = this->mesh();
const word limiterFieldName(type() + "Limiter(" + phi.name() + ')');
if (this->mesh().cache("limiter"))
{
auto* fldptr = mesh.getObjectPtr<surfaceScalarField>(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<surfaceScalarField>::New
(
limiterFieldName,
limiterField
);
}
else
{
auto tlimiterField = surfaceScalarField::New
(
limiterFieldName,
mesh,
dimless
);
calcLimiter(phi, tlimiterField.ref());
return tlimiterField;
}
}
// ************************************************************************* //