diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/pointLinear/pointLinear.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/pointLinear/pointLinear.C
new file mode 100644
index 0000000000..bd7d2aaab9
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/pointLinear/pointLinear.C
@@ -0,0 +1,111 @@
+/*---------------------------------------------------------------------------*\
+ ========= |
+ \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
+ \\ / O peration |
+ \\ / A nd | Copyright (C) 2004-2011 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 "pointLinear.H"
+#include "fvMesh.H"
+#include "volPointInterpolation.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+template
+Foam::tmp >
+Foam::pointLinear::
+correction
+(
+ const GeometricField& vf
+) const
+{
+ const fvMesh& mesh = this->mesh();
+
+ GeometricField pvf =
+ volPointInterpolation::New(mesh).interpolate(vf);
+
+ tmp > tsfCorr =
+ linearInterpolate(vf);
+
+ Field& sfCorr = tsfCorr().internalField();
+
+ const pointField& points = mesh.points();
+ const pointField& C = mesh.C().internalField();
+ const faceList& faces = mesh.faces();
+ const scalarField& w = mesh.weights().internalField();
+ const labelList& owner = mesh.owner();
+ const labelList& neighbour = mesh.neighbour();
+
+ forAll(sfCorr, facei)
+ {
+ point pi =
+ w[owner[facei]]*C[owner[facei]]
+ + (1.0 - w[owner[facei]])*C[neighbour[facei]];
+
+ scalar at = triangle
+ (
+ pi,
+ points[faces[facei][0]],
+ points[faces[facei][faces[facei].size()-1]]
+ ).mag();
+
+ scalar sumAt = at;
+ Type sumPsip = at*(1.0/3.0)*
+ (
+ sfCorr[facei]
+ + pvf[faces[facei][0]]
+ + pvf[faces[facei][faces[facei].size()-1]]
+ );
+
+ for (label pointi=1; pointi
+ (
+ pi,
+ points[faces[facei][pointi]],
+ points[faces[facei][pointi-1]]
+ ).mag();
+
+ sumAt += at;
+ sumPsip += at*(1.0/3.0)*
+ (
+ sfCorr[facei]
+ + pvf[faces[facei][pointi]]
+ + pvf[faces[facei][pointi-1]]
+ );
+
+ }
+
+ sfCorr[facei] = sumPsip/sumAt - sfCorr[facei];
+ }
+
+ tsfCorr().boundaryField() = pTraits::zero;
+
+ return tsfCorr;
+}
+
+
+namespace Foam
+{
+ makeSurfaceInterpolationScheme(pointLinear);
+}
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/pointLinear/pointLinear.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/pointLinear/pointLinear.H
new file mode 100644
index 0000000000..ef72c7316c
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/pointLinear/pointLinear.H
@@ -0,0 +1,129 @@
+/*---------------------------------------------------------------------------*\
+ ========= |
+ \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
+ \\ / O peration |
+ \\ / A nd | Copyright (C) 2011-2011 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
+ pointLinear
+
+Description
+ Face-point interpolation scheme class derived from linear and
+ returns linear weighting factors but also applies an explicit correction.
+
+ Uses volPointInterpolation to obtain the field values at the face-points.
+
+SourceFiles
+ pointLinear.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef pointLinear_H
+#define pointLinear_H
+
+#include "linear.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+ Class pointLinear Declaration
+\*---------------------------------------------------------------------------*/
+
+template
+class pointLinear
+:
+ public linear
+{
+ // Private Member Functions
+
+ //- Disallow default bitwise copy construct
+ pointLinear(const pointLinear&);
+
+ //- Disallow default bitwise assignment
+ void operator=(const pointLinear&);
+
+
+public:
+
+ //- Runtime type information
+ TypeName("pointLinear");
+
+
+ // Constructors
+
+ //- Construct from mesh
+ pointLinear(const fvMesh& mesh)
+ :
+ linear(mesh)
+ {}
+
+
+ //- Construct from mesh and Istream
+ pointLinear
+ (
+ const fvMesh& mesh,
+ Istream&
+ )
+ :
+ linear(mesh)
+ {}
+
+
+ //- Construct from mesh, faceFlux and Istream
+ pointLinear
+ (
+ const fvMesh& mesh,
+ const surfaceScalarField&,
+ Istream&
+ )
+ :
+ linear(mesh)
+ {}
+
+
+ // Member Functions
+
+ //- 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& vf
+ ) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //