Commit 929ed3ca authored by Andrew Heather's avatar Andrew Heather
Browse files

ENH: New deferred correction scheme

The characteristics of the base scheme are recovered by applying an
explicit correction to the upwind scheme weights.

Usage
    Example of the \c deferredCorrection scheme applied to the \c linear
    scheme:
    \verbatim
    divSchemes
    {
        .
        .
        div(phi,U)      Gauss deferredCorrection linear;
        .
        .
    }
    \endverbatim

Based on a generalised form of a deferred correction linear scheme
supplied by CFD Software E+F GmbH
parent 90a9cd06
......@@ -310,6 +310,7 @@ $(schemes)/linearPureUpwindFit/linearPureUpwindFit.C
$(schemes)/linearUpwind/linearUpwind.C
$(schemes)/linearUpwind/linearUpwindV.C
$(schemes)/LUST/LUST.C
$(schemes)/deferredCorrection/deferredCorrection.C
limitedSchemes = $(surfaceInterpolation)/limitedSchemes
$(limitedSchemes)/limitedSurfaceInterpolationScheme/limitedSurfaceInterpolationSchemes.C
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
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 "deferredCorrection.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Foam::deferredCorrection<Type>::correction
(
const GeometricField<Type, fvPatchField, volMesh>& vf
) const
{
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsfCorr
(
new GeometricField<Type, fvsPatchField, surfaceMesh>
(
"deferredCorrection::correction(" + vf.name() + ')',
tbaseScheme_().interpolate(vf)
)
);
GeometricField<Type, fvsPatchField, surfaceMesh>& sfCorr = tsfCorr.ref();
// Interpolate using the upwind weights to avoid circular reference to
// [this] explicit correction
sfCorr -= upwind<Type>::interpolate(vf, upwind<Type>::weights());
/*
auto& sfCorrBf = sfCorr.boundaryFieldRef();
for (auto& pf : sfCorrBf)
{
if (!pf.coupled())
{
pf = pTraits<Type>::zero;
}
}
*/
return tsfCorr;
}
namespace Foam
{
//makeSurfaceInterpolationScheme(deferredCorrection)
makeSurfaceInterpolationTypeScheme(deferredCorrection, scalar)
makeSurfaceInterpolationTypeScheme(deferredCorrection, vector)
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
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/>.
Class
Foam::deferredCorrection
Description
Deferred correction interpolation scheme wrapper around a run-time
selectable base scheme.
The characteristics of the base scheme are recovered by applying an
explicit correction to the upwind scheme weights.
Usage
Example of the \c deferredCorrection scheme applied to the \c linear
scheme:
\verbatim
divSchemes
{
.
.
div(phi,U) Gauss deferredCorrection linear;
.
.
}
\endverbatim
SourceFiles
deferredCorrection.C
SeeAlso
Foam::upwind
\*---------------------------------------------------------------------------*/
#ifndef deferredCorrection_H
#define deferredCorrection_H
#include "upwind.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class deferredCorrection Declaration
\*---------------------------------------------------------------------------*/
template<class Type>
class deferredCorrection
:
public upwind<Type>
{
// Private data
//- Base scheme
tmp<surfaceInterpolationScheme<Type>> tbaseScheme_;
// Private Member Functions
//- Disallow default bitwise copy construct
deferredCorrection(const deferredCorrection&) = delete;
//- Disallow default bitwise assignment
void operator=(const deferredCorrection&) = delete;
public:
//- Runtime type information
TypeName("deferredCorrection");
// Constructors
//- Construct from Istream.
// The name of the flux field is read from the Istream and looked-up
// from the mesh objectRegistry
deferredCorrection
(
const fvMesh& mesh,
Istream& is
)
:
upwind<Type>(mesh, is),
tbaseScheme_
(
surfaceInterpolationScheme<Type>::New(mesh, is)
)
{}
//- Construct from faceFlux and Istream
deferredCorrection
(
const fvMesh& mesh,
const surfaceScalarField& faceFlux,
Istream& is
)
:
upwind<Type>(mesh, faceFlux, is),
tbaseScheme_
(
surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
)
{}
// Member Functions
using upwind<Type>::weights;
//- Return true if this scheme uses an explicit correction
virtual bool corrected() const
{
return true;
}
//- Return the explicit correction to the face-interpolate
virtual tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
correction
(
const GeometricField<Type, fvPatchField, volMesh>&
) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
Markdown is supported
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment