diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files
index e7023dea21..c89ad9ea5a 100644
--- a/src/finiteVolume/Make/files
+++ b/src/finiteVolume/Make/files
@@ -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
diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/deferredCorrection/deferredCorrection.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/deferredCorrection/deferredCorrection.C
new file mode 100644
index 0000000000..4d0975eab5
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/deferredCorrection/deferredCorrection.C
@@ -0,0 +1,71 @@
+/*---------------------------------------------------------------------------*\
+ ========= |
+ \\ / 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 .
+
+\*---------------------------------------------------------------------------*/
+
+#include "deferredCorrection.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+template
+Foam::tmp>
+Foam::deferredCorrection::correction
+(
+ const GeometricField& vf
+) const
+{
+ tmp> tsfCorr
+ (
+ new GeometricField
+ (
+ "deferredCorrection::correction(" + vf.name() + ')',
+ tbaseScheme_().interpolate(vf)
+ )
+ );
+ GeometricField& sfCorr = tsfCorr.ref();
+
+ // Interpolate using the upwind weights to avoid circular reference to
+ // [this] explicit correction
+ sfCorr -= upwind::interpolate(vf, upwind::weights());
+/*
+ auto& sfCorrBf = sfCorr.boundaryFieldRef();
+ for (auto& pf : sfCorrBf)
+ {
+ if (!pf.coupled())
+ {
+ pf = pTraits::zero;
+ }
+ }
+*/
+ return tsfCorr;
+}
+
+
+namespace Foam
+{
+ //makeSurfaceInterpolationScheme(deferredCorrection)
+ makeSurfaceInterpolationTypeScheme(deferredCorrection, scalar)
+ makeSurfaceInterpolationTypeScheme(deferredCorrection, vector)
+}
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/deferredCorrection/deferredCorrection.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/deferredCorrection/deferredCorrection.H
new file mode 100644
index 0000000000..af69763d86
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/deferredCorrection/deferredCorrection.H
@@ -0,0 +1,157 @@
+/*---------------------------------------------------------------------------*\
+ ========= |
+ \\ / 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 .
+
+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 deferredCorrection
+:
+ public upwind
+{
+ // Private data
+
+ //- Base scheme
+ tmp> 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(mesh, is),
+ tbaseScheme_
+ (
+ surfaceInterpolationScheme::New(mesh, is)
+ )
+ {}
+
+ //- Construct from faceFlux and Istream
+ deferredCorrection
+ (
+ const fvMesh& mesh,
+ const surfaceScalarField& faceFlux,
+ Istream& is
+ )
+ :
+ upwind(mesh, faceFlux, is),
+ tbaseScheme_
+ (
+ surfaceInterpolationScheme::New(mesh, faceFlux, is)
+ )
+ {}
+
+
+ // Member Functions
+
+ using upwind::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>
+ correction
+ (
+ const GeometricField&
+ ) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //