From 5ac31aabc8bddb8dacc01a4d85db65aa62a2bdbd Mon Sep 17 00:00:00 2001
From: mattijs <mattijs>
Date: Thu, 12 Sep 2013 15:37:53 +0100
Subject: [PATCH] ENH: FieldMapper: extend with unmapped checking flag

---
 .../decomposePar/dimFieldDecomposer.H         |   3 +-
 .../decomposePar/pointFieldDecomposer.C       |   7 +-
 .../decomposePar/pointFieldDecomposer.H       |  10 +-
 .../fields/Fields/Field/FieldMapper.H         |   6 +-
 .../fields/Fields/Field/directFieldMapper.H   | 104 ++++++++++++
 .../directPointPatchFieldMapper.H             | 111 ++++++++++++
 .../pointMesh/pointMeshMapper/pointMapper.H   |   9 +-
 .../pointMeshMapper/pointPatchMapper.C        |  16 +-
 .../pointMeshMapper/pointPatchMapper.H        |   9 +-
 .../mapPolyMesh/cellMapper/cellMapper.H       |  11 +-
 .../mapPolyMesh/faceMapper/faceMapper.H       |   7 +-
 src/dynamicMesh/fvMeshAdder/fvMeshAdder.H     |  49 +-----
 .../fvMeshAdder/fvMeshAdderTemplates.C        |   4 +-
 .../fixedFluxPressureFvPatchScalarField.C     |   7 +-
 .../fvPatchField/directFvPatchFieldMapper.H   | 111 ++++++++++++
 .../fvMesh/fvMeshMapper/fvPatchMapper.C       |   9 +
 .../fvMesh/fvMeshMapper/fvPatchMapper.H       |   9 +-
 .../fvMesh/fvMeshMapper/fvSurfaceMapper.H     |   8 +-
 .../fvMesh/fvMeshSubset/fvMeshSubset.H        |  86 ----------
 .../fvMeshSubset/fvMeshSubsetInterpolate.C    |  91 ++++++----
 .../singleCellFvMesh/singleCellFvMesh.H       |  55 +++---
 .../singleCellFvMeshInterpolate.C             |  81 +++++----
 .../decompose/decompose/fvFieldDecomposer.H   |  20 ++-
 .../fvFieldDecomposerDecomposeFields.C        | 159 ++++++++++++------
 .../reconstruct/fvFieldReconstructor.H        |   7 +-
 .../reconstruct/pointFieldReconstructor.H     |   7 +-
 .../meshToMesh/meshToMesh.H                   |  42 -----
 .../meshToMesh/meshToMeshInterpolate.C        |  38 ++++-
 28 files changed, 714 insertions(+), 362 deletions(-)
 create mode 100644 src/OpenFOAM/fields/Fields/Field/directFieldMapper.H
 create mode 100644 src/OpenFOAM/fields/pointPatchFields/pointPatchField/directPointPatchFieldMapper.H
 create mode 100644 src/finiteVolume/fields/fvPatchFields/fvPatchField/directFvPatchFieldMapper.H

diff --git a/applications/utilities/parallelProcessing/decomposePar/dimFieldDecomposer.H b/applications/utilities/parallelProcessing/decomposePar/dimFieldDecomposer.H
index ac3700d8cd8..5de57330139 100644
--- a/applications/utilities/parallelProcessing/decomposePar/dimFieldDecomposer.H
+++ b/applications/utilities/parallelProcessing/decomposePar/dimFieldDecomposer.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -37,7 +37,6 @@ SourceFiles
 #define dimFieldDecomposer_H
 
 #include "fvMesh.H"
-#include "fvPatchFieldMapper.H"
 #include "surfaceFields.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
diff --git a/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.C b/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.C
index c978cd4f0a5..3934fdc021d 100644
--- a/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.C
+++ b/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -39,7 +39,8 @@ Foam::pointFieldDecomposer::patchFieldDecomposer::patchFieldDecomposer
         completeMeshPatch,
         procMeshPatch
     ),
-    directAddressing_(procMeshPatch.size(), -1)
+    directAddressing_(procMeshPatch.size(), -1),
+    hasUnmapped_(false)
 {
     // Create the inverse-addressing of the patch point labels.
     labelList pointMap(completeMeshPatch.boundaryMesh().mesh().size(), -1);
@@ -64,6 +65,8 @@ Foam::pointFieldDecomposer::patchFieldDecomposer::patchFieldDecomposer
     // Check that all the patch point addresses are set
     if (directAddressing_.size() && min(directAddressing_) < 0)
     {
+        hasUnmapped_ = true;
+
         FatalErrorIn
         (
             "pointFieldDecomposer::patchFieldDecomposer()"
diff --git a/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.H b/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.H
index efe328833bc..7214b9419af 100644
--- a/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.H
+++ b/applications/utilities/parallelProcessing/decomposePar/pointFieldDecomposer.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -63,6 +63,9 @@ public:
 
                 labelList directAddressing_;
 
+                //- Does map contain any unmapped values
+                bool hasUnmapped_;
+
         public:
 
             // Constructors
@@ -88,6 +91,11 @@ public:
                     return true;
                 }
 
+                bool hasUnmapped() const
+                {
+                    return hasUnmapped_;
+                }
+
                 const labelUList& directAddressing() const
                 {
                     return directAddressing_;
diff --git a/src/OpenFOAM/fields/Fields/Field/FieldMapper.H b/src/OpenFOAM/fields/Fields/Field/FieldMapper.H
index 319c8f0575e..f463f97c576 100644
--- a/src/OpenFOAM/fields/Fields/Field/FieldMapper.H
+++ b/src/OpenFOAM/fields/Fields/Field/FieldMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -64,6 +64,10 @@ public:
 
         virtual bool direct() const = 0;
 
+        //- Are there unmapped values? I.e. do all size() elements get
+        //  get value
+        virtual bool hasUnmapped() const = 0;
+
         virtual const labelUList& directAddressing() const
         {
             FatalErrorIn("FieldMapper::directAddressing() const")
diff --git a/src/OpenFOAM/fields/Fields/Field/directFieldMapper.H b/src/OpenFOAM/fields/Fields/Field/directFieldMapper.H
new file mode 100644
index 00000000000..34c8935c445
--- /dev/null
+++ b/src/OpenFOAM/fields/Fields/Field/directFieldMapper.H
@@ -0,0 +1,104 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2013 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::directFieldMapper
+
+Description
+    FieldMapper with direct mapping.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef directFieldMapper_H
+#define directFieldMapper_H
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                           Class directFieldMapper Declaration
+\*---------------------------------------------------------------------------*/
+
+class directFieldMapper
+:
+    public FieldMapper
+{
+    const labelUList& directAddressing_;
+
+    bool hasUnmapped_;
+
+public:
+
+    // Constructors
+
+        //- Construct given addressing
+        patchFieldSubset(const labelUList& directAddressing)
+        :
+            directAddressing_(directAddressing),
+            hasUnmapped_(false)
+        {
+            if (directAddressing_.size() && min(directAddressing_) < 0)
+            {
+                hasUnmapped_ = true;
+            }
+        }
+
+    //- Destructor
+    virtual ~directFieldMapper()
+    {}
+
+
+    // Member Functions
+
+        label size() const
+        {
+            return directAddressing_.size();
+        }
+
+        bool direct() const
+        {
+            return true;
+        }
+
+        bool hasUnmapped() const
+        {
+            return hasUnmapped_;
+        }
+
+        const labelUList& directAddressing() const
+        {
+            return directAddressing_;
+        }
+};
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/OpenFOAM/fields/pointPatchFields/pointPatchField/directPointPatchFieldMapper.H b/src/OpenFOAM/fields/pointPatchFields/pointPatchField/directPointPatchFieldMapper.H
new file mode 100644
index 00000000000..4990061219a
--- /dev/null
+++ b/src/OpenFOAM/fields/pointPatchFields/pointPatchField/directPointPatchFieldMapper.H
@@ -0,0 +1,111 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2013 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::directPointPatchFieldMapper
+
+Description
+    direct pointPatchFieldMapper
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef directPointPatchFieldMapper_H
+#define directPointPatchFieldMapper_H
+
+#include "pointPatchFieldMapper.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                       Class directPointPatchFieldMapper Declaration
+\*---------------------------------------------------------------------------*/
+
+class directPointPatchFieldMapper
+:
+    public pointPatchFieldMapper
+{
+
+    //- Addressing from new back to old
+    const labelUList& directAddressing_;
+
+    //- Does map contain any unmapped values
+    bool hasUnmapped_;
+
+
+public:
+
+    // Constructors
+
+        //- Construct given addressing
+        directPointPatchFieldMapper(const labelUList& directAddressing)
+        :
+            directAddressing_(directAddressing),
+            hasUnmapped_(false)
+        {
+            if (directAddressing_.size() && min(directAddressing_) < 0)
+            {
+                hasUnmapped_ = true;
+            }
+        }
+
+    //- Destructor
+    virtual ~directPointPatchFieldMapper()
+    {}
+
+
+    // Member Functions
+
+        label size() const
+        {
+            return directAddressing_.size();
+        }
+
+        bool direct() const
+        {
+            return true;
+        }
+
+        bool hasUnmapped() const
+        {
+            return hasUnmapped_;
+        }
+
+        const labelUList& directAddressing() const
+        {
+            return directAddressing_;
+        }
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointMapper.H b/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointMapper.H
index 4900c7e0749..4435edf3d73 100644
--- a/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointMapper.H
+++ b/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -129,6 +129,13 @@ public:
             return direct_;
         }
 
+        //- Are there unmapped values? I.e. do all size() elements get
+        //  get value
+        virtual bool hasUnmapped() const
+        {
+            return insertedObjects();
+        }
+
         //- Return direct addressing
         virtual const labelUList& directAddressing() const;
 
diff --git a/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.C b/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.C
index 542455775fd..624897bc098 100644
--- a/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.C
+++ b/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -47,6 +47,8 @@ void Foam::pointPatchMapper::calcAddressing() const
             << abort(FatalError);
     }
 
+    hasUnmapped_ = false;
+
     if (direct())
     {
         // Direct mapping.
@@ -57,7 +59,7 @@ void Foam::pointPatchMapper::calcAddressing() const
         {
             if (addr[i] < 0)
             {
-                addr[i] = 0;
+                hasUnmapped_ = true;
             }
         }
     }
@@ -87,9 +89,11 @@ void Foam::pointPatchMapper::calcAddressing() const
             }
             else
             {
-                // Inserted point. Map from point0 (arbitrary choice)
-                addr[i] = labelList(1, label(0));
-                w[i] = scalarList(1, 1.0);
+                // Inserted point.
+                ///// Map from point0 (arbitrary choice)
+                //addr[i] = labelList(1, label(0));
+                //w[i] = scalarList(1, 1.0);
+                hasUnmapped_ = true;
             }
         }
     }
@@ -101,6 +105,7 @@ void Foam::pointPatchMapper::clearOut()
     deleteDemandDrivenData(directAddrPtr_);
     deleteDemandDrivenData(interpolationAddrPtr_);
     deleteDemandDrivenData(weightsPtr_);
+    hasUnmapped_ = false;
 }
 
 
@@ -124,6 +129,7 @@ Foam::pointPatchMapper::pointPatchMapper
       ? mpm_.oldPatchNMeshPoints()[patch_.index()]
       : 0
     ),
+    hasUnmapped_(false),
     directAddrPtr_(NULL),
     interpolationAddrPtr_(NULL),
     weightsPtr_(NULL)
diff --git a/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.H b/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.H
index 575d6f6bc94..54d89cbdefe 100644
--- a/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.H
+++ b/src/OpenFOAM/meshes/pointMesh/pointMeshMapper/pointPatchMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -74,6 +74,8 @@ class pointPatchMapper
 
     // Demand-driven private data
 
+        mutable bool hasUnmapped_;
+
         //- Direct addressing (only one for of addressing is used)
         mutable labelList* directAddrPtr_;
 
@@ -130,6 +132,11 @@ public:
             return patch_.size();
         }
 
+        virtual bool hasUnmapped() const
+        {
+            return hasUnmapped_;
+        }
+
         //- Return size of field before mapping
         virtual label sizeBeforeMapping() const
         {
diff --git a/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/cellMapper/cellMapper.H b/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/cellMapper/cellMapper.H
index 67c82e43537..03da8e0821b 100644
--- a/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/cellMapper/cellMapper.H
+++ b/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/cellMapper/cellMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -130,6 +130,11 @@ public:
             return direct_;
         }
 
+        virtual bool hasUnmapped() const
+        {
+            return insertedObjects();
+        }
+
         //- Return direct addressing
         virtual const labelUList& directAddressing() const;
 
@@ -140,13 +145,13 @@ public:
         virtual const scalarListList& weights() const;
 
         //- Are there any inserted cells
-        bool insertedObjects() const
+        virtual bool insertedObjects() const
         {
             return insertedCells_;
         }
 
         //- Return list of inserted cells
-        const labelList& insertedObjectLabels() const;
+        const virtual labelList& insertedObjectLabels() const;
 };
 
 
diff --git a/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/faceMapper/faceMapper.H b/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/faceMapper/faceMapper.H
index ad48462250e..199481ea710 100644
--- a/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/faceMapper/faceMapper.H
+++ b/src/OpenFOAM/meshes/polyMesh/mapPolyMesh/faceMapper/faceMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -134,6 +134,11 @@ public:
             return direct_;
         }
 
+        virtual bool hasUnmapped() const
+        {
+            return insertedObjects();
+        }
+
         //- Return direct addressing
         virtual const labelUList& directAddressing() const;
 
diff --git a/src/dynamicMesh/fvMeshAdder/fvMeshAdder.H b/src/dynamicMesh/fvMeshAdder/fvMeshAdder.H
index f246a1887f1..1247e68623f 100644
--- a/src/dynamicMesh/fvMeshAdder/fvMeshAdder.H
+++ b/src/dynamicMesh/fvMeshAdder/fvMeshAdder.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -67,53 +67,6 @@ class fvMeshAdder
 
 private:
 
-
-    // Private class
-
-    class directFvPatchFieldMapper
-    :
-        public fvPatchFieldMapper
-    {
-        // Private data
-
-            const labelList& directAddressing_;
-
-    public:
-
-        // Constructors
-
-            //- Construct from components
-            directFvPatchFieldMapper(const labelList& directAddressing)
-            :
-                fvPatchFieldMapper(),
-                directAddressing_(directAddressing)
-            {}
-
-
-        //- Destructor
-        virtual ~directFvPatchFieldMapper()
-        {}
-
-
-        // Member Functions
-
-            label size() const
-            {
-                return directAddressing_.size();
-            }
-
-            bool direct() const
-            {
-                return true;
-            }
-
-            const labelUList& directAddressing() const
-            {
-                return directAddressing_;
-            }
-    };
-
-
     // Private Member Functions
 
         //- Calculate map from new patch faces to old patch faces. -1 where
diff --git a/src/dynamicMesh/fvMeshAdder/fvMeshAdderTemplates.C b/src/dynamicMesh/fvMeshAdder/fvMeshAdderTemplates.C
index bb01e2fa833..d5afadb9b39 100644
--- a/src/dynamicMesh/fvMeshAdder/fvMeshAdderTemplates.C
+++ b/src/dynamicMesh/fvMeshAdder/fvMeshAdderTemplates.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -23,10 +23,10 @@ License
 
 \*---------------------------------------------------------------------------*/
 
-
 #include "volFields.H"
 #include "surfaceFields.H"
 #include "emptyFvPatchField.H"
+#include "directFvPatchFieldMapper.H"
 
 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
 
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/fixedFluxPressure/fixedFluxPressureFvPatchScalarField.C b/src/finiteVolume/fields/fvPatchFields/derived/fixedFluxPressure/fixedFluxPressureFvPatchScalarField.C
index 3894a5c3e45..7120a099199 100644
--- a/src/finiteVolume/fields/fvPatchFields/derived/fixedFluxPressure/fixedFluxPressureFvPatchScalarField.C
+++ b/src/finiteVolume/fields/fvPatchFields/derived/fixedFluxPressure/fixedFluxPressureFvPatchScalarField.C
@@ -53,8 +53,13 @@ Foam::fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
     fixedGradientFvPatchScalarField(p, iF),
     curTimeIndex_(-1)
 {
+    patchType() = ptf.patchType();
+
     // Map value. Set unmapped values and overwrite with mapped ptf
-    fvPatchField<scalar>::operator=(patchInternalField());
+    if (&iF && iF.size())
+    {
+        fvPatchField<scalar>::operator=(patchInternalField());
+    }
     map(ptf, mapper);
     // Map gradient. Set unmapped values and overwrite with mapped ptf
     gradient() = 0.0;
diff --git a/src/finiteVolume/fields/fvPatchFields/fvPatchField/directFvPatchFieldMapper.H b/src/finiteVolume/fields/fvPatchFields/fvPatchField/directFvPatchFieldMapper.H
new file mode 100644
index 00000000000..7a87d1ad534
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/fvPatchField/directFvPatchFieldMapper.H
@@ -0,0 +1,111 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2013 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::directFvPatchFieldMapper
+
+Description
+    direct fvPatchFieldMapper
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef directFvPatchFieldMapper_H
+#define directFvPatchFieldMapper_H
+
+#include "fvPatchFieldMapper.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                       Class directFvPatchFieldMapper Declaration
+\*---------------------------------------------------------------------------*/
+
+class directFvPatchFieldMapper
+:
+    public fvPatchFieldMapper
+{
+
+    //- Addressing from new back to old
+    const labelUList& directAddressing_;
+
+    //- Does map contain any unmapped values
+    bool hasUnmapped_;
+
+
+public:
+
+    // Constructors
+
+        //- Construct given addressing
+        directFvPatchFieldMapper(const labelUList& directAddressing)
+        :
+            directAddressing_(directAddressing),
+            hasUnmapped_(false)
+        {
+            if (directAddressing_.size() && min(directAddressing_) < 0)
+            {
+                hasUnmapped_ = true;
+            }
+        }
+
+    //- Destructor
+    virtual ~directFvPatchFieldMapper()
+    {}
+
+
+    // Member Functions
+
+        label size() const
+        {
+            return directAddressing_.size();
+        }
+
+        bool direct() const
+        {
+            return true;
+        }
+
+        bool hasUnmapped() const
+        {
+            return hasUnmapped_;
+        }
+
+        const labelUList& directAddressing() const
+        {
+            return directAddressing_;
+        }
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.C b/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.C
index 6a3e337357c..ada0749df30 100644
--- a/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.C
+++ b/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.C
@@ -55,6 +55,8 @@ void Foam::fvPatchMapper::calcAddressing() const
     const label oldPatchEnd =
         oldPatchStart + faceMap_.oldPatchSizes()[patch_.index()];
 
+    hasUnmapped_ = false;
+
     // Assemble the maps: slice to patch
     if (direct())
     {
@@ -84,6 +86,7 @@ void Foam::fvPatchMapper::calcAddressing() const
             {
                 //addr[faceI] = 0;
                 addr[faceI] = -1;
+                hasUnmapped_ = true;
             }
         }
 
@@ -175,6 +178,10 @@ void Foam::fvPatchMapper::calcAddressing() const
                 {
                     newWeights /= sum(newWeights);
                 }
+                else
+                {
+                    hasUnmapped_ = true;
+                }
 
                 // Reset addressing and weights
                 curAddr = newAddr;
@@ -206,6 +213,7 @@ void Foam::fvPatchMapper::clearOut()
     deleteDemandDrivenData(directAddrPtr_);
     deleteDemandDrivenData(interpolationAddrPtr_);
     deleteDemandDrivenData(weightsPtr_);
+    hasUnmapped_ = false;
 }
 
 
@@ -221,6 +229,7 @@ Foam::fvPatchMapper::fvPatchMapper
     patch_(patch),
     faceMap_(faceMap),
     sizeBeforeMapping_(faceMap.oldPatchSizes()[patch_.index()]),
+    hasUnmapped_(false),
     directAddrPtr_(NULL),
     interpolationAddrPtr_(NULL),
     weightsPtr_(NULL)
diff --git a/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.H b/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.H
index aa169fa082a..bcfd49566b7 100644
--- a/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.H
+++ b/src/finiteVolume/fvMesh/fvMeshMapper/fvPatchMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -71,6 +71,8 @@ class fvPatchMapper
 
     // Demand-driven private data
 
+        mutable bool hasUnmapped_;
+
         //- Direct addressing (only one for of addressing is used)
         mutable labelList* directAddrPtr_;
 
@@ -138,6 +140,11 @@ public:
             return faceMap_.direct();
         }
 
+        virtual bool hasUnmapped() const
+        {
+            return hasUnmapped_;
+        }
+
         //- Return direct addressing
         virtual const labelUList& directAddressing() const;
 
diff --git a/src/finiteVolume/fvMesh/fvMeshMapper/fvSurfaceMapper.H b/src/finiteVolume/fvMesh/fvMeshMapper/fvSurfaceMapper.H
index c144bb99603..bf47b1f55c4 100644
--- a/src/finiteVolume/fvMesh/fvMeshMapper/fvSurfaceMapper.H
+++ b/src/finiteVolume/fvMesh/fvMeshMapper/fvSurfaceMapper.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -133,6 +133,12 @@ public:
             return faceMap_.direct();
         }
 
+        //- Has unmapped elements
+        virtual bool hasUnmapped() const
+        {
+            return insertedObjects();
+        }
+
         //- Return direct addressing
         virtual const labelUList& directAddressing() const;
 
diff --git a/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubset.H b/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubset.H
index 2406da5db52..5f2a281d065 100644
--- a/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubset.H
+++ b/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubset.H
@@ -56,8 +56,6 @@ SourceFiles
 
 #include "fvMesh.H"
 #include "pointMesh.H"
-#include "fvPatchFieldMapper.H"
-#include "pointPatchFieldMapper.H"
 #include "GeometricField.H"
 #include "HashSet.H"
 #include "surfaceMesh.H"
@@ -74,90 +72,6 @@ namespace Foam
 class fvMeshSubset
 {
 
-public:
-
-    //- Patch-field subset interpolation class
-    class patchFieldSubset
-    :
-        public fvPatchFieldMapper
-    {
-        const labelList& directAddressing_;
-
-    public:
-
-        // Constructors
-
-            //- Construct given addressing
-            patchFieldSubset(const labelList& directAddressing)
-            :
-                directAddressing_(directAddressing)
-            {}
-
-        //- Destructor
-        virtual ~patchFieldSubset()
-        {}
-
-
-        // Member Functions
-
-            label size() const
-            {
-                return directAddressing_.size();
-            }
-
-            bool direct() const
-            {
-                return true;
-            }
-
-            const labelUList& directAddressing() const
-            {
-                return directAddressing_;
-            }
-    };
-
-
-    //- Patch-field subset interpolation class
-    class pointPatchFieldSubset
-    :
-        public pointPatchFieldMapper
-    {
-        const labelList& directAddressing_;
-
-    public:
-
-        // Constructors
-
-            //- Construct given addressing
-            pointPatchFieldSubset(const labelList& directAddressing)
-            :
-                directAddressing_(directAddressing)
-            {}
-
-        //- Destructor
-        virtual ~pointPatchFieldSubset()
-        {}
-
-
-        // Member Functions
-
-            label size() const
-            {
-                return directAddressing_.size();
-            }
-
-            bool direct() const
-            {
-                return true;
-            }
-
-            const labelUList& directAddressing() const
-            {
-                return directAddressing_;
-            }
-    };
-
-
 private:
 
     // Private data
diff --git a/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubsetInterpolate.C b/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubsetInterpolate.C
index 130b08072a8..78811948187 100644
--- a/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubsetInterpolate.C
+++ b/src/finiteVolume/fvMesh/fvMeshSubset/fvMeshSubsetInterpolate.C
@@ -27,6 +27,8 @@ License
 #include "emptyFvsPatchField.H"
 #include "emptyPointPatchField.H"
 #include "emptyFvPatchFields.H"
+#include "directFvPatchFieldMapper.H"
+#include "directPointPatchFieldMapper.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
@@ -70,8 +72,9 @@ tmp<GeometricField<Type, fvPatchField, volMesh> > fvMeshSubset::interpolate
             patchFields.set
             (
                 patchI,
-                new calculatedFvPatchField<Type>
+                fvPatchField<Type>::New
                 (
+                    calculatedFvPatchField<Type>::typeName,
                     sMesh.boundary()[patchI],
                     DimensionedField<Type, volMesh>::null()
                 )
@@ -142,7 +145,7 @@ tmp<GeometricField<Type, fvPatchField, volMesh> > fvMeshSubset::interpolate
                     vf.boundaryField()[patchMap[patchI]],
                     subPatch,
                     resF.dimensionedInternalField(),
-                    patchFieldSubset(directAddressing)
+                    directFvPatchFieldMapper(directAddressing)
                 )
             );
         }
@@ -203,8 +206,9 @@ tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > fvMeshSubset::interpolate
             patchFields.set
             (
                 patchI,
-                new calculatedFvsPatchField<Type>
+                fvsPatchField<Type>::New
                 (
+                    calculatedFvsPatchField<Type>::typeName,
                     sMesh.boundary()[patchI],
                     DimensionedField<Type, surfaceMesh>::null()
                 )
@@ -285,7 +289,7 @@ tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > fvMeshSubset::interpolate
                     vf.boundaryField()[patchMap[patchI]],
                     subPatch,
                     resF.dimensionedInternalField(),
-                    patchFieldSubset(directAddressing)
+                    directFvPatchFieldMapper(directAddressing)
                 )
             );
 
@@ -348,10 +352,7 @@ fvMeshSubset::interpolate
     const labelList& pointMap
 )
 {
-    // Create and map the internal-field values
-    Field<Type> internalField(vf.internalField(), pointMap);
-
-    // Create and map the patch field values
+    // 1. Create the complete field with dummy patch fields
     PtrList<pointPatchField<Type> > patchFields(patchMap.size());
 
     forAll(patchFields, patchI)
@@ -372,6 +373,54 @@ fvMeshSubset::interpolate
             );
         }
         else
+        {
+            patchFields.set
+            (
+                patchI,
+                pointPatchField<Type>::New
+                (
+                    calculatedPointPatchField<Type>::typeName,
+                    sMesh.boundary()[patchI],
+                    DimensionedField<Type, pointMesh>::null()
+                )
+            );
+        }
+    }
+
+    // Create the complete field from the pieces
+    tmp<GeometricField<Type, pointPatchField, pointMesh> > tresF
+    (
+        new GeometricField<Type, pointPatchField, pointMesh>
+        (
+            IOobject
+            (
+                "subset"+vf.name(),
+                sMesh.time().timeName(),
+                sMesh.thisDb(),
+                IOobject::NO_READ,
+                IOobject::NO_WRITE
+            ),
+            sMesh,
+            vf.dimensions(),
+            Field<Type>(vf.internalField(), pointMap),
+            patchFields
+        )
+    );
+    GeometricField<Type, pointPatchField, pointMesh>& resF = tresF();
+
+
+    // 2. Change the pointPatchFields to the correct type using a mapper
+    //  constructor (with reference to the now correct internal field)
+
+    typename GeometricField<Type, pointPatchField, pointMesh>::
+        GeometricBoundaryField& bf = resF.boundaryField();
+
+    forAll(bf, patchI)
+    {
+        // Set the first one by hand as it corresponds to the
+        // exposed internal faces.  Additional interpolation can be put here
+        // as necessary.
+        if (patchMap[patchI] != -1)
         {
             // Construct addressing
             const pointPatch& basePatch =
@@ -406,40 +455,20 @@ fvMeshSubset::interpolate
                 }
             }
 
-            patchFields.set
+            bf.set
             (
                 patchI,
                 pointPatchField<Type>::New
                 (
                     vf.boundaryField()[patchMap[patchI]],
                     subPatch,
-                    DimensionedField<Type, pointMesh>::null(),
-                    pointPatchFieldSubset(directAddressing)
+                    resF.dimensionedInternalField(),
+                    directPointPatchFieldMapper(directAddressing)
                 )
             );
         }
     }
 
-    // Create the complete field from the pieces
-    tmp<GeometricField<Type, pointPatchField, pointMesh> > tresF
-    (
-        new GeometricField<Type, pointPatchField, pointMesh>
-        (
-            IOobject
-            (
-                "subset"+vf.name(),
-                vf.time().timeName(),
-                sMesh.thisDb(),
-                IOobject::NO_READ,
-                IOobject::NO_WRITE
-            ),
-            sMesh,
-            vf.dimensions(),
-            internalField,
-            patchFields
-        )
-    );
-
     return tresF;
 }
 
diff --git a/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMesh.H b/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMesh.H
index 7ee6c37b687..036e4f24c93 100644
--- a/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMesh.H
+++ b/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMesh.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -87,39 +87,6 @@ class singleCellFvMesh
 
 public:
 
-        //- Patch field mapper class for non-agglomerated meshes
-        class directPatchFieldMapper
-        :
-            public fvPatchFieldMapper
-        {
-            // Private data
-
-                const labelUList& directAddressing_;
-
-        public:
-
-                //- Construct given addressing
-                directPatchFieldMapper(const labelUList& directAddressing)
-                :
-                    directAddressing_(directAddressing)
-                {}
-
-                virtual label size() const
-                {
-                    return directAddressing_.size();
-                }
-
-                virtual bool direct() const
-                {
-                    return true;
-                }
-
-                virtual const labelUList& directAddressing() const
-                {
-                    return directAddressing_;
-                }
-        };
-
         //- Patch field mapper class for agglomerated meshes
         class agglomPatchFieldMapper
         :
@@ -129,6 +96,7 @@ public:
 
                 const labelListList& addressing_;
                 const scalarListList& weights_;
+                bool hasUnmapped_;
 
         public:
 
@@ -140,8 +108,18 @@ public:
                 )
                 :
                     addressing_(addressing),
-                    weights_(weights)
-                {}
+                    weights_(weights),
+                    hasUnmapped_(false)
+                {
+                    forAll(addressing_, i)
+                    {
+                        if (addressing_[i].empty())
+                        {
+                            hasUnmapped_ = true;
+                            break;
+                        }
+                    }
+                }
 
                 virtual label size() const
                 {
@@ -153,6 +131,11 @@ public:
                     return false;
                 }
 
+                bool hasUnmapped() const
+                {
+                    return hasUnmapped_;
+                }
+
                 virtual const labelListList& addressing() const
                 {
                     return addressing_;
diff --git a/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMeshInterpolate.C b/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMeshInterpolate.C
index 33e8316ac7f..93ef66101f4 100644
--- a/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMeshInterpolate.C
+++ b/src/finiteVolume/fvMesh/singleCellFvMesh/singleCellFvMeshInterpolate.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -24,6 +24,8 @@ License
 \*---------------------------------------------------------------------------*/
 
 #include "singleCellFvMesh.H"
+#include "calculatedFvPatchFields.H"
+#include "directFvPatchFieldMapper.H"
 #include "Time.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@@ -39,12 +41,51 @@ tmp<GeometricField<Type, fvPatchField, volMesh> > singleCellFvMesh::interpolate
     const GeometricField<Type, fvPatchField, volMesh>& vf
 ) const
 {
-    // Create internal-field values
-    Field<Type> internalField(1, gAverage(vf));
-
-    // Create and map the patch field values
+    // 1. Create the complete field with dummy patch fields
     PtrList<fvPatchField<Type> > patchFields(vf.boundaryField().size());
 
+    forAll(patchFields, patchI)
+    {
+        patchFields.set
+        (
+            patchI,
+            fvPatchField<Type>::New
+            (
+                calculatedFvPatchField<Type>::typeName,
+                boundary()[patchI],
+                DimensionedField<Type, volMesh>::null()
+            )
+        );
+    }
+
+    // Create the complete field from the pieces
+    tmp<GeometricField<Type, fvPatchField, volMesh> > tresF
+    (
+        new GeometricField<Type, fvPatchField, volMesh>
+        (
+            IOobject
+            (
+                vf.name(),
+                time().timeName(),
+                *this,
+                IOobject::NO_READ,
+                IOobject::NO_WRITE
+            ),
+            *this,
+            vf.dimensions(),
+            Field<Type>(1, gAverage(vf)),
+            patchFields
+        )
+    );
+    GeometricField<Type, fvPatchField, volMesh>& resF = tresF();
+
+
+    // 2. Change the fvPatchFields to the correct type using a mapper
+    //  constructor (with reference to the now correct internal field)
+
+    typename GeometricField<Type, fvPatchField, volMesh>::
+        GeometricBoundaryField& bf = resF.boundaryField();
+
     if (agglomerate())
     {
         forAll(vf.boundaryField(), patchI)
@@ -67,14 +108,14 @@ tmp<GeometricField<Type, fvPatchField, volMesh> > singleCellFvMesh::interpolate
                 );
             }
 
-            patchFields.set
+            bf.set
             (
                 patchI,
                 fvPatchField<Type>::New
                 (
                     vf.boundaryField()[patchI],
                     boundary()[patchI],
-                    DimensionedField<Type, volMesh>::null(),
+                    resF.dimensionedInternalField(),
                     agglomPatchFieldMapper(coarseToFine, coarseWeights)
                 )
             );
@@ -86,40 +127,20 @@ tmp<GeometricField<Type, fvPatchField, volMesh> > singleCellFvMesh::interpolate
         {
             labelList map(identity(vf.boundaryField()[patchI].size()));
 
-            patchFields.set
+            bf.set
             (
                 patchI,
                 fvPatchField<Type>::New
                 (
                     vf.boundaryField()[patchI],
                     boundary()[patchI],
-                    DimensionedField<Type, volMesh>::null(),
-                    directPatchFieldMapper(map)
+                    resF.dimensionedInternalField(),
+                    directFvPatchFieldMapper(map)
                 )
             );
         }
     }
 
-    // Create the complete field from the pieces
-    tmp<GeometricField<Type, fvPatchField, volMesh> > tresF
-    (
-        new GeometricField<Type, fvPatchField, volMesh>
-        (
-            IOobject
-            (
-                vf.name(),
-                time().timeName(),
-                *this,
-                IOobject::NO_READ,
-                IOobject::NO_WRITE
-            ),
-            *this,
-            vf.dimensions(),
-            internalField,
-            patchFields
-        )
-    );
-
     return tresF;
 }
 
diff --git a/src/parallel/decompose/decompose/fvFieldDecomposer.H b/src/parallel/decompose/decompose/fvFieldDecomposer.H
index 899b45b4b75..0bf83369e62 100644
--- a/src/parallel/decompose/decompose/fvFieldDecomposer.H
+++ b/src/parallel/decompose/decompose/fvFieldDecomposer.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -88,6 +88,12 @@ public:
                     return true;
                 }
 
+                //- Are there unmapped values
+                bool hasUnmapped() const
+                {
+                    return false;
+                }
+
                 const labelUList& directAddressing() const
                 {
                     return directAddressing_;
@@ -128,6 +134,12 @@ public:
                     return true;
                 }
 
+                //- Are there unmapped values
+                bool hasUnmapped() const
+                {
+                    return false;
+                }
+
                 const labelUList& directAddressing() const
                 {
                     return directAddressing_;
@@ -165,6 +177,12 @@ public:
                     return false;
                 }
 
+                //- Are there unmapped values
+                bool hasUnmapped() const
+                {
+                    return false;
+                }
+
                 const labelListList& addressing() const
                 {
                     return addressing_;
diff --git a/src/parallel/decompose/decompose/fvFieldDecomposerDecomposeFields.C b/src/parallel/decompose/decompose/fvFieldDecomposerDecomposeFields.C
index ca9d0e31b5e..c72ceedfb89 100644
--- a/src/parallel/decompose/decompose/fvFieldDecomposerDecomposeFields.C
+++ b/src/parallel/decompose/decompose/fvFieldDecomposerDecomposeFields.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -40,37 +40,76 @@ Foam::fvFieldDecomposer::decomposeField
     const bool allowUnknownPatchFields
 ) const
 {
-    // Create and map the internal field values
-    Field<Type> internalField(field.internalField(), cellAddressing_);
-
-    // Create and map the patch field values
+    // 1. Create the complete field with dummy patch fields
     PtrList<fvPatchField<Type> > patchFields(boundaryAddressing_.size());
 
     forAll(boundaryAddressing_, patchi)
+    {
+        patchFields.set
+        (
+            patchi,
+            fvPatchField<Type>::New
+            (
+                calculatedFvPatchField<Type>::typeName,
+                procMesh_.boundary()[patchi],
+                DimensionedField<Type, volMesh>::null()
+            )
+        );
+    }
+
+    // Create the field for the processor
+    tmp<GeometricField<Type, fvPatchField, volMesh> > tresF
+    (
+        new GeometricField<Type, fvPatchField, volMesh>
+        (
+            IOobject
+            (
+                field.name(),
+                procMesh_.time().timeName(),
+                procMesh_,
+                IOobject::NO_READ,
+                IOobject::NO_WRITE
+            ),
+            procMesh_,
+            field.dimensions(),
+            Field<Type>(field.internalField(), cellAddressing_),
+            patchFields
+        )
+    );
+    GeometricField<Type, fvPatchField, volMesh>& resF = tresF();
+
+
+    // 2. Change the fvPatchFields to the correct type using a mapper
+    //  constructor (with reference to the now correct internal field)
+
+    typename GeometricField<Type, fvPatchField, volMesh>::
+        GeometricBoundaryField& bf = resF.boundaryField();
+
+    forAll(bf, patchi)
     {
         if (patchFieldDecomposerPtrs_[patchi])
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 fvPatchField<Type>::New
                 (
                     field.boundaryField()[boundaryAddressing_[patchi]],
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, volMesh>::null(),
+                    resF.dimensionedInternalField(),
                     *patchFieldDecomposerPtrs_[patchi]
                 )
             );
         }
         else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 new processorCyclicFvPatchField<Type>
                 (
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, volMesh>::null(),
+                    resF.dimensionedInternalField(),
                     Field<Type>
                     (
                         field.internalField(),
@@ -81,13 +120,13 @@ Foam::fvFieldDecomposer::decomposeField
         }
         else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 new processorFvPatchField<Type>
                 (
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, volMesh>::null(),
+                    resF.dimensionedInternalField(),
                     Field<Type>
                     (
                         field.internalField(),
@@ -98,13 +137,13 @@ Foam::fvFieldDecomposer::decomposeField
         }
         else if (allowUnknownPatchFields)
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 new emptyFvPatchField<Type>
                 (
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, volMesh>::null()
+                    resF.dimensionedInternalField()
                 )
             );
         }
@@ -116,24 +155,7 @@ Foam::fvFieldDecomposer::decomposeField
     }
 
     // Create the field for the processor
-    return tmp<GeometricField<Type, fvPatchField, volMesh> >
-    (
-        new GeometricField<Type, fvPatchField, volMesh>
-        (
-            IOobject
-            (
-                field.name(),
-                procMesh_.time().timeName(),
-                procMesh_,
-                IOobject::NO_READ,
-                IOobject::NO_WRITE
-            ),
-            procMesh_,
-            field.dimensions(),
-            internalField,
-            patchFields
-        )
-    );
+    return tresF;
 }
 
 
@@ -188,34 +210,76 @@ Foam::fvFieldDecomposer::decomposeField
         }
     }
 
-    // Create and map the patch field values
+
+    // 1. Create the complete field with dummy patch fields
     PtrList<fvsPatchField<Type> > patchFields(boundaryAddressing_.size());
 
+    forAll(boundaryAddressing_, patchi)
+    {
+        patchFields.set
+        (
+            patchi,
+            fvsPatchField<Type>::New
+            (
+                calculatedFvsPatchField<Type>::typeName,
+                procMesh_.boundary()[patchi],
+                DimensionedField<Type, surfaceMesh>::null()
+            )
+        );
+    }
+
+    tmp<GeometricField<Type, fvsPatchField, surfaceMesh> > tresF
+    (
+        new GeometricField<Type, fvsPatchField, surfaceMesh>
+        (
+            IOobject
+            (
+                field.name(),
+                procMesh_.time().timeName(),
+                procMesh_,
+                IOobject::NO_READ,
+                IOobject::NO_WRITE
+            ),
+            procMesh_,
+            field.dimensions(),
+            Field<Type>(field.internalField(), mapAddr),
+            patchFields
+        )
+    );
+    GeometricField<Type, fvsPatchField, surfaceMesh>& resF = tresF();
+
+
+    // 2. Change the fvsPatchFields to the correct type using a mapper
+    //  constructor (with reference to the now correct internal field)
+
+    typename GeometricField<Type, fvsPatchField, surfaceMesh>::
+        GeometricBoundaryField& bf = resF.boundaryField();
+
     forAll(boundaryAddressing_, patchi)
     {
         if (patchFieldDecomposerPtrs_[patchi])
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 fvsPatchField<Type>::New
                 (
                     field.boundaryField()[boundaryAddressing_[patchi]],
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, surfaceMesh>::null(),
+                    resF.dimensionedInternalField(),
                     *patchFieldDecomposerPtrs_[patchi]
                 )
             );
         }
         else if (isA<processorCyclicFvPatch>(procMesh_.boundary()[patchi]))
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 new processorCyclicFvsPatchField<Type>
                 (
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, surfaceMesh>::null(),
+                    resF.dimensionedInternalField(),
                     Field<Type>
                     (
                         allFaceField,
@@ -226,13 +290,13 @@ Foam::fvFieldDecomposer::decomposeField
         }
         else if (isA<processorFvPatch>(procMesh_.boundary()[patchi]))
         {
-            patchFields.set
+            bf.set
             (
                 patchi,
                 new processorFvsPatchField<Type>
                 (
                     procMesh_.boundary()[patchi],
-                    DimensionedField<Type, surfaceMesh>::null(),
+                    resF.dimensionedInternalField(),
                     Field<Type>
                     (
                         allFaceField,
@@ -249,24 +313,7 @@ Foam::fvFieldDecomposer::decomposeField
     }
 
     // Create the field for the processor
-    return tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
-    (
-        new GeometricField<Type, fvsPatchField, surfaceMesh>
-        (
-            IOobject
-            (
-                field.name(),
-                procMesh_.time().timeName(),
-                procMesh_,
-                IOobject::NO_READ,
-                IOobject::NO_WRITE
-            ),
-            procMesh_,
-            field.dimensions(),
-            internalField,
-            patchFields
-        )
-    );
+    return tresF;
 }
 
 
diff --git a/src/parallel/reconstruct/reconstruct/fvFieldReconstructor.H b/src/parallel/reconstruct/reconstruct/fvFieldReconstructor.H
index 43469d24b4f..0ee1decbb4c 100644
--- a/src/parallel/reconstruct/reconstruct/fvFieldReconstructor.H
+++ b/src/parallel/reconstruct/reconstruct/fvFieldReconstructor.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -116,6 +116,11 @@ public:
                     return true;
                 }
 
+                bool hasUnmapped() const
+                {
+                    return false;
+                }
+
                 const labelUList& directAddressing() const
                 {
                     return labelUList::null();
diff --git a/src/parallel/reconstruct/reconstruct/pointFieldReconstructor.H b/src/parallel/reconstruct/reconstruct/pointFieldReconstructor.H
index 47ad5b23c2d..c15b7111da8 100644
--- a/src/parallel/reconstruct/reconstruct/pointFieldReconstructor.H
+++ b/src/parallel/reconstruct/reconstruct/pointFieldReconstructor.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -112,6 +112,11 @@ public:
                     return true;
                 }
 
+                bool hasUnmapped() const
+                {
+                    return false;
+                }
+
                 const labelUList& directAddressing() const
                 {
                     return labelUList::null();
diff --git a/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMesh.H b/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMesh.H
index b99c0629ccb..bfedab7433d 100644
--- a/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMesh.H
+++ b/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMesh.H
@@ -175,48 +175,6 @@ public:
     ~meshToMesh();
 
 
-    //- Patch-field interpolation class
-    class patchFieldInterpolator
-    :
-        public fvPatchFieldMapper
-    {
-        const labelList& directAddressing_;
-
-    public:
-
-        // Constructors
-
-            //- Construct given addressing
-            patchFieldInterpolator(const labelList& addr)
-            :
-                directAddressing_(addr)
-            {}
-
-
-        //- Destructor
-        virtual ~patchFieldInterpolator()
-        {}
-
-
-        // Member Functions
-
-            label size() const
-            {
-                return directAddressing_.size();
-            }
-
-            bool direct() const
-            {
-                return true;
-            }
-
-            const labelList& directAddressing() const
-            {
-                return directAddressing_;
-            }
-    };
-
-
     // Member Functions
 
         // Access
diff --git a/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMeshInterpolate.C b/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMeshInterpolate.C
index 372194e3225..42a5269cd2d 100644
--- a/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMeshInterpolate.C
+++ b/src/sampling/meshToMeshInterpolation/meshToMesh/meshToMeshInterpolate.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2011-2013 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -28,6 +28,7 @@ License
 #include "interpolationCellPoint.H"
 #include "SubField.H"
 #include "mixedFvPatchField.H"
+#include "directFvPatchFieldMapper.H"
 
 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
 
@@ -408,7 +409,7 @@ Foam::meshToMesh::interpolate
             << exit(FatalError);
     }
 
-    // Create and map the patch field values
+    // 1. Create the complete field with dummy patch fields
     PtrList<fvPatchField<Type> > patchFields
     (
         boundaryAddressing_.size()
@@ -421,13 +422,9 @@ Foam::meshToMesh::interpolate
             patchI,
             fvPatchField<Type>::New
             (
-                fromVf.boundaryField()[patchI],
+                calculatedFvPatchField<Type>::typeName,
                 toMesh_.boundary()[patchI],
-                DimensionedField<Type, volMesh>::null(),
-                patchFieldInterpolator
-                (
-                    boundaryAddressing_[patchI]
-                )
+                DimensionedField<Type, volMesh>::null()
             )
         );
     }
@@ -452,6 +449,31 @@ Foam::meshToMesh::interpolate
             patchFields
         )
     );
+    GeometricField<Type, fvPatchField, volMesh>& toF = ttoF();
+
+    // 2. Change the fvPatchFields to the correct type using a mapper
+    //  constructor (with reference to the now correct internal field)
+
+    typename GeometricField<Type, fvPatchField, volMesh>::
+        GeometricBoundaryField& bf = toF.boundaryField();
+
+    forAll(boundaryAddressing_, patchI)
+    {
+        bf.set
+        (
+            patchI,
+            fvPatchField<Type>::New
+            (
+                fromVf.boundaryField()[patchI],
+                toMesh_.boundary()[patchI],
+                toF.dimensionedInternalField(),
+                directFvPatchFieldMapper
+                (
+                    boundaryAddressing_[patchI]
+                )
+            )
+        );
+    }
 
     return ttoF;
 }
-- 
GitLab