diff --git a/applications/utilities/mesh/generation/extrude/extrudeToRegionMesh/extrudeToRegionMesh.C b/applications/utilities/mesh/generation/extrude/extrudeToRegionMesh/extrudeToRegionMesh.C
index 0e84353209f248b1660cd0e7a9bf3c00deeb4180..44bd1b43674b69aa7d16b2bcccc3e36828582eb6 100644
--- a/applications/utilities/mesh/generation/extrude/extrudeToRegionMesh/extrudeToRegionMesh.C
+++ b/applications/utilities/mesh/generation/extrude/extrudeToRegionMesh/extrudeToRegionMesh.C
@@ -492,15 +492,8 @@ labelListList globalEdgeFaces
 
     forAll(edgeFaces, edgeI)
     {
-        const labelList& eFaces = edgeFaces[edgeI];
-
         // Store pp face and processor as unique tag.
-        labelList& globalEFaces = globalEdgeFaces[edgeI];
-        globalEFaces.setSize(eFaces.size());
-        forAll(eFaces, i)
-        {
-            globalEFaces[i] = globalFaces.toGlobal(eFaces[i]);
-        }
+        globalEdgeFaces[edgeI] = globalFaces.toGlobal(edgeFaces[edgeI]);
     }
 
     // Synchronise across coupled edges.
diff --git a/applications/utilities/mesh/manipulation/checkMesh/checkTools.C b/applications/utilities/mesh/manipulation/checkMesh/checkTools.C
index cf981adb581bd3150bfbcd881cd4f091a3543c2b..76e4fc5b18cdcfa602cdd962d64f16634284445b 100644
--- a/applications/utilities/mesh/manipulation/checkMesh/checkTools.C
+++ b/applications/utilities/mesh/manipulation/checkMesh/checkTools.C
@@ -416,11 +416,7 @@ void Foam::mergeAndWrite
 
         // Get renumbered local data
         pointField myPoints(mesh.points(), setPointIDs);
-        labelList myIDs(setPointIDs.size());
-        forAll(setPointIDs, i)
-        {
-            myIDs[i] = globalNumbering.toGlobal(setPointIDs[i]);
-        }
+        labelList myIDs(globalNumbering.toGlobal(setPointIDs));
 
         if (Pstream::master())
         {
diff --git a/applications/utilities/preProcessing/viewFactorsGen/viewFactorsGen.C b/applications/utilities/preProcessing/viewFactorsGen/viewFactorsGen.C
index 3670eceab1177999421b3ca25bd3fba2252cdebb..5bae07c86af572f1a09fca644a88b135080beb16 100644
--- a/applications/utilities/preProcessing/viewFactorsGen/viewFactorsGen.C
+++ b/applications/utilities/preProcessing/viewFactorsGen/viewFactorsGen.C
@@ -822,9 +822,9 @@ int main(int argc, char *argv[])
     {
         const Map<label>& localToCompactMap = compactMap[procI];
 
-        forAllConstIter(Map<label>, localToCompactMap, iter)
+        forAllConstIters(localToCompactMap, iter)
         {
-            compactToGlobal[iter()] = globalNumbering.toGlobal
+            compactToGlobal[*iter] = globalNumbering.toGlobal
             (
                 procI,
                 iter.key()
diff --git a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGProcAgglomerations/GAMGProcAgglomeration/GAMGProcAgglomeration.C b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGProcAgglomerations/GAMGProcAgglomeration/GAMGProcAgglomeration.C
index 96849e559c11dd0e0dbd0053db97e896af79f803..1bb2057f24e7dd087408c5eb1039860e94b78367 100644
--- a/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGProcAgglomerations/GAMGProcAgglomeration/GAMGProcAgglomeration.C
+++ b/src/OpenFOAM/matrices/lduMatrix/solvers/GAMG/GAMGProcAgglomerations/GAMGProcAgglomeration/GAMGProcAgglomeration.C
@@ -133,12 +133,14 @@ Foam::labelListList Foam::GAMGProcAgglomeration::globalCellCells
         Pstream::parRun()
     );
 
-    labelList globalIndices(addr.size());
-    forAll(globalIndices, celli)
-    {
-        globalIndices[celli] = globalNumbering.toGlobal(myProcID, celli);
-    }
-
+    labelList globalIndices
+    (
+        identity
+        (
+            globalNumbering.localSize(myProcID),
+            globalNumbering.localStart(myProcID)
+        )
+    );
 
     // Get the interface cells
     PtrList<labelList> nbrGlobalCells(interfaces.size());
diff --git a/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndex.H b/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndex.H
index f89d245192e9826700a214b5eae19c2359256982..ff387c388106d86cedc3fc97853e5e1a0f077f96 100644
--- a/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndex.H
+++ b/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndex.H
@@ -85,9 +85,9 @@ public:
         inline globalIndex
         (
             const label localSize,
-            const int tag,
-            const label comm,
-            const bool parallel     // use parallel comms
+            const int tag,          //!< message tag
+            const label comm,       //!< communicator
+            const bool parallel     //!< use parallel comms
         );
 
         //- Copy construct from list of labels
@@ -97,7 +97,7 @@ public:
         inline explicit globalIndex(labelList&& offsets);
 
         //- Construct from Istream
-        globalIndex(Istream& is);
+        explicit globalIndex(Istream& is);
 
 
     // Member Functions
@@ -123,9 +123,9 @@ public:
         void reset
         (
             const label localSize,
-            const int tag,
-            const label comm,
-            const bool parallel     // use parallel comms
+            const int tag,          //!< message tag
+            const label comm,       //!< communicator
+            const bool parallel     //!< use parallel comms
         );
 
 
@@ -140,12 +140,18 @@ public:
             //- Return start/size range of local processor data
             inline labelRange range() const;
 
-            //- From local to global
-            inline label toGlobal(const label i) const;
-
             //- Is on local processor
             inline bool isLocal(const label i) const;
 
+            //- From local to global index
+            inline label toGlobal(const label i) const;
+
+            //- From local to global index
+            inline labelList toGlobal(const labelUList& labels) const;
+
+            //- From local to global index (inplace)
+            inline void inplaceToGlobal(labelList& labels) const;
+
             //- From global to local on current processor.
             //  FatalError if not on local processor.
             inline label toLocal(const label i) const;
@@ -168,11 +174,26 @@ public:
             //- Return start/size range of proci data
             inline labelRange range(const label proci) const;
 
+            //- Is on processor proci
+            inline bool isLocal(const label proci, const label i) const;
+
             //- From local to global on proci
             inline label toGlobal(const label proci, const label i) const;
 
-            //- Is on processor proci
-            inline bool isLocal(const label proci, const label i) const;
+            //- From local to global on proci
+            inline labelList toGlobal
+            (
+                const label proci,
+                const labelUList& labels
+            ) const;
+
+            //- From local to global index on proci (inplace)
+            inline void inplaceToGlobal
+            (
+                const label proci,
+                labelList& labels
+            ) const;
+
 
             //- From global to local on proci
             inline label toLocal(const label proci, const label i) const;
diff --git a/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndexI.H b/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndexI.H
index 3adfb13ab09f4d2bedf968e85347f25fa6761392..b2ae4c403037c0ae134ac0329c2238b7128d14a5 100644
--- a/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndexI.H
+++ b/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalIndexI.H
@@ -3,7 +3,7 @@
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
     \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
-     \\/     M anipulation  |
+     \\/     M anipulation  | Copyright (C) 2018 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -130,6 +130,18 @@ inline Foam::label Foam::globalIndex::size() const
 }
 
 
+inline bool Foam::globalIndex::isLocal(const label proci, const label i) const
+{
+    return i >= offsets_[proci] && i < offsets_[proci+1];
+}
+
+
+inline bool Foam::globalIndex::isLocal(const label i) const
+{
+    return isLocal(Pstream::myProcNo(), i);
+}
+
+
 inline Foam::label Foam::globalIndex::toGlobal
 (
     const label proci,
@@ -146,15 +158,46 @@ inline Foam::label Foam::globalIndex::toGlobal(const label i) const
 }
 
 
-inline bool Foam::globalIndex::isLocal(const label proci, const label i) const
+inline Foam::labelList Foam::globalIndex::toGlobal
+(
+    const label proci,
+    const labelUList& labels
+) const
 {
-    return i >= offsets_[proci] && i < offsets_[proci+1];
+    labelList result(labels);
+    inplaceToGlobal(proci, result);
+
+    return result;
 }
 
 
-inline bool Foam::globalIndex::isLocal(const label i) const
+inline Foam::labelList Foam::globalIndex::toGlobal
+(
+    const labelUList& labels
+) const
 {
-    return isLocal(Pstream::myProcNo(), i);
+    return toGlobal(Pstream::myProcNo(), labels);
+}
+
+
+inline void Foam::globalIndex::inplaceToGlobal
+(
+    const label proci,
+    labelList& labels
+) const
+{
+    const label off = offsets_[proci];
+
+    for (label& val : labels)
+    {
+        val += off;
+    }
+}
+
+
+inline void Foam::globalIndex::inplaceToGlobal(labelList& labels) const
+{
+    inplaceToGlobal(Pstream::myProcNo(), labels);
 }
 
 
diff --git a/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalMeshData.C b/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalMeshData.C
index feb7781a19622a64e9a66a796c2e6c9fa0b3c299..9810cb6cfb46ca5c5a8af0446fcfbe7a4b0165c9 100644
--- a/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalMeshData.C
+++ b/src/OpenFOAM/meshes/polyMesh/globalMeshData/globalMeshData.C
@@ -689,12 +689,7 @@ void Foam::globalMeshData::calcGlobalPointEdges
     forAll(pointEdges, pointi)
     {
         const labelList& pEdges = pointEdges[pointi];
-        labelList& globalPEdges = globalPointEdges[pointi];
-        globalPEdges.setSize(pEdges.size());
-        forAll(pEdges, i)
-        {
-            globalPEdges[i] = globalEdgeNumbers.toGlobal(pEdges[i]);
-        }
+        globalPointEdges[pointi] = globalEdgeNumbers.toGlobal(pEdges);
 
         labelPairList& globalPPoints = globalPointPoints[pointi];
         globalPPoints.setSize(pEdges.size());
diff --git a/src/mesh/snappyHexMesh/meshRefinement/meshRefinement.C b/src/mesh/snappyHexMesh/meshRefinement/meshRefinement.C
index 2f7fb70d989b203a34750e4206d5e5702b2f2d81..7b49cfa7b55c30ea3b97dd3c13820cb683a0169b 100644
--- a/src/mesh/snappyHexMesh/meshRefinement/meshRefinement.C
+++ b/src/mesh/snappyHexMesh/meshRefinement/meshRefinement.C
@@ -269,14 +269,7 @@ void Foam::meshRefinement::updateIntersections(const labelList& changedFaces)
     bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
 
     {
-        label nMasterFaces = 0;
-        forAll(isMasterFace, facei)
-        {
-            if (isMasterFace.test(facei))
-            {
-                ++nMasterFaces;
-            }
-        }
+        label nMasterFaces = isMasterFace.count();
         reduce(nMasterFaces, sumOp<label>());
 
         label nChangedFaces = 0;
@@ -2659,11 +2652,10 @@ Foam::bitSet Foam::meshRefinement::getMasterPoints
 {
     const globalIndex globalPoints(meshPoints.size());
 
-    labelList myPoints(meshPoints.size());
-    forAll(meshPoints, pointi)
-    {
-        myPoints[pointi] = globalPoints.toGlobal(pointi);
-    }
+    labelList myPoints
+    (
+        identity(globalPoints.localSize(), globalPoints.localStart())
+    );
 
     syncTools::syncPointList
     (
@@ -2696,11 +2688,10 @@ Foam::bitSet Foam::meshRefinement::getMasterEdges
 {
     const globalIndex globalEdges(meshEdges.size());
 
-    labelList myEdges(meshEdges.size());
-    forAll(meshEdges, edgei)
-    {
-        myEdges[edgei] = globalEdges.toGlobal(edgei);
-    }
+    labelList myEdges
+    (
+        identity(globalEdges.localSize(), globalEdges.localStart())
+    );
 
     syncTools::syncEdgeList
     (
@@ -2741,24 +2732,10 @@ const
 
     {
         bitSet isMasterFace(syncTools::getMasterFaces(mesh_));
-        label nMasterFaces = 0;
-        forAll(isMasterFace, i)
-        {
-            if (isMasterFace.test(i))
-            {
-                ++nMasterFaces;
-            }
-        }
+        label nMasterFaces = isMasterFace.count();
 
         bitSet isMeshMasterPoint(syncTools::getMasterPoints(mesh_));
-        label nMasterPoints = 0;
-        forAll(isMeshMasterPoint, i)
-        {
-            if (isMeshMasterPoint.test(i))
-            {
-                ++nMasterPoints;
-            }
-        }
+        label nMasterPoints = isMeshMasterPoint.count();
 
         Info<< msg.c_str()
             << " : cells:" << pData.nTotalCells()
diff --git a/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C b/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C
index 2a76bc90dce768bf4032b8b56961393a25d4c5bb..871de94683c7c357a5126ec342a6d64696015a8a 100644
--- a/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C
+++ b/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolation.C
@@ -966,10 +966,7 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::update
 
         for (labelList& addressing : tgtAddress_)
         {
-            for (label& addr : addressing)
-            {
-                addr = globalSrcFaces.toGlobal(addr);
-            }
+            globalSrcFaces.inplaceToGlobal(addressing);
         }
 
         // Send data back to originating procs. Note that contributions
diff --git a/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolationParallelOps.C b/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolationParallelOps.C
index 64f803c43fceafce4aadc4baf26b52a13aaad498..dc1b1559f7d8d857235095d482f542475ba0ab80 100644
--- a/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolationParallelOps.C
+++ b/src/meshTools/AMIInterpolation/AMIInterpolation/AMIInterpolationParallelOps.C
@@ -134,12 +134,6 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::distributePatches
 
         if (domain != Pstream::myProcNo() && sendElems.size())
         {
-            labelList globalElems(sendElems.size());
-            forAll(sendElems, i)
-            {
-                globalElems[i] = gi.toGlobal(sendElems[i]);
-            }
-
             faceList subFaces(UIndirectList<face>(pp, sendElems));
             primitivePatch subPatch
             (
@@ -156,7 +150,7 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::distributePatches
             UOPstream toDomain(domain, pBufs);
             toDomain
                 << subPatch.localFaces() << subPatch.localPoints()
-                << globalElems;
+                << gi.toGlobal(sendElems);
         }
     }
 
@@ -186,12 +180,7 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::distributePatches
 
         faces[Pstream::myProcNo()] = subPatch.localFaces();
         points[Pstream::myProcNo()] = subPatch.localPoints();
-
-        faceIDs[Pstream::myProcNo()].setSize(sendElems.size());
-        forAll(sendElems, i)
-        {
-            faceIDs[Pstream::myProcNo()][i] = gi.toGlobal(sendElems[i]);
-        }
+        faceIDs[Pstream::myProcNo()] = gi.toGlobal(sendElems);
     }
 
     // Consume
diff --git a/src/meshTools/regionSplit2D/regionSplit2D.C b/src/meshTools/regionSplit2D/regionSplit2D.C
index aca03a82f52bbef98c1a1751f40f3b81ca913ef4..319a0ae7d63f5161b729a99a04502a1b4f04bb5a 100644
--- a/src/meshTools/regionSplit2D/regionSplit2D.C
+++ b/src/meshTools/regionSplit2D/regionSplit2D.C
@@ -111,10 +111,9 @@ Foam::regionSplit2D::regionSplit2D
 
     // In-place renumber the local regionI to global (compact) regionI
     globalIndex giCompact(compactRegionI);
-    forAllIter(Map<label>, regionToCompactAddr, iter)
+    forAllIters(regionToCompactAddr, iter)
     {
-        label compactRegionI = iter();
-        iter() = giCompact.toGlobal(compactRegionI);
+        *iter = giCompact.toGlobal(*iter);
     }
 
     // Ensure regionToCompactAddr consistent across all processors
diff --git a/src/overset/cellCellStencil/cellCellStencil/cellCellStencil.C b/src/overset/cellCellStencil/cellCellStencil/cellCellStencil.C
index afb63ccaac2e9f91fe183a9f4c4172fc02308fd5..0a8954272a778ace3adfa5cb138fb9a32f72750b 100644
--- a/src/overset/cellCellStencil/cellCellStencil/cellCellStencil.C
+++ b/src/overset/cellCellStencil/cellCellStencil/cellCellStencil.C
@@ -200,11 +200,7 @@ void Foam::cellCellStencil::globalCellCells
 
     // 1. Determine global cell number on other side of coupled patches
 
-    labelList globalCellIDs(mesh.nCells());
-    forAll(globalCellIDs, celli)
-    {
-        globalCellIDs[celli] = gi.toGlobal(celli);
-    }
+    labelList globalCellIDs(identity(gi.localSize(), gi.localStart()));
 
     labelList nbrGlobalCellIDs;
     syncTools::swapBoundaryCellList
diff --git a/src/overset/oversetPolyPatch/oversetGAMGInterface.C b/src/overset/oversetPolyPatch/oversetGAMGInterface.C
index cae275f9574ea882b3edd1992a49b4e9dcd7a1ac..a70c7acf0c9b0277ff1518613949fafb235ff70e 100644
--- a/src/overset/oversetPolyPatch/oversetGAMGInterface.C
+++ b/src/overset/oversetPolyPatch/oversetGAMGInterface.C
@@ -134,12 +134,8 @@ Foam::oversetGAMGInterface::oversetGAMGInterface
         label nCoarseCells = max(restrictMap)+1;
         globalIndex globalNumbering(nCoarseCells);
 
-        labelList globalCoarseIDs(restrictMap.size());
-        forAll(restrictMap, fineCelli)
-        {
-            globalCoarseIDs[fineCelli] =
-                globalNumbering.toGlobal(restrictMap[fineCelli]);
-        }
+        labelList globalCoarseIDs(globalNumbering.toGlobal(restrictMap));
+
         fineMap.distribute(globalCoarseIDs);
 
         //Pout<< this->name()
diff --git a/src/sampling/meshToMesh/meshToMesh.C b/src/sampling/meshToMesh/meshToMesh.C
index b141f3ada8b854a77c130660f8346bdec2ea3db2..7f429a96d99c49c1c615cd98148a547498f80951 100644
--- a/src/sampling/meshToMesh/meshToMesh.C
+++ b/src/sampling/meshToMesh/meshToMesh.C
@@ -524,24 +524,19 @@ void Foam::meshToMesh::calculate(const word& methodName, const bool normalise)
 
         calcAddressing(methodName, srcRegion_, newTgt);
 
-        // per source cell the target cell address in newTgt mesh
-        forAll(srcToTgtCellAddr_, i)
+        // Per source cell the target cell address in newTgt mesh
+        for (labelList& addressing : srcToTgtCellAddr_)
         {
-            labelList& addressing = srcToTgtCellAddr_[i];
-            forAll(addressing, addrI)
+            for (label& addr : addressing)
             {
-                addressing[addrI] = newTgtCellIDs[addressing[addrI]];
+                addr = newTgtCellIDs[addr];
             }
         }
 
-        // convert target addresses in newTgtMesh into global cell numbering
-        forAll(tgtToSrcCellAddr_, i)
+        // Convert target addresses in newTgtMesh into global cell numbering
+        for (labelList& addressing : tgtToSrcCellAddr_)
         {
-            labelList& addressing = tgtToSrcCellAddr_[i];
-            forAll(addressing, addrI)
-            {
-                addressing[addrI] = globalSrcCells.toGlobal(addressing[addrI]);
-            }
+            globalSrcCells.inplaceToGlobal(addressing);
         }
 
         // set up as a reverse distribute
diff --git a/src/sampling/meshToMesh/meshToMeshParallelOps.C b/src/sampling/meshToMesh/meshToMeshParallelOps.C
index b391a6132f1222198cfb46c8675a1b653727ca5b..d60b386ab6a32708c216d809f71e622df0f93812 100644
--- a/src/sampling/meshToMesh/meshToMeshParallelOps.C
+++ b/src/sampling/meshToMesh/meshToMeshParallelOps.C
@@ -455,18 +455,17 @@ void Foam::meshToMesh::distributeCells
             }
 
             // tgt cells into global numbering
-            labelList globalElems(sendElems.size());
-            forAll(sendElems, i)
+            labelList globalElems(globalI.toGlobal(sendElems));
+
+            if (debug > 1)
             {
-                if (debug > 1)
+                forAll(sendElems, i)
                 {
                     Pout<< "tgtProc:" << Pstream::myProcNo()
                         << " sending tgt cell " << sendElems[i]
-                        << "[" << globalI.toGlobal(sendElems[i]) << "]"
+                        << "[" << globalElems[i] << "]"
                         << " to srcProc " << domain << endl;
                 }
-
-                globalElems[i] = globalI.toGlobal(sendElems[i]);
             }
 
             // pass data
diff --git a/src/thermophysicalModels/radiation/radiationModels/viewFactor/viewFactor.C b/src/thermophysicalModels/radiation/radiationModels/viewFactor/viewFactor.C
index d4f7607958bd4cc4b770fd5b2ecb8b00712276ef..959a8f559308c49c9cee379ff551e53c8fe6d57a 100644
--- a/src/thermophysicalModels/radiation/radiationModels/viewFactor/viewFactor.C
+++ b/src/thermophysicalModels/radiation/radiationModels/viewFactor/viewFactor.C
@@ -476,20 +476,13 @@ void Foam::radiation::viewFactor::calculate()
     map_->distribute(compactCoarseHo);
 
     // Distribute local global ID
-    labelList compactGlobalIds(map_->constructSize(), 0.0);
-
-    labelList localGlobalIds(nLocalCoarseFaces_);
-
-    for(label k = 0; k < nLocalCoarseFaces_; k++)
-    {
-        localGlobalIds[k] = globalNumbering.toGlobal(Pstream::myProcNo(), k);
-    }
+    labelList compactGlobalIds(map_->constructSize(), Zero);
 
     SubList<label>
     (
         compactGlobalIds,
         nLocalCoarseFaces_
-    ) = localGlobalIds;
+    ) = identity(globalNumbering.localSize(), globalNumbering.localStart());
 
     map_->distribute(compactGlobalIds);