diff --git a/applications/solvers/heatTransfer/chtMultiRegionFoam/coupledDerivedFvPatchFields/turbulentTemperatureRadCoupledMixed/turbulentTemperatureRadCoupledMixedFvPatchScalarField.C b/applications/solvers/heatTransfer/chtMultiRegionFoam/coupledDerivedFvPatchFields/turbulentTemperatureRadCoupledMixed/turbulentTemperatureRadCoupledMixedFvPatchScalarField.C
index e9bfaf8dfd57d7b37527876c39522a6b81116d8a..a137719a44c3ec5fb2c2396d90314e74e27d7723 100644
--- a/applications/solvers/heatTransfer/chtMultiRegionFoam/coupledDerivedFvPatchFields/turbulentTemperatureRadCoupledMixed/turbulentTemperatureRadCoupledMixedFvPatchScalarField.C
+++ b/applications/solvers/heatTransfer/chtMultiRegionFoam/coupledDerivedFvPatchFields/turbulentTemperatureRadCoupledMixed/turbulentTemperatureRadCoupledMixedFvPatchScalarField.C
@@ -191,15 +191,7 @@ void turbulentTemperatureRadCoupledMixedFvPatchScalarField::updateCoeffs()
     if (QrNbrName_ != "none")
     {
         QrNbr = nbrPatch.lookupPatchField<volScalarField, scalar>(QrNbrName_);
-        mapDistribute::distribute
-        (
-            Pstream::defaultCommsType,
-            distMap.schedule(),
-            distMap.constructSize(),
-            distMap.subMap(),           // what to send
-            distMap.constructMap(),     // what to receive
-            QrNbr
-        );
+        mpp.map().distribute(QrNbr);
     }
 
     scalarField alpha(KDeltaNbr - (Qr + QrNbr)/Tp);
diff --git a/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict b/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict
index 1a27956dbc5e54f91526642cdd89a47b2e2cb616..8ee2c4b74c09687f516c86e4bcf066945f99b55d 100644
--- a/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict
+++ b/applications/utilities/mesh/generation/snappyHexMesh/snappyHexMeshDict
@@ -270,7 +270,7 @@ addLayersControls
     //- If points get not extruded do nGrow layers of connected faces that are
     //  also not grown. This helps convergence of the layer addition process
     //  close to features.
-    // Note: changed(corrected) w.r.t 16x! (didn't do anything in 16x)
+    // Note: changed(corrected) w.r.t 17x! (didn't do anything in 17x)
     nGrow 0;
 
 
@@ -301,7 +301,7 @@ addLayersControls
     maxThicknessToMedialRatio 0.3;
 
     // Angle used to pick up medial axis points
-    // Note: changed(corrected) w.r.t 16x! 90 degrees corresponds to 130 in 16x.
+    // Note: changed(corrected) w.r.t 17x! 90 degrees corresponds to 130 in 17x.
     minMedianAxisAngle 90;
 
     // Create buffer region for new layer terminations
diff --git a/src/OpenFOAM/db/Time/Time.C b/src/OpenFOAM/db/Time/Time.C
index f7566f5f935f163ea9bfe9f849b63d7175b5dd9b..af3b5544d869ae79e6c1d3993825bf2bc49a5631 100644
--- a/src/OpenFOAM/db/Time/Time.C
+++ b/src/OpenFOAM/db/Time/Time.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
+    \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -86,18 +86,25 @@ void Foam::Time::adjustDeltaT()
             (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
         );
 
-        label nStepsToNextWrite = label(timeToNextWrite/deltaT_ - SMALL) + 1;
-        scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
+        scalar nSteps = timeToNextWrite/deltaT_ - SMALL;
 
-        // Control the increase of the time step to within a factor of 2
-        // and the decrease within a factor of 5.
-        if (newDeltaT >= deltaT_)
+        // For tiny deltaT the label can overflow!
+        if (nSteps < labelMax)
         {
-            deltaT_ = min(newDeltaT, 2.0*deltaT_);
-        }
-        else
-        {
-            deltaT_ = max(newDeltaT, 0.2*deltaT_);
+            label nStepsToNextWrite = label(nSteps) + 1;
+
+            scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
+
+            // Control the increase of the time step to within a factor of 2
+            // and the decrease within a factor of 5.
+            if (newDeltaT >= deltaT_)
+            {
+                deltaT_ = min(newDeltaT, 2.0*deltaT_);
+            }
+            else
+            {
+                deltaT_ = max(newDeltaT, 0.2*deltaT_);
+            }
         }
     }
 }
diff --git a/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.C b/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.C
index 34b2b091de9b35583c6184e5f303116ee6e138d6..b26b225690fb3c5e9e80c286cbac290e62a655ac 100644
--- a/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.C
+++ b/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
+    \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -30,9 +30,24 @@ License
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
-defineTypeNameAndDebug(Foam::coupledPolyPatch, 0);
+namespace Foam
+{
+    defineTypeNameAndDebug(coupledPolyPatch, 0);
+
+    scalar coupledPolyPatch::matchTol = 1E-3;
 
-Foam::scalar Foam::coupledPolyPatch::matchTol = 1E-3;
+    template<>
+    const char* NamedEnum<coupledPolyPatch::transformType, 4>::names[] =
+    {
+        "unknown",
+        "rotational",
+        "translational",
+        "noOrdering"
+    };
+
+    const NamedEnum<coupledPolyPatch::transformType, 4>
+        coupledPolyPatch::transformTypeNames;
+}
 
 
 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
@@ -204,12 +219,14 @@ void Foam::coupledPolyPatch::calcTransformTensors
     const vectorField& nf,
     const vectorField& nr,
     const scalarField& smallDist,
-    const scalar absTol
+    const scalar absTol,
+    const transformType transform
 ) const
 {
     if (debug)
     {
         Pout<< "coupledPolyPatch::calcTransformTensors : " << name() << endl
+            << "    transform:" << transformTypeNames[transform] << nl
             << "    (half)size:" << Cf.size() << nl
             << "    absTol:" << absTol << nl
             << "    smallDist min:" << min(smallDist) << nl
@@ -242,9 +259,16 @@ void Foam::coupledPolyPatch::calcTransformTensors
             Pout<< "    error:" << error << endl;
         }
 
-        if (sum(mag(nf & nr)) < Cf.size()-error)
+        if
+        (
+            transform == ROTATIONAL
+         || (
+                transform != TRANSLATIONAL
+             && (sum(mag(nf & nr)) < Cf.size()-error)
+            )
+        )
         {
-            // Rotation, no separation
+            // Type is rotation or unknown and normals not aligned
 
             // Assume per-face differing transformation, correct later
 
@@ -284,7 +308,7 @@ void Foam::coupledPolyPatch::calcTransformTensors
         }
         else
         {
-            // No rotation, possible separation
+            // Translational or (unknown and normals aligned)
 
             forwardT_.setSize(0);
             reverseT_.setSize(0);
diff --git a/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.H b/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.H
index 31fc9bcb1baf1099805c10127a3c9e1154228293..682e9962b62c3aeafa2c35a0576b69d2e7cafe41 100644
--- a/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.H
+++ b/src/OpenFOAM/meshes/polyMesh/polyPatches/basic/coupled/coupledPolyPatch.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
+    \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -53,6 +53,20 @@ class coupledPolyPatch
 :
     public polyPatch
 {
+public:
+
+    enum transformType
+    {
+        UNKNOWN,            // unspecified; automatic ordering
+        ROTATIONAL,         // rotation along coordinate axis
+        TRANSLATIONAL,      // translation
+        NOORDERING          // unspecified, no automatic ordering
+    };
+    static const NamedEnum<transformType, 4> transformTypeNames;
+
+
+private:
+
     // Private data
 
         //- offset (distance) vector from one side of the couple to the other
@@ -82,6 +96,8 @@ protected:
         //- Calculate the transformation tensors
         //  smallDist : matching distance per face
         //  absTol    : absolute error in normal
+        //  if transformType = unknown it first tries rotational, then
+        //  translational transform
         void calcTransformTensors
         (
             const vectorField& Cf,
@@ -89,7 +105,8 @@ protected:
             const vectorField& nf,
             const vectorField& nr,
             const scalarField& smallDist,
-            const scalar absTol = matchTol
+            const scalar absTol = matchTol,
+            const transformType = UNKNOWN
         ) const;
 
         //- Initialise the calculation of the patch geometry
diff --git a/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.C b/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.C
index a5039132f9f47c92e7c6171e378c9445ebe02953..03032fbb53780bd860f8b2039d7c3f40ab5cd65b 100644
--- a/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.C
+++ b/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.C
@@ -44,24 +44,8 @@ namespace Foam
 
     addToRunTimeSelectionTable(polyPatch, cyclicPolyPatch, word);
     addToRunTimeSelectionTable(polyPatch, cyclicPolyPatch, dictionary);
-
-    template<>
-    const char* Foam::NamedEnum
-    <
-        Foam::cyclicPolyPatch::transformType,
-        4
-    >::names[] =
-    {
-        "unknown",
-        "rotational",
-        "translational",
-        "noOrdering"
-    };
 }
 
-const Foam::NamedEnum<Foam::cyclicPolyPatch::transformType, 4>
-    Foam::cyclicPolyPatch::transformTypeNames;
-
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
@@ -271,7 +255,9 @@ void Foam::cyclicPolyPatch::calcTransforms
             static_cast<const pointField&>(half1Ctrs),
             half0Normals,
             half1Normals,
-            half0Tols
+            half0Tols,
+            matchTol,
+            transform_
         );
 
         if (transform_ == ROTATIONAL && !parallel() && forwardT().size() > 1)
diff --git a/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.H b/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.H
index dd0d14db86e6e53a3d036da906f5bcaee2ba2b34..68a008bdedd606dde76b16bf6cf162e180be8a4a 100644
--- a/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.H
+++ b/src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/cyclic/cyclicPolyPatch.H
@@ -64,20 +64,6 @@ class cyclicPolyPatch
 :
     public coupledPolyPatch
 {
-public:
-
-    enum transformType
-    {
-        UNKNOWN,            // unspecified; automatic ordering
-        ROTATIONAL,         // rotation along coordinate axis
-        TRANSLATIONAL,      // translation
-        NOORDERING          // unspecified, no automatic ordering
-    };
-    static const NamedEnum<transformType, 4> transformTypeNames;
-
-
-private:
-
     // Private data
 
         //- Name of other half
diff --git a/src/meshTools/directMapped/directMappedPolyPatch/directMappedPatchBase.C b/src/meshTools/directMapped/directMappedPolyPatch/directMappedPatchBase.C
index fee6f86ff465211bb95d767251c4c5e85073942f..544245faf57bcbb3e984dbb1308e516279e24b23 100644
--- a/src/meshTools/directMapped/directMappedPolyPatch/directMappedPatchBase.C
+++ b/src/meshTools/directMapped/directMappedPolyPatch/directMappedPatchBase.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
+    \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -709,16 +709,17 @@ Foam::directMappedPatchBase::directMappedPatchBase
     }
     else
     {
-        FatalErrorIn
+        FatalIOErrorIn
         (
             "directMappedPatchBase::directMappedPatchBase\n"
             "(\n"
             "    const polyPatch& pp,\n"
             "    const dictionary& dict\n"
-            ")\n"
+            ")\n",
+            dict
         )   << "Please supply the offsetMode as one of "
             << NamedEnum<offsetMode, 3>::words()
-            << exit(FatalError);
+            << exit(FatalIOError);
     }
 }
 
diff --git a/src/sampling/probes/patchProbes.C b/src/sampling/probes/patchProbes.C
index 7e545c9f96aa0315a726bd2cd981fe42108c97d6..d7ea026dce52f0c350dbd91b317dcf4d9c6728a6 100644
--- a/src/sampling/probes/patchProbes.C
+++ b/src/sampling/probes/patchProbes.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2010-2010 OpenCFD Ltd.
+    \\  /    A nd           | Copyright (C) 2010-2011 OpenCFD Ltd.
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -53,23 +53,15 @@ void Foam::patchProbes::findElements(const fvMesh& mesh)
     {
         const vector& sample = operator[](probeI);
         label faceI = meshSearchEngine.findNearestBoundaryFace(sample);
-        if (faceI == -1)
-        {
-            nearest[probeI].second().first() = Foam::sqr(GREAT);
-            nearest[probeI].second().second() = Pstream::myProcNo();
-        }
-        else
-        {
-            const point& fc = mesh.faceCentres()[faceI];
-            nearest[probeI].first() = pointIndexHit
-            (
-                true,
-                fc,
-                faceI
-            );
-            nearest[probeI].second().first() = magSqr(fc-sample);
-            nearest[probeI].second().second() = Pstream::myProcNo();
-        }
+        const point& fc = mesh.faceCentres()[faceI];
+        nearest[probeI].first() = pointIndexHit
+        (
+            true,
+            fc,
+            faceI
+        );
+        nearest[probeI].second().first() = magSqr(fc-sample);
+        nearest[probeI].second().second() = Pstream::myProcNo();
     }
 
 
@@ -92,27 +84,16 @@ void Foam::patchProbes::findElements(const fvMesh& mesh)
         }
     }
 
-
-
     // Check if all patchProbes have been found.
     forAll(nearest, sampleI)
     {
-        label localI = nearest[sampleI].first().index();
-
-        if (localI == -1)
+        label localI = -1;
+        if (nearest[sampleI].second().second() == Pstream::myProcNo())
         {
-             if (Pstream::master())
-             {
-                WarningIn("patchProbes::findElements()")
-                    << "Did not find location "
-                    <<  nearest[sampleI].second().first()
-                    << " in any cell. Skipping location." << endl;
-             }
-        }
-        else
-        {
-            elementList_[sampleI] = localI;
+            localI = nearest[sampleI].first().index();
         }
+
+        elementList_[sampleI] = localI;
     }
 }
 
diff --git a/src/thermophysicalModels/radiationModels/derivedFvPatchFields/radiationCoupledBase/radiationCoupledBase.C b/src/thermophysicalModels/radiationModels/derivedFvPatchFields/radiationCoupledBase/radiationCoupledBase.C
index eab4624c295b675920d3537e24a91d4f05109b1c..f8d034ca34963952df97d50bd540e15c05e87379 100644
--- a/src/thermophysicalModels/radiationModels/derivedFvPatchFields/radiationCoupledBase/radiationCoupledBase.C
+++ b/src/thermophysicalModels/radiationModels/derivedFvPatchFields/radiationCoupledBase/radiationCoupledBase.C
@@ -81,17 +81,18 @@ Foam::radiationCoupledBase::radiationCoupledBase
         {
             if (!isA<directMappedPatchBase>(patch_.patch()))
             {
-                FatalErrorIn
+                FatalIOErrorIn
                 (
                     "radiationCoupledBase::radiationCoupledBase\n"
                     "(\n"
                     "    const fvPatch& p,\n"
                     "    const dictionary& dict\n"
-                    ")\n"
+                    ")\n",
+                    dict
                 )   << "\n    patch type '" << patch_.type()
                     << "' not type '" << directMappedPatchBase::typeName << "'"
                     << "\n    for patch " << patch_.name()
-                    << exit(FatalError);
+                    << exit(FatalIOError);
             }
 
             emissivity_ = scalarField(patch_.size(), 0.0);
@@ -102,16 +103,17 @@ Foam::radiationCoupledBase::radiationCoupledBase
         {
             if(!dict.found("emissivity"))
             {
-                FatalErrorIn
+                FatalIOErrorIn
                 (
                     "radiationCoupledBase::radiationCoupledBase\n"
                     "(\n"
                     "    const fvPatch& p,\n"
                     "    const dictionary& dict\n"
-                    ")\n"
+                    ")\n",
+                    dict
                 )   << "\n    emissivity key does not exist for patch "
                     << patch_.name()
-                    << exit(FatalError);
+                    << exit(FatalIOError);
             }
             else
             {
diff --git a/tutorials/incompressible/simpleSRFFoam/mixer/constant/polyMesh/blockMeshDict b/tutorials/incompressible/simpleSRFFoam/mixer/constant/polyMesh/blockMeshDict
index 8decdddd8c7f612f2c5333bf819e418fc1b37134..204134e5b4e9b55f8b17ac3e773f2dcaf7986104 100644
--- a/tutorials/incompressible/simpleSRFFoam/mixer/constant/polyMesh/blockMeshDict
+++ b/tutorials/incompressible/simpleSRFFoam/mixer/constant/polyMesh/blockMeshDict
@@ -70,45 +70,80 @@ edges
       arc  23 19 (       0.034        0.094        2.000 )
 );
 
-patches
+boundary
 (
-      patch inlet
-      (
+    inlet
+    {
+        type patch;
+        faces
+        (
             (13 12 21 16)
             (14 13 16 18)
             (15 14 18 20)
             (17 22 12 13)
             (23 19 14 15)
-      )
-      patch outlet
-      (
+        );
+    }
+    outlet
+    {
+        type patch;
+        faces
+        (
             (1 4 9 0)
             (2 6 4 1)
             (3 8 6 2)
             (5 1 0 10)
             (11 3 2 7)
-      )
-      wall innerWall
-      (
+        );
+    }
+    innerWall
+    {
+        type wall;
+        faces
+        (
             (2 1 13 14)
             (5 10 22 17)
             (5 17 13 1)
             (11 7 19 23)
             (7 2 14 19)
-      )
-      wall outerWall
-      (
+        );
+    }
+    outerWall
+    {
+        type wall;
+        faces
+        (
             (4 16 21 9)
             (6 18 16 4)
             (8 20 18 6)
-      )
-      cyclic cyclic
-      (
+        );
+    }
+    cyclic_half0
+    {
+        type cyclic;
+        neighbourPatch cyclic_half1;
+        transform rotational;
+        rotationAxis (0 0 1);
+        rotationCentre (0 0 0);
+        faces
+        (
             (0 9 21 12)
             (10 0 12 22)
+        );
+    }
+    cyclic_half1
+    {
+        type cyclic;
+        neighbourPatch cyclic_half0;
+        transform rotational;
+        rotationAxis (0 0 1);
+        rotationCentre (0 0 0);
+        faces
+        (
             (3 15 20 8)
             (11 23 15 3)
-      )
+        );
+    }
 );
 
 mergeMatchPairs