diff --git a/src/finiteVolume/Make/files b/src/finiteVolume/Make/files
index 521ac3647b7fbc507ebdaff4914d244ef6fdaccb..a99d3e3a3196dd0a5791fdc876c18bde5b85dedd 100644
--- a/src/finiteVolume/Make/files
+++ b/src/finiteVolume/Make/files
@@ -183,6 +183,7 @@ $(derivedFvPatchFields)/turbulentIntensityKineticEnergyInlet/turbulentIntensityK
 $(derivedFvPatchFields)/uniformDensityHydrostaticPressure/uniformDensityHydrostaticPressureFvPatchScalarField.C
 $(derivedFvPatchFields)/uniformFixedGradient/uniformFixedGradientFvPatchFields.C
 $(derivedFvPatchFields)/uniformFixedValue/uniformFixedValueFvPatchFields.C
+$(derivedFvPatchFields)/uniformInletOutlet/uniformInletOutletFvPatchFields.C
 $(derivedFvPatchFields)/uniformJump/uniformJumpFvPatchFields.C
 $(derivedFvPatchFields)/uniformJumpAMI/uniformJumpAMIFvPatchFields.C
 $(derivedFvPatchFields)/uniformTotalPressure/uniformTotalPressureFvPatchScalarField.C
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchField.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchField.C
new file mode 100644
index 0000000000000000000000000000000000000000..7b5f58a2d0cf1171587f0a77a1639c57c9d41457
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchField.C
@@ -0,0 +1,216 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#include "uniformInletOutletFvPatchField.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class Type>
+Foam::uniformInletOutletFvPatchField<Type>::uniformInletOutletFvPatchField
+(
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF
+)
+:
+    mixedFvPatchField<Type>(p, iF),
+    phiName_("phi")
+{
+    this->refValue() = pTraits<Type>::zero;
+    this->refGrad() = pTraits<Type>::zero;
+    this->valueFraction() = 0.0;
+}
+
+
+template<class Type>
+Foam::uniformInletOutletFvPatchField<Type>::uniformInletOutletFvPatchField
+(
+    const uniformInletOutletFvPatchField<Type>& ptf,
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF,
+    const fvPatchFieldMapper& mapper
+)
+:
+    mixedFvPatchField<Type>(p, iF),
+    phiName_(ptf.phiName_),
+    uniformInletValue_(ptf.uniformInletValue_, false)
+{
+    this->patchType() = ptf.patchType();
+
+    // For safety re-evaluate
+    const scalar t = this->db().time().timeOutputValue();
+    this->refValue() = uniformInletValue_->value(t);
+    this->refGrad() = pTraits<Type>::zero;
+    this->valueFraction() = 0.0;
+
+    // Map value (unmapped get refValue)
+    if (&iF && iF.size())
+    {
+        fvPatchField<Type>::operator=(this->refValue());
+    }
+    this->map(ptf, mapper);
+
+}
+
+
+template<class Type>
+Foam::uniformInletOutletFvPatchField<Type>::uniformInletOutletFvPatchField
+(
+    const fvPatch& p,
+    const DimensionedField<Type, volMesh>& iF,
+    const dictionary& dict
+)
+:
+    mixedFvPatchField<Type>(p, iF),
+    phiName_(dict.lookupOrDefault<word>("phi", "phi")),
+    uniformInletValue_(DataEntry<Type>::New("uniformInletValue", dict))
+{
+    const scalar t = this->db().time().timeOutputValue();
+    this->refValue() = uniformInletValue_->value(t);
+
+    if (dict.found("value"))
+    {
+        fvPatchField<Type>::operator=
+        (
+            Field<Type>("value", dict, p.size())
+        );
+    }
+    else
+    {
+        fvPatchField<Type>::operator=(this->refValue());
+    }
+
+    this->refGrad() = pTraits<Type>::zero;
+    this->valueFraction() = 0.0;
+}
+
+
+template<class Type>
+Foam::uniformInletOutletFvPatchField<Type>::uniformInletOutletFvPatchField
+(
+    const uniformInletOutletFvPatchField<Type>& ptf
+)
+:
+    mixedFvPatchField<Type>(ptf),
+    phiName_(ptf.phiName_),
+    uniformInletValue_(ptf.uniformInletValue_, false)
+{}
+
+
+template<class Type>
+Foam::uniformInletOutletFvPatchField<Type>::uniformInletOutletFvPatchField
+(
+    const uniformInletOutletFvPatchField<Type>& ptf,
+    const DimensionedField<Type, volMesh>& iF
+)
+:
+    mixedFvPatchField<Type>(ptf, iF),
+    phiName_(ptf.phiName_),
+    uniformInletValue_(ptf.uniformInletValue_, false)
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::uniformInletOutletFvPatchField<Type>::updateCoeffs()
+{
+    if (this->updated())
+    {
+        return;
+    }
+
+    const Field<scalar>& phip =
+        this->patch().template lookupPatchField<surfaceScalarField, scalar>
+        (
+            phiName_
+        );
+
+    this->valueFraction() = 1.0 - pos(phip);
+
+    mixedFvPatchField<Type>::updateCoeffs();
+}
+
+
+template<class Type>
+void Foam::uniformInletOutletFvPatchField<Type>::write(Ostream& os) const
+{
+    fvPatchField<Type>::write(os);
+    if (phiName_ != "phi")
+    {
+        os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
+    }
+    this->uniformInletValue_->writeData(os);
+    this->writeEntry("value", os);
+}
+
+
+// * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::uniformInletOutletFvPatchField<Type>::autoMap
+(
+    const fvPatchFieldMapper& m
+)
+{
+    mixedFvPatchField<Type>::autoMap(m);
+
+    // Override
+    const scalar t = this->db().time().timeOutputValue();
+    this->refValue() = uniformInletValue_->value(t);
+}
+
+
+template<class Type>
+void Foam::uniformInletOutletFvPatchField<Type>::rmap
+(
+    const fvPatchField<Type>& ptf,
+    const labelList& addr
+)
+{
+    mixedFvPatchField<Type>::rmap(ptf, addr);
+
+    // Override
+    const scalar t = this->db().time().timeOutputValue();
+    this->refValue() = uniformInletValue_->value(t);
+}
+
+
+// * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
+
+template<class Type>
+void Foam::uniformInletOutletFvPatchField<Type>::operator=
+(
+    const fvPatchField<Type>& ptf
+)
+{
+    fvPatchField<Type>::operator=
+    (
+        this->valueFraction()*this->refValue()
+        + (1 - this->valueFraction())*ptf
+    );
+}
+
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchField.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchField.H
new file mode 100644
index 0000000000000000000000000000000000000000..cc9e2d0a9df2fa705ba05f2ca67ccd4a849a0897
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchField.H
@@ -0,0 +1,214 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  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::uniformInletOutletFvPatchField
+
+Group
+    grpOutletBoundaryConditions
+
+Description
+    Variant of inletOutlet boundary condition with uniform inletValue.
+
+    \heading Patch usage
+
+    \table
+        Property     | Description             | Required    | Default value
+        phi          | flux field name         | no          | phi
+        uniformInletValue   | inlet value for reverse flow | yes    |
+    \endtable
+
+    Example of the boundary condition specification:
+    \verbatim
+    myPatch
+    {
+        type                uniformInletOutlet;
+        phi                 phi;
+        uniformInletValue   0;
+        value               uniform 0;
+    }
+    \endverbatim
+
+    The mode of operation is determined by the sign of the flux across the
+    patch faces.
+
+Note
+    Sign conventions:
+    - positive flux (out of domain): apply zero-gradient condition
+    - negative flux (into of domain): apply the user-specified fixed value
+
+SeeAlso
+    Foam::inletOutletFvPatchField
+
+SourceFiles
+    uniformInletOutletFvPatchField.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef uniformInletOutletFvPatchField_H
+#define uniformInletOutletFvPatchField_H
+
+#include "mixedFvPatchField.H"
+#include "DataEntry.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                   Class uniformInletOutletFvPatchField Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class Type>
+class uniformInletOutletFvPatchField
+:
+    public mixedFvPatchField<Type>
+{
+
+protected:
+
+    // Protected data
+
+        //- Name of flux field
+        word phiName_;
+
+        //- Value
+        autoPtr<DataEntry<Type> > uniformInletValue_;
+
+
+public:
+
+    //- Runtime type information
+    TypeName("uniformInletOutlet");
+
+
+    // Constructors
+
+        //- Construct from patch and internal field
+        uniformInletOutletFvPatchField
+        (
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&
+        );
+
+        //- Construct from patch, internal field and dictionary
+        uniformInletOutletFvPatchField
+        (
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&,
+            const dictionary&
+        );
+
+        //- Construct by mapping given uniformInletOutletFvPatchField
+        //  onto a new patch
+        uniformInletOutletFvPatchField
+        (
+            const uniformInletOutletFvPatchField<Type>&,
+            const fvPatch&,
+            const DimensionedField<Type, volMesh>&,
+            const fvPatchFieldMapper&
+        );
+
+        //- Construct as copy
+        uniformInletOutletFvPatchField
+        (
+            const uniformInletOutletFvPatchField<Type>&
+        );
+
+        //- Construct and return a clone
+        virtual tmp<fvPatchField<Type> > clone() const
+        {
+            return tmp<fvPatchField<Type> >
+            (
+                new uniformInletOutletFvPatchField<Type>(*this)
+            );
+        }
+
+        //- Construct as copy setting internal field reference
+        uniformInletOutletFvPatchField
+        (
+            const uniformInletOutletFvPatchField<Type>&,
+            const DimensionedField<Type, volMesh>&
+        );
+
+        //- Construct and return a clone setting internal field reference
+        virtual tmp<fvPatchField<Type> > clone
+        (
+            const DimensionedField<Type, volMesh>& iF
+        ) const
+        {
+            return tmp<fvPatchField<Type> >
+            (
+                new uniformInletOutletFvPatchField<Type>(*this, iF)
+            );
+        }
+
+
+    // Mapping functions
+
+        //- Map (and resize as needed) from self given a mapping object
+        virtual void autoMap
+        (
+            const fvPatchFieldMapper&
+        );
+
+        //- Reverse map the given fvPatchField onto this fvPatchField
+        virtual void rmap
+        (
+            const fvPatchField<Type>&,
+            const labelList&
+        );
+
+
+    // Member functions
+
+        //- Update the coefficients associated with the patch field
+        virtual void updateCoeffs();
+
+        //- Write
+        virtual void write(Ostream&) const;
+
+
+    // Member operators
+
+        virtual void operator=(const fvPatchField<Type>& pvf);
+
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+#   include "uniformInletOutletFvPatchField.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFields.C b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFields.C
new file mode 100644
index 0000000000000000000000000000000000000000..c34bd1caff938a32db0211c798fd75b6b44e88a2
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFields.C
@@ -0,0 +1,44 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#include "volFields.H"
+#include "surfaceFields.H"
+#include "uniformInletOutletFvPatchFields.H"
+#include "addToRunTimeSelectionTable.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
+
+makePatchFields(uniformInletOutlet);
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFields.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFields.H
new file mode 100644
index 0000000000000000000000000000000000000000..940749210b878fade3163240c0da3cf13ebc6f5d
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFields.H
@@ -0,0 +1,49 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef uniformInletOutletFvPatchFields_H
+#define uniformInletOutletFvPatchFields_H
+
+#include "uniformInletOutletFvPatchField.H"
+#include "fieldTypes.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+makePatchTypeFieldTypedefs(uniformInletOutlet);
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFieldsFwd.H b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFieldsFwd.H
new file mode 100644
index 0000000000000000000000000000000000000000..71f272a479a08a62f01ab6df3abefc5d0cc9c5ed
--- /dev/null
+++ b/src/finiteVolume/fields/fvPatchFields/derived/uniformInletOutlet/uniformInletOutletFvPatchFieldsFwd.H
@@ -0,0 +1,50 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  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/>.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef uniformInletOutletFvPatchFieldsFwd_H
+#define uniformInletOutletFvPatchFieldsFwd_H
+
+#include "fieldTypes.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+template<class Type> class uniformInletOutletFvPatchField;
+
+makePatchTypeFieldTypedefs(uniformInletOutlet);
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/applications/utilities/mesh/manipulation/orientFaceZone/patchFaceOrientation.C b/src/mesh/autoMesh/autoHexMesh/meshRefinement/patchFaceOrientation.C
similarity index 100%
rename from applications/utilities/mesh/manipulation/orientFaceZone/patchFaceOrientation.C
rename to src/mesh/autoMesh/autoHexMesh/meshRefinement/patchFaceOrientation.C
diff --git a/applications/utilities/mesh/manipulation/orientFaceZone/patchFaceOrientation.H b/src/mesh/autoMesh/autoHexMesh/meshRefinement/patchFaceOrientation.H
similarity index 100%
rename from applications/utilities/mesh/manipulation/orientFaceZone/patchFaceOrientation.H
rename to src/mesh/autoMesh/autoHexMesh/meshRefinement/patchFaceOrientation.H
diff --git a/applications/utilities/mesh/manipulation/orientFaceZone/patchFaceOrientationI.H b/src/mesh/autoMesh/autoHexMesh/meshRefinement/patchFaceOrientationI.H
similarity index 100%
rename from applications/utilities/mesh/manipulation/orientFaceZone/patchFaceOrientationI.H
rename to src/mesh/autoMesh/autoHexMesh/meshRefinement/patchFaceOrientationI.H