diff --git a/applications/utilities/parallelProcessing/redistributeMeshPar/redistributeMeshPar.C b/applications/utilities/parallelProcessing/redistributeMeshPar/redistributeMeshPar.C
index 7bddb26d7cb710947919f2545011bde6d5b88c10..a76cbb1393452173e1a34e5394efc1ca58fa428b 100644
--- a/applications/utilities/parallelProcessing/redistributeMeshPar/redistributeMeshPar.C
+++ b/applications/utilities/parallelProcessing/redistributeMeshPar/redistributeMeshPar.C
@@ -53,6 +53,7 @@ Description
 #include "fvMeshDistribute.H"
 #include "mapDistributePolyMesh.H"
 #include "IOobjectList.H"
+#include "globalIndex.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
@@ -107,7 +108,9 @@ autoPtr<fvMesh> createMesh
     fvMesh& mesh = meshPtr();
 
 
-    // Determine patches.
+    // Sync patches
+    // ~~~~~~~~~~~~
+
     if (Pstream::master())
     {
         // Send patches
@@ -118,14 +121,14 @@ autoPtr<fvMesh> createMesh
             slave++
         )
         {
-            OPstream toSlave(Pstream::blocking, slave);
+            OPstream toSlave(Pstream::scheduled, slave);
             toSlave << mesh.boundaryMesh();
         }
     }
     else
     {
         // Receive patches
-        IPstream fromMaster(Pstream::blocking, Pstream::masterNo());
+        IPstream fromMaster(Pstream::scheduled, Pstream::masterNo());
         PtrList<entry> patchEntries(fromMaster);
 
         if (haveMesh)
@@ -224,6 +227,58 @@ autoPtr<fvMesh> createMesh
         }
     }
 
+
+    // Determine zones
+    // ~~~~~~~~~~~~~~~
+
+    wordList pointZoneNames(mesh.pointZones().names());
+    Pstream::scatter(pointZoneNames);
+    wordList faceZoneNames(mesh.faceZones().names());
+    Pstream::scatter(faceZoneNames);
+    wordList cellZoneNames(mesh.cellZones().names());
+    Pstream::scatter(cellZoneNames);
+
+    if (!haveMesh)
+    {
+        // Add the zones
+        List<pointZone*> pz(pointZoneNames.size());
+        forAll(pointZoneNames, i)
+        {
+            pz[i] = new pointZone
+            (
+                pointZoneNames[i],
+                labelList(0),
+                i,
+                mesh.pointZones()
+            );
+        }
+        List<faceZone*> fz(faceZoneNames.size());
+        forAll(faceZoneNames, i)
+        {
+            fz[i] = new faceZone
+            (
+                faceZoneNames[i],
+                labelList(0),
+                boolList(0),
+                i,
+                mesh.faceZones()
+            );
+        }
+        List<cellZone*> cz(cellZoneNames.size());
+        forAll(cellZoneNames, i)
+        {
+            cz[i] = new cellZone
+            (
+                cellZoneNames[i],
+                labelList(0),
+                i,
+                mesh.cellZones()
+            );
+        }
+        mesh.addZones(pz, fz, cz);        
+    }
+
+
     if (!haveMesh)
     {
         // We created a dummy mesh file above. Delete it.
@@ -236,6 +291,21 @@ autoPtr<fvMesh> createMesh
     mesh.clearOut();
     mesh.globalData();
 
+
+    // Do some checks.
+
+    // Check if the boundary definition is unique
+    mesh.boundaryMesh().checkDefinition(true);
+    // Check if the boundary processor patches are correct
+    mesh.boundaryMesh().checkParallelSync(true);
+    // Check names of zones are equal
+    mesh.cellZones().checkDefinition(true);
+    mesh.cellZones().checkParallelSync(true);
+    mesh.faceZones().checkDefinition(true);
+    mesh.faceZones().checkParallelSync(true);
+    mesh.pointZones().checkDefinition(true);
+    mesh.pointZones().checkParallelSync(true);
+
     return meshPtr;
 }
 
@@ -292,6 +362,59 @@ void printMeshData(Ostream& os, const polyMesh& mesh)
         << "          face zones:       " << mesh.faceZones().size() << nl
         << "          cell zones:       " << mesh.cellZones().size() << nl;
 }
+void printMeshData(const polyMesh& mesh)
+{
+    // Collect all data on master
+
+    globalIndex globalCells(mesh.nCells());
+    labelListList patchNeiProcNo(Pstream::nProcs());
+    labelListList patchSize(Pstream::nProcs());
+    const labelList& pPatches = mesh.globalData().processorPatches();
+    patchNeiProcNo[Pstream::myProcNo()].setSize(pPatches.size());
+    patchSize[Pstream::myProcNo()].setSize(pPatches.size());
+    forAll(pPatches, i)
+    {
+        const processorPolyPatch& ppp = refCast<const processorPolyPatch>
+        (
+            mesh.boundaryMesh()[pPatches[i]]
+        );
+        patchNeiProcNo[Pstream::myProcNo()][i] = ppp.neighbProcNo();
+        patchSize[Pstream::myProcNo()][i] = ppp.size();
+    }
+    Pstream::gatherList(patchNeiProcNo);
+    Pstream::gatherList(patchSize);
+
+
+    // Print stats
+
+    globalIndex globalBoundaryFaces(mesh.nFaces()-mesh.nInternalFaces());
+
+    for (label procI = 0; procI < Pstream::nProcs(); procI++)
+    {
+        Info<< endl
+            << "Processor " << procI << nl
+            << "    Number of cells = " << globalCells.localSize(procI)
+            << endl;
+
+        label nProcFaces = 0;
+
+        const labelList& nei = patchNeiProcNo[procI];
+        
+        forAll(patchNeiProcNo[procI], i)
+        {
+            Info<< "    Number of faces shared with processor "
+                << patchNeiProcNo[procI][i] << " = " << patchSize[procI][i]
+                << endl;
+
+            nProcFaces += patchSize[procI][i];
+        }
+
+        Info<< "    Number of processor patches = " << nei.size() << nl
+            << "    Number of processor faces = " << nProcFaces << nl
+            << "    Number of boundary faces = "
+            << globalBoundaryFaces.localSize(procI) << endl;
+    }
+}
 
 
 // Debugging: write volScalarField with decomposition for post processing.
@@ -507,6 +630,7 @@ void compareFields
 int main(int argc, char *argv[])
 {
 #   include "addRegionOption.H"
+#   include "addOverwriteOption.H"
     argList::addOption
     (
         "mergeTol",
@@ -539,6 +663,8 @@ int main(int argc, char *argv[])
     }
     Info<< "Using mesh subdirectory " << meshSubDir << nl << endl;
 
+    const bool overwrite = args.optionFound("overwrite");
+
 
     // Get time instance directory. Since not all processors have meshes
     // just use the master one everywhere.
@@ -573,9 +699,11 @@ int main(int argc, char *argv[])
     );
     fvMesh& mesh = meshPtr();
 
-    Pout<< "Read mesh:" << endl;
-    printMeshData(Pout, mesh);
-    Pout<< endl;
+    //Pout<< "Read mesh:" << endl;
+    //printMeshData(Pout, mesh);
+    //Pout<< endl;
+
+
 
     IOdictionary decompositionDict
     (
@@ -618,7 +746,10 @@ int main(int argc, char *argv[])
     }
 
     // Dump decomposition to volScalarField
-    writeDecomposition("decomposition", mesh, finalDecomp);
+    if (!overwrite)
+    {
+        writeDecomposition("decomposition", mesh, finalDecomp);
+    }
 
 
     // Create 0 sized mesh to do all the generation of zero sized
@@ -796,12 +927,20 @@ int main(int argc, char *argv[])
 
 
     // Print a bit
-    Pout<< "After distribution mesh:" << endl;
-    printMeshData(Pout, mesh);
-    Pout<< endl;
+    // Print some statistics
+    Info<< "After distribution:" << endl;
+    printMeshData(mesh);
+
 
-    runTime++;
-    Pout<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
+    if (!overwrite)
+    {
+        runTime++;
+    }
+    else
+    {
+        mesh.setInstance(masterInstDir);
+    }
+    Info<< "Writing redistributed mesh to " << runTime.timeName() << nl << endl;
     mesh.write();
 
 
diff --git a/tutorials/incompressible/simpleFoam/windTurbineTerrain/Allrun b/tutorials/incompressible/simpleFoam/windTurbineTerrain/Allrun
index c69808608c7d4af7432036135d885e3d6e553a99..48a73959b0288e498d0e1ef54dc9d7ceece0d04c 100755
--- a/tutorials/incompressible/simpleFoam/windTurbineTerrain/Allrun
+++ b/tutorials/incompressible/simpleFoam/windTurbineTerrain/Allrun
@@ -18,13 +18,22 @@ runApplication decomposePar
 cp system/decomposeParDict-par system/decomposeParDict
 runParallel snappyHexMesh 2 -overwrite
 
-# Add wildcard entries for meshes patches since not preserved
+# Add wildcard entries for meshed patches since not preserved
 # by decomposePar. Notice -literalRE option to add wildcard itself
 # without evaluation.
 runParallel changeDictionary 2 -literalRE
 
-runParallel setSet 2 -batch makeZones
-runParallel windSimpleFoam 2
+cp system/decomposeParDict-4proc system/decomposeParDict
+runParallel redistributeMeshPar 4 -overwrite
+runParallel renumberMesh 4 -overwrite
+
+# Add wildcard entries for meshes patches since not preserved
+# by decomposePar. Notice -literalRE option to add wildcard itself
+# without evaluation.
+#runParallel changeDictionary 4 -literalRE
+
+runParallel setSet 4 -batch makeZones
+runParallel windSimpleFoam 4
 
 runApplication reconstructParMesh -constant
 runApplication reconstructPar
diff --git a/tutorials/incompressible/simpleFoam/windTurbineTerrain/system/decomposeParDict-4proc b/tutorials/incompressible/simpleFoam/windTurbineTerrain/system/decomposeParDict-4proc
new file mode 100644
index 0000000000000000000000000000000000000000..d81e9121240067cf798e5ee552fe70e2e46f6348
--- /dev/null
+++ b/tutorials/incompressible/simpleFoam/windTurbineTerrain/system/decomposeParDict-4proc
@@ -0,0 +1,22 @@
+/*--------------------------------*- C++ -*----------------------------------*\
+| =========                 |                                                 |
+| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
+|  \\    /   O peration     | Version:  dev                                   |
+|   \\  /    A nd           | Web:      www.OpenFOAM.com                      |
+|    \\/     M anipulation  |                                                 |
+\*---------------------------------------------------------------------------*/
+FoamFile
+{
+    version     2.0;
+    format      ascii;
+    class       dictionary;
+    object      decomposeParDict;
+}
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+numberOfSubdomains 4;
+
+method          ptscotch;
+
+// ************************************************************************* //