diff --git a/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict b/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict
index 074d5b8ee3405f62c8dd80b1d068a22a54d5a0c6..25405d2ff3e63ecc8972f225d6921bc44192de94 100644
--- a/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict
+++ b/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict
@@ -154,6 +154,7 @@ castellatedMeshControls
             //faceZone sphere;
             //cellZone sphere;
             //cellZoneInside inside;  //outside/insidePoint
+
             //- Optional specification of what to do with faceZone faces:
             //      internal : keep them as internal faces (default)
             //      baffle   : create baffles from them. This gives more
diff --git a/applications/utilities/mesh/manipulation/topoSet/topoSetDict b/applications/utilities/mesh/manipulation/topoSet/topoSetDict
index 4c2180f21d489ff5b916f8bce73628f23443558c..1f58d59ddee437f838f92bd7905cd044ec58232c 100644
--- a/applications/utilities/mesh/manipulation/topoSet/topoSetDict
+++ b/applications/utilities/mesh/manipulation/topoSet/topoSetDict
@@ -238,6 +238,14 @@ FoamFile
 //        cos     0.01;       // Tolerance (max cos of angle)
 //    }
 //
+//    // Walk on faces in faceSet, starting from face nearest given position
+//    source  regionToFace;
+//    sourceInfo
+//    {
+//        set         f0;
+//        nearPoint   (0.1 0.1 0.005);
+//    }
+//
 //
 //
 // pointSet
diff --git a/applications/utilities/preProcessing/setFields/setFields.C b/applications/utilities/preProcessing/setFields/setFields.C
index dfde514992d607d36f62079b8c21e2f9b70f5f30..41c09b6053c7cb08aaa4c62f13fea1ad772dbb74 100644
--- a/applications/utilities/preProcessing/setFields/setFields.C
+++ b/applications/utilities/preProcessing/setFields/setFields.C
@@ -110,6 +110,9 @@ bool setCellFieldType
             "(const fvMesh& mesh, const labelList& selectedCells,"
             "Istream& fieldValueStream)"
         ) << "Field " << fieldName << " not found" << endl;
+
+        // Consume value
+        (void)pTraits<Type>(fieldValueStream);
     }
 
     return true;
@@ -286,6 +289,9 @@ bool setFaceFieldType
             "(const fvMesh& mesh, const labelList& selectedFaces,"
             "Istream& fieldValueStream)"
         ) << "Field " << fieldName << " not found" << endl;
+
+        // Consume value
+        (void)pTraits<Type>(fieldValueStream);
     }
 
     return true;
diff --git a/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C b/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C
index 5acc8dfa4c4925ef0998cc515394d4cbf14f018c..ad5a974ed70fd338dcd6add55b6c2a8429df6ea0 100644
--- a/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C
+++ b/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C
@@ -149,13 +149,8 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::checkPatches
     const scalar maxBoundsError = 0.05;
 
     // check bounds of source and target
-    boundBox bbSrc(srcPatch.points(), srcPatch.meshPoints());
-    reduce(bbSrc.min(), minOp<point>());
-    reduce(bbSrc.max(), maxOp<point>());
-
-    boundBox bbTgt(tgtPatch.points(), tgtPatch.meshPoints());
-    reduce(bbTgt.min(), minOp<point>());
-    reduce(bbTgt.max(), maxOp<point>());
+    boundBox bbSrc(srcPatch.points(), srcPatch.meshPoints(), true);
+    boundBox bbTgt(tgtPatch.points(), tgtPatch.meshPoints(), true);
 
     boundBox bbTgtInf(bbTgt);
     bbTgtInf.inflate(maxBoundsError);
diff --git a/src/meshTools/Make/files b/src/meshTools/Make/files
index 3329c2b28d1115a290a7cc0b4b11537fadc7cd11..4abc57669ebf212139abd71a21710d3b9fcbfc00 100644
--- a/src/meshTools/Make/files
+++ b/src/meshTools/Make/files
@@ -110,6 +110,8 @@ $(faceSources)/patchToFace/patchToFace.C
 $(faceSources)/boundaryToFace/boundaryToFace.C
 $(faceSources)/zoneToFace/zoneToFace.C
 $(faceSources)/boxToFace/boxToFace.C
+$(faceSources)/regionToFace/regionToFace.C
+$(faceSources)/regionToFace/patchEdgeFaceRegion.C
 
 pointSources = sets/pointSources
 $(pointSources)/labelToPoint/labelToPoint.C
diff --git a/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegion.C b/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegion.C
new file mode 100644
index 0000000000000000000000000000000000000000..91e9b4c331919fd0c0a6e6e4097d0801c828aacd
--- /dev/null
+++ b/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegion.C
@@ -0,0 +1,50 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 OpenFOAM Foundation
+     \\/     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 "patchEdgeFaceRegion.H"
+
+// * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
+
+Foam::Ostream& Foam::operator<<
+(
+    Foam::Ostream& os,
+    const Foam::patchEdgeFaceRegion& wDist
+)
+{
+    return os << wDist.region_;
+}
+
+
+Foam::Istream& Foam::operator>>
+(
+    Foam::Istream& is,
+    Foam::patchEdgeFaceRegion& wDist
+)
+{
+    return is >> wDist.region_;
+}
+
+
+// ************************************************************************* //
diff --git a/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegion.H b/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegion.H
new file mode 100644
index 0000000000000000000000000000000000000000..6c0b65f2d7525a2c6cdeb490978dac48dd88c824
--- /dev/null
+++ b/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegion.H
@@ -0,0 +1,193 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 OpenFOAM Foundation
+     \\/     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::patchEdgeFaceRegion
+
+Description
+    Transport of region for use in PatchEdgeFaceWave.
+
+    Set element to -2 to denote blocked.
+
+SourceFiles
+    patchEdgeFaceRegionI.H
+    patchEdgeFaceRegion.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef patchEdgeFaceRegion_H
+#define patchEdgeFaceRegion_H
+
+#include "point.H"
+#include "label.H"
+#include "scalar.H"
+#include "tensor.H"
+#include "indirectPrimitivePatch.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// Forward declaration of classes
+class polyPatch;
+class polyMesh;
+
+/*---------------------------------------------------------------------------*\
+                           Class patchEdgeFaceRegion Declaration
+\*---------------------------------------------------------------------------*/
+
+class patchEdgeFaceRegion
+{
+    // Private data
+
+        //- region
+        label region_;
+
+    // Private Member Functions
+
+        //- Combine current with w2. Update region_ if w2 has smaller
+        //  quantities and returns true.
+        template<class TrackingData>
+        inline bool update
+        (
+            const patchEdgeFaceRegion& w2,
+            const scalar tol,
+            TrackingData& td
+        );
+
+
+public:
+
+    // Constructors
+
+        //- Construct null
+        inline patchEdgeFaceRegion();
+
+        //- Construct from origin, distance
+        inline patchEdgeFaceRegion(const label);
+
+
+    // Member Functions
+
+        // Access
+
+            inline label region() const;
+
+
+        // Needed by meshWave
+
+            //- Check whether origin has been changed at all or
+            //  still contains original (invalid) value.
+            template<class TrackingData>
+            inline bool valid(TrackingData& td) const;
+
+            //- Apply rotation matrix
+            template<class TrackingData>
+            inline void transform
+            (
+                const polyMesh& mesh,
+                const indirectPrimitivePatch& patch,
+                const tensor& rotTensor,
+                const scalar tol,
+                TrackingData& td
+            );
+
+            //- Influence of face on edge
+            template<class TrackingData>
+            inline bool updateEdge
+            (
+                const polyMesh& mesh,
+                const indirectPrimitivePatch& patch,
+                const label edgeI,
+                const label faceI,
+                const patchEdgeFaceRegion& faceInfo,
+                const scalar tol,
+                TrackingData& td
+            );
+
+            //- New information for edge (from e.g. coupled edge)
+            template<class TrackingData>
+            inline bool updateEdge
+            (
+                const polyMesh& mesh,
+                const indirectPrimitivePatch& patch,
+                const patchEdgeFaceRegion& edgeInfo,
+                const bool sameOrientation,
+                const scalar tol,
+                TrackingData& td
+            );
+
+            //- Influence of edge on face.
+            template<class TrackingData>
+            inline bool updateFace
+            (
+                const polyMesh& mesh,
+                const indirectPrimitivePatch& patch,
+                const label faceI,
+                const label edgeI,
+                const patchEdgeFaceRegion& edgeInfo,
+                const scalar tol,
+                TrackingData& td
+            );
+
+            //- Same (like operator==)
+            template<class TrackingData>
+            inline bool equal(const patchEdgeFaceRegion&, TrackingData&) const;
+
+
+    // Member Operators
+
+        // Needed for List IO
+        inline bool operator==(const patchEdgeFaceRegion&) const;
+        inline bool operator!=(const patchEdgeFaceRegion&) const;
+
+
+    // IOstream Operators
+
+        friend Ostream& operator<<(Ostream&, const patchEdgeFaceRegion&);
+        friend Istream& operator>>(Istream&, patchEdgeFaceRegion&);
+};
+
+
+//- Data associated with patchEdgeFaceRegion type are contiguous
+template<>
+inline bool contiguous<patchEdgeFaceRegion>()
+{
+    return true;
+}
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#include "patchEdgeFaceRegionI.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegionI.H b/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegionI.H
new file mode 100644
index 0000000000000000000000000000000000000000..60a6d13eb3f2797a2884162008f6841f697f4093
--- /dev/null
+++ b/src/meshTools/sets/faceSources/regionToFace/patchEdgeFaceRegionI.H
@@ -0,0 +1,209 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 OpenFOAM Foundation
+     \\/     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 "polyMesh.H"
+#include "transform.H"
+
+// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
+
+// Update this with w2 if w2 nearer to pt.
+template<class TrackingData>
+inline bool Foam::patchEdgeFaceRegion::update
+(
+    const patchEdgeFaceRegion& w2,
+    const scalar tol,
+    TrackingData& td
+)
+{
+    if (!w2.valid(td))
+    {
+        FatalErrorIn("patchEdgeFaceRegion::update(..)")
+            << "problem." << abort(FatalError);
+    }
+
+    if (w2.region_ == -2 || region_ == -2)
+    {
+
+Pout<< "update : " << *this << "   w2:" << w2 << " return FALSE" << endl;
+
+        // Blocked edge/face
+        return false;
+    }
+
+    if (!valid(td))
+    {
+        // current not yet set so use any value
+        label oldRegion = region_;
+        operator=(w2);
+Pout<< "update : " << *this << " was:" << oldRegion
+    << "   w2:" << w2 << " return TRUE" << endl;
+        return true;
+    }
+    else
+    {
+        if (w2.region_ < region_)
+        {
+            label oldRegion = region_;
+            operator=(w2);
+Pout<< "update : " << *this << " was:" << oldRegion
+    << "   w2:" << w2 << " return TRUE" << endl;
+        return true;
+            return true;
+        }
+        else
+        {
+Pout<< "update : " << *this
+    << "   w2:" << w2 << " return FALSE" << endl;
+            return false;
+        }
+    }
+}
+
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+// Null constructor
+inline Foam::patchEdgeFaceRegion::patchEdgeFaceRegion()
+:
+    region_(-1)
+{}
+
+
+// Construct from origin, distance
+inline Foam::patchEdgeFaceRegion::patchEdgeFaceRegion
+(
+    const label region
+)
+:
+    region_(region)
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+inline Foam::label Foam::patchEdgeFaceRegion::region() const
+{
+    return region_;
+}
+
+
+template<class TrackingData>
+inline bool Foam::patchEdgeFaceRegion::valid(TrackingData& td) const
+{
+    return region_ != -1;
+}
+
+
+template<class TrackingData>
+inline void Foam::patchEdgeFaceRegion::transform
+(
+    const polyMesh& mesh,
+    const indirectPrimitivePatch& patch,
+    const tensor& rotTensor,
+    const scalar tol,
+    TrackingData& td
+)
+{}
+
+
+template<class TrackingData>
+inline bool Foam::patchEdgeFaceRegion::updateEdge
+(
+    const polyMesh& mesh,
+    const indirectPrimitivePatch& patch,
+    const label edgeI,
+    const label faceI,
+    const patchEdgeFaceRegion& faceInfo,
+    const scalar tol,
+    TrackingData& td
+)
+{
+    return update(faceInfo, tol, td);
+}
+
+
+template<class TrackingData>
+inline bool Foam::patchEdgeFaceRegion::updateEdge
+(
+    const polyMesh& mesh,
+    const indirectPrimitivePatch& patch,
+    const patchEdgeFaceRegion& edgeInfo,
+    const bool sameOrientation,
+    const scalar tol,
+    TrackingData& td
+)
+{
+    return update(edgeInfo, tol, td);
+}
+
+
+template<class TrackingData>
+inline bool Foam::patchEdgeFaceRegion::updateFace
+(
+    const polyMesh& mesh,
+    const indirectPrimitivePatch& patch,
+    const label faceI,
+    const label edgeI,
+    const patchEdgeFaceRegion& edgeInfo,
+    const scalar tol,
+    TrackingData& td
+)
+{
+    return update(edgeInfo, tol, td);
+}
+
+
+template <class TrackingData>
+inline bool Foam::patchEdgeFaceRegion::equal
+(
+    const patchEdgeFaceRegion& rhs,
+    TrackingData& td
+) const
+{
+    return operator==(rhs);
+}
+
+
+// * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
+
+inline bool Foam::patchEdgeFaceRegion::operator==
+(
+    const Foam::patchEdgeFaceRegion& rhs
+) const
+{
+    return region() == rhs.region();
+}
+
+
+inline bool Foam::patchEdgeFaceRegion::operator!=
+(
+    const Foam::patchEdgeFaceRegion& rhs
+) const
+{
+    return !(*this == rhs);
+}
+
+
+// ************************************************************************* //
diff --git a/src/meshTools/sets/faceSources/regionToFace/regionToFace.C b/src/meshTools/sets/faceSources/regionToFace/regionToFace.C
new file mode 100644
index 0000000000000000000000000000000000000000..5c5f538876becc47012933bfb83376e11eb1d45f
--- /dev/null
+++ b/src/meshTools/sets/faceSources/regionToFace/regionToFace.C
@@ -0,0 +1,254 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 OpenFOAM Foundation
+     \\/     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 "regionToFace.H"
+#include "polyMesh.H"
+#include "faceSet.H"
+#include "mappedPatchBase.H"
+#include "indirectPrimitivePatch.H"
+#include "PatchTools.H"
+#include "addToRunTimeSelectionTable.H"
+#include "PatchEdgeFaceWave.H"
+#include "patchEdgeFaceRegion.H"
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+defineTypeNameAndDebug(regionToFace, 0);
+
+addToRunTimeSelectionTable(topoSetSource, regionToFace, word);
+
+addToRunTimeSelectionTable(topoSetSource, regionToFace, istream);
+
+}
+
+
+Foam::topoSetSource::addToUsageTable Foam::regionToFace::usage_
+(
+    regionToFace::typeName,
+    "\n    Usage: regionToFace <faceSet> (x y z)\n\n"
+    "    Select all faces in the connected region of the faceSet"
+    " starting from the point.\n"
+);
+
+
+// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
+
+void Foam::regionToFace::markZone
+(
+    const indirectPrimitivePatch& patch,
+    const label procI,
+    const label faceI,
+    const label zoneI,
+    labelList& faceZone
+) const
+{
+    // Data on all edges and faces
+    List<patchEdgeFaceRegion> allEdgeInfo(patch.nEdges());
+    List<patchEdgeFaceRegion> allFaceInfo(patch.size());
+
+    DynamicList<label> changedEdges;
+    DynamicList<patchEdgeFaceRegion> changedInfo;
+
+    if (Pstream::myProcNo() == procI)
+    {
+        const labelList& fEdges = patch.faceEdges()[faceI];
+        forAll(fEdges, i)
+        {
+            changedEdges.append(fEdges[i]);
+            changedInfo.append(zoneI);
+        }
+    }
+
+    // Walk
+    PatchEdgeFaceWave
+    <
+        indirectPrimitivePatch,
+        patchEdgeFaceRegion
+    > calc
+    (
+        mesh_,
+        patch,
+        changedEdges,
+        changedInfo,
+        allEdgeInfo,
+        allFaceInfo,
+        returnReduce(patch.nEdges(), sumOp<label>())
+    );
+
+    forAll(allFaceInfo, faceI)
+    {
+        if (allFaceInfo[faceI].region() == zoneI)
+        {
+            faceZone[faceI] = zoneI;
+        }
+    }
+}
+
+
+void Foam::regionToFace::combine(topoSet& set, const bool add) const
+{
+    Info<< "    Loading subset " << setName_ << " to delimit search region."
+        << endl;
+    faceSet subSet(mesh_, setName_);
+
+    indirectPrimitivePatch patch
+    (
+        IndirectList<face>(mesh_.faces(), subSet.toc()),
+        mesh_.points()
+    );
+
+    mappedPatchBase::nearInfo ni
+    (
+        pointIndexHit(false, vector::zero, -1),
+        Tuple2<scalar, label>
+        (
+            sqr(GREAT),
+            Pstream::myProcNo()
+        )
+    );
+
+    forAll(patch, i)
+    {
+        const point& fc = patch.faceCentres()[i];
+        scalar d2 = magSqr(fc-nearPoint_);
+
+        if (!ni.first().hit() || d2 < ni.second().first())
+        {
+            ni.second().first() = d2;
+            ni.first().setHit();
+            ni.first().setPoint(fc);
+            ni.first().setIndex(i);
+        }
+    }
+
+    // Globally reduce
+    combineReduce(ni, mappedPatchBase::nearestEqOp());
+
+    Info<< "    Found nearest face at " << ni.first().rawPoint()
+        << " on processor " << ni.second().second()
+        << " face " << ni.first().index()
+        << " distance " << Foam::sqrt(ni.second().first()) << endl;
+
+    labelList faceRegion(patch.size(), -1);
+    markZone
+    (
+        patch,
+        ni.second().second(),   // procI
+        ni.first().index(),     // start face
+        0,                      // currentZone
+        faceRegion
+    );
+
+    forAll(faceRegion, faceI)
+    {
+        if (faceRegion[faceI] == 0)
+        {
+            addOrDelete(set, patch.addressing()[faceI], add);
+        }
+    }
+}
+
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+// Construct from components
+Foam::regionToFace::regionToFace
+(
+    const polyMesh& mesh,
+    const word& setName,
+    const point& nearPoint
+)
+:
+    topoSetSource(mesh),
+    setName_(setName),
+    nearPoint_(nearPoint)
+{}
+
+
+// Construct from dictionary
+Foam::regionToFace::regionToFace
+(
+    const polyMesh& mesh,
+    const dictionary& dict
+)
+:
+    topoSetSource(mesh),
+    setName_(dict.lookup("set")),
+    nearPoint_(dict.lookup("nearPoint"))
+{}
+
+
+// Construct from Istream
+Foam::regionToFace::regionToFace
+(
+    const polyMesh& mesh,
+    Istream& is
+)
+:
+    topoSetSource(mesh),
+    setName_(checkIs(is)),
+    nearPoint_(checkIs(is))
+{}
+
+
+// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
+
+Foam::regionToFace::~regionToFace()
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+void Foam::regionToFace::applyToSet
+(
+    const topoSetSource::setAction action,
+    topoSet& set
+) const
+{
+    if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
+    {
+        Info<< "    Adding all faces of connected region of set "
+            << setName_
+            << " starting from point "
+            << nearPoint_ << " ..." << endl;
+
+        combine(set, true);
+    }
+    else if (action == topoSetSource::DELETE)
+    {
+        Info<< "    Removing all cells of connected region of set "
+            << setName_
+            << " starting from point "
+            << nearPoint_ << " ..." << endl;
+
+        combine(set, false);
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/src/meshTools/sets/faceSources/regionToFace/regionToFace.H b/src/meshTools/sets/faceSources/regionToFace/regionToFace.H
new file mode 100644
index 0000000000000000000000000000000000000000..1dc1d48069ff74a2311ed47de996498a5f971c7f
--- /dev/null
+++ b/src/meshTools/sets/faceSources/regionToFace/regionToFace.H
@@ -0,0 +1,139 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2012 OpenFOAM Foundation
+     \\/     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::regionToFace
+
+Description
+    A topoSetSource to select faces belonging to topological connected region
+    (that contains given point)
+
+SourceFiles
+    regionToFace.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef regionToFace_H
+#define regionToFace_H
+
+#include "topoSetSource.H"
+#include "PackedBoolList.H"
+#include "indirectPrimitivePatch.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                        Class regionToFace Declaration
+\*---------------------------------------------------------------------------*/
+
+class regionToFace
+:
+    public topoSetSource
+{
+    // Private data
+
+        //- Add usage string
+        static addToUsageTable usage_;
+
+        //- Name of set to use
+        word setName_;
+
+        //- Coordinate that is nearest/on connected region
+        point nearPoint_;
+
+    // Private Member Functions
+
+        //- Walk edge-face-edge
+        void markZone
+        (
+            const indirectPrimitivePatch& patch,
+            const label procI,
+            const label faceI,
+            const label zoneI,
+            labelList& faceZone
+        ) const;
+
+        void combine(topoSet& set, const bool add) const;
+
+public:
+
+    //- Runtime type information
+    TypeName("regionToFace");
+
+    // Constructors
+
+        //- Construct from components
+        regionToFace
+        (
+            const polyMesh& mesh,
+            const word& setName,
+            const point& nearPoint
+        );
+
+        //- Construct from dictionary
+        regionToFace
+        (
+            const polyMesh& mesh,
+            const dictionary& dict
+        );
+
+        //- Construct from Istream
+        regionToFace
+        (
+            const polyMesh& mesh,
+            Istream&
+        );
+
+
+    //- Destructor
+    virtual ~regionToFace();
+
+
+    // Member Functions
+
+        virtual sourceType setType() const
+        {
+            return FACESETSOURCE;
+        }
+
+        virtual void applyToSet
+        (
+            const topoSetSource::setAction action,
+            topoSet&
+        ) const;
+
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //