diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files
index 36b69f77f089d7a0955c1a142007e6008eed04fa..b0cb2bc25b20e82d180823eb595a38040ea785b0 100644
--- a/src/finiteVolume/Make/files
+++ b/src/finiteVolume/Make/files
@@ -363,6 +363,7 @@ $(schemes)/localBlended/localBlended.C
 $(schemes)/limiterBlended/limiterBlended.C
 $(schemes)/CoBlended/CoBlended.C
 $(schemes)/cellCoBlended/cellCoBlended.C
+$(schemes)/zoneBlended/zoneBlended.C
 $(schemes)/localMax/localMax.C
 $(schemes)/localMin/localMin.C
 
diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlended.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlended.C
new file mode 100644
index 0000000000000000000000000000000000000000..f8ac1722f7369048b6635b7bddb275d3110e284a
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlended.C
@@ -0,0 +1,39 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 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 "fvMesh.H"
+#include "zoneBlended.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+    makeSurfaceInterpolationScheme(zoneBlended);
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlended.H b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlended.H
new file mode 100644
index 0000000000000000000000000000000000000000..49a2cb6c5eab2bdffe324348f0744c39b452f45b
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlended.H
@@ -0,0 +1,370 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 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/>.
+
+Class
+    Foam::zoneBlended
+
+Group
+    grpFvSurfaceInterpolationSchemes
+
+Description
+    Multi-faceZone based blending differencing scheme.
+
+    Schemes are set in dictonary format according to:
+
+    \verbatim
+    divSchemes
+    {
+        .
+        .
+        div(phi,U)      Gauss zoneBlended
+        {
+            default         defaultScheme;
+            faceZone1       scheme1;
+            faceZone2       scheme2;
+            ...
+            faceZoneN       schemeN;
+        }
+        .
+        .
+    }
+    \endverbatim
+
+    The \c default entry specifies the background scheme; additional schemes
+    can be set per \c faceZone, e.g. \c scheme1 is applied to \c facZone1, 
+    \c scheme2 is applied to \c facZone2 etc. 
+
+
+Usage
+    Example of the \c zoneBlended scheme to use \c linearUpwind as the
+    background scheme and \c upwind in \c faceZone1:
+
+    \verbatim
+    divSchemes
+    {
+        .
+        .
+        div(phi,U)      Gauss zoneBlended 
+        {
+            default         linearUpwind grad(U); 
+            faceZone1       upwind;
+        };
+        .
+        .
+    }
+    \endverbatim
+
+SourceFiles
+    zoneBlended.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef Foam_zoneBlended_H
+#define Foam_zoneBlended_H
+
+#include "surfaceInterpolationScheme.H"
+#include "blendedSchemeBase.H"
+#include "surfaceInterpolate.H"
+#include "UIndirectList.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                         Class zoneBlended Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class zoneBlended
+:
+    public surfaceInterpolationScheme<Type>,
+    public blendedSchemeBase<Type>
+{
+    using SurfaceField = GeometricField<Type, fvsPatchField, surfaceMesh>;
+    using VolumeField = GeometricField<Type, fvPatchField, volMesh>;
+
+    // Private data
+
+        //- Face zones
+        wordList zoneNames_;
+
+        //- Schemes
+        //  Note: 0 index used to describe default/background scheme
+        List<tmp<surfaceInterpolationScheme<Type>>> schemePtrs_;
+
+        //- Corrected flag - true if any of the schemes has corrected() set
+        bool corrected_;
+
+
+    // Private Member Functions
+
+        //- Set the lists of face zone names and schemes
+        void setSchemes(const fvMesh& mesh, const dictionary& dict);
+
+        //- Set the lists of face zone names and schemes
+        void setSchemes
+        (
+            const fvMesh& mesh,
+            const surfaceScalarField& faceFlux,
+            const dictionary& dict
+        );
+
+        //- Retrieve a scheme from the list
+        const surfaceInterpolationScheme<Type>& scheme
+        (
+            const label schemei
+        ) const
+        {
+            return schemePtrs_[schemei]();
+        }
+
+        //- Set destination values from source values for a face zone
+        template<class FieldType>
+        void setFaceZoneValues
+        (
+            FieldType& dest,
+            const FieldType& src,
+            const faceZone& fz
+        ) const;
+
+        //- Set destination values to zero for a face zone
+        template<class FieldType>
+        void zeroFaceZoneValues(FieldType& dest, const faceZone& fz) const;
+
+        //- No copy construct
+        zoneBlended(const zoneBlended&) = delete;
+
+        //- No copy assignment
+        void operator=(const zoneBlended&) = delete;
+
+
+public:
+
+    //- Runtime type information
+    TypeName("zoneBlended");
+
+
+    // Constructors
+
+        //- Construct from mesh and Istream.
+        //  The name of the flux field is read from the Istream and looked-up
+        //  from the mesh objectRegistry
+        zoneBlended(const fvMesh& mesh, Istream& is)
+        :
+            surfaceInterpolationScheme<Type>(mesh),
+            zoneNames_(),
+            schemePtrs_(),
+            corrected_(false)
+        {
+            const dictionary dict(is);
+
+            setSchemes(mesh, dict);
+        }
+
+
+        //- Construct from mesh, faceFlux and Istream
+        zoneBlended
+        (
+            const fvMesh& mesh,
+            const surfaceScalarField& faceFlux,
+            Istream& is
+        )
+        :
+            surfaceInterpolationScheme<Type>(mesh),
+            zoneNames_(),
+            schemePtrs_(),
+            corrected_(false)
+        {
+            const dictionary dict(is);
+
+            setSchemes(mesh, faceFlux, dict);
+        }
+
+
+    // Member Functions
+
+        //- Return the face-based blending factor
+        virtual tmp<surfaceScalarField> blendingFactor
+        (
+            const VolumeField& vf
+        ) const
+        {
+            const auto& mesh = vf.mesh();
+            auto tbf = surfaceScalarField::New("blendingFactor", mesh, dimless);
+            auto& bf = tbf.ref();
+            auto& bbf = bf.boundaryFieldRef();
+            bf = 0.0;
+
+            const auto& pbm = mesh.boundaryMesh();
+            const auto& zones = mesh.faceZones();
+
+            // Use blending factor to show different zones
+            for (label zonei=1; zonei<zoneNames_.size(); ++zonei)
+            {
+                const word& name = zoneNames_[zonei];
+                const auto& fz = zones[name];
+
+                for (const label facei : fz)
+                {
+                    if (mesh.isInternalFace(facei))
+                    {
+                        bf[facei] = zonei;
+                    }
+                    else
+                    {
+                        const labelPair pf = pbm.whichPatchFace(facei);
+                        auto& pbf = bbf[pf.first()];
+                        if (pbf.size())
+                        {
+                            pbf[pf.second()] = zonei;
+                        }
+                    }
+                }
+            }
+
+            return tbf;
+        }
+
+
+        //- Return the interpolation weighting factors
+        tmp<surfaceScalarField> weights(const VolumeField& vf) const
+        {
+            const auto& mesh = vf.mesh();
+            auto tweights =
+                surfaceScalarField::New("weights", vf.mesh(), dimless);
+            auto& weights = tweights.ref();
+
+            // Set default scheme weights
+            weights = this->scheme(0).weights(vf);
+
+            // Set face zone weights
+            const auto& zones = mesh.faceZones();
+
+            for (label schemei=1; schemei<schemePtrs_.size(); ++schemei)
+            {
+                const auto& scheme = this->scheme(schemei);
+
+                auto tschemeWeights = scheme.weights(vf);
+                const auto& schemeWeights = tschemeWeights();
+                const auto& fz = zones[zoneNames_[schemei]];
+
+                setFaceZoneValues(weights, schemeWeights, fz);
+            }
+
+            return tweights;
+        }
+
+
+        //- Return the face-interpolate of the given cell field
+        //  with explicit correction
+        tmp<SurfaceField> interpolate(const VolumeField& vf) const
+        {
+            return
+                surfaceInterpolationScheme<Type>::interpolate
+                (
+                    vf,
+                    weights(vf)
+                );
+        }
+
+
+        //- Return true if this scheme uses an explicit correction
+        virtual bool corrected() const
+        {
+            return corrected_;
+        }
+
+
+        //- Return the explicit correction to the face-interpolate
+        //- for the given field
+        virtual tmp<SurfaceField> correction
+        (
+            const VolumeField& vf
+        ) const
+        {
+            const auto& mesh = vf.mesh();
+            auto tcorr =
+                SurfaceField::New("correction", vf.mesh(), vf.dimensions());
+            auto& corr = tcorr.ref();
+            corr = dimensioned<Type>(vf.dimensions(), Zero);
+
+            // Set default scheme correction
+            const auto& scheme0 = this->scheme(0);
+            if (scheme0.corrected())
+            {
+                corr = scheme0.correction(vf);
+            }
+
+            // Only the default scheme exists - can exit early
+            if (schemePtrs_.size() == 1) return tcorr;
+
+            // Set correction field in faceZones
+            const auto& zones = mesh.faceZones();
+
+            for (label schemei=1; schemei<schemePtrs_.size(); ++schemei)
+            {
+                const auto& scheme = this->scheme(schemei);
+
+                if (scheme.corrected())
+                {
+                    auto tschemeCorr = scheme.correction(vf);
+                    const auto& schemeCorr = tschemeCorr();
+                    const auto& fz = zones[zoneNames_[schemei]];
+
+                    setFaceZoneValues(corr, schemeCorr, fz);
+                }
+                else
+                {
+                    if (scheme0.corrected())
+                    {
+                        // Remove correction from base scheme face zone faces
+                        const auto& fz = zones[zoneNames_[schemei]];
+                        zeroFaceZoneValues(corr, fz);
+                    }
+                }
+            }
+
+            return tcorr;
+        }
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "zoneBlendedTemplates.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlendedTemplates.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlendedTemplates.C
new file mode 100644
index 0000000000000000000000000000000000000000..687d01847e0d788a2ae653ed176c24b1230e7ddb
--- /dev/null
+++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/zoneBlended/zoneBlendedTemplates.C
@@ -0,0 +1,203 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 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/>.
+
+\*---------------------------------------------------------------------------*/
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::zoneBlended<Type>::setSchemes
+(
+    const fvMesh& mesh,
+    const dictionary& dict
+)
+{
+    zoneNames_.resize(dict.size());
+    schemePtrs_.resize(dict.size());
+
+    schemePtrs_[0] =
+        surfaceInterpolationScheme<Type>::New(mesh, dict.lookup("default"));
+
+    zoneNames_[0] = "default";
+
+    corrected_ = scheme(0).corrected();
+
+    label schemei = 1;
+
+    for (const auto& e : dict)
+    {
+        if (e.isDict())
+        {
+            FatalIOErrorInFunction(dict)
+                << "Entries must be given in zoneName-scheme pairs"
+                << abort(FatalIOError);
+        }
+
+        const word& key = e.keyword();
+
+        if (key == "default")
+        {
+            // Background scheme - already handled at index 0
+        }
+        else
+        {
+            zoneNames_[schemei] = key;
+
+            schemePtrs_[schemei] =
+                surfaceInterpolationScheme<Type>::New(mesh, e.stream());
+
+            corrected_ = corrected_ || scheme(schemei).corrected();
+
+            ++schemei;
+        }
+    }
+}
+
+
+template<class Type>
+void Foam::zoneBlended<Type>::setSchemes
+(
+    const fvMesh& mesh,
+    const surfaceScalarField& faceFlux,
+    const dictionary& dict
+)
+{
+    zoneNames_.resize(dict.size());
+    schemePtrs_.resize(dict.size());
+
+    schemePtrs_[0] =
+        surfaceInterpolationScheme<Type>::New
+        (
+            mesh,
+            faceFlux,
+            dict.lookup("default")
+        );
+
+    zoneNames_[0] = "default";
+
+    corrected_ = scheme(0).corrected();
+
+    label schemei = 1;
+
+    for (const auto& e : dict)
+    {
+        if (e.isDict())
+        {
+            FatalIOErrorInFunction(dict)
+                << "Entries must be given in faceZoneName-scheme pairs"
+                << abort(FatalIOError);
+        }
+
+        const word& key = e.keyword();
+
+        if (key == "default")
+        {
+            // Background scheme - already handled at index 0
+        }
+        else
+        {
+            zoneNames_[schemei] = key;
+
+            schemePtrs_[schemei] =
+                surfaceInterpolationScheme<Type>::New
+                (
+                    mesh,
+                    faceFlux,
+                    e.stream()
+                );
+
+            corrected_ = corrected_ || scheme(schemei).corrected();
+
+            ++schemei;
+        }
+    }
+}
+
+
+template<class Type>
+template<class FieldType>
+void Foam::zoneBlended<Type>::setFaceZoneValues
+(
+    FieldType& dest,
+    const FieldType& src,
+    const faceZone& fz
+) const
+{
+    const auto& mesh = dest.mesh();
+    const auto& pbm = mesh.boundaryMesh();
+    const auto& srcBf = src.boundaryField();
+    auto& destBf = dest.boundaryFieldRef();
+
+    for (const label facei : fz)
+    {
+        if (mesh.isInternalFace(facei))
+        {
+            dest[facei] = src[facei];
+        }
+        else
+        {
+            const labelPair pf = pbm.whichPatchFace(facei);
+            auto& pdest = destBf[pf.first()];
+            if (pdest.size())
+            {
+                pdest[pf.second()] = srcBf[pf.first()][pf.second()];
+            }
+        }
+    }
+}
+
+
+template<class Type>
+template<class FieldType>
+void Foam::zoneBlended<Type>::zeroFaceZoneValues
+(
+    FieldType& dest,
+    const faceZone& fz
+) const
+{
+    const auto& mesh = dest.mesh();
+    const auto& pbm = mesh.boundaryMesh();
+    auto& destBf = dest.boundaryFieldRef();
+
+    for (const label facei : fz)
+    {
+        if (mesh.isInternalFace(facei))
+        {
+            dest[facei] = pTraits<Type>::zero;
+        }
+        else
+        {
+            const labelPair pf = pbm.whichPatchFace(facei);
+            auto& pdest = destBf[pf.first()];
+            if (pdest.size())
+            {
+                pdest[pf.second()] = pTraits<Type>::zero;
+            }
+        }
+    }
+}
+
+
+// ************************************************************************* //