diff --git a/src/fvOptions/Make/files b/src/fvOptions/Make/files
index 4d172dbf076029c0e44ba3f74900300199c73cdf..1692b88e3ea92e8848ec16f092ca3d1a2cd538cb 100644
--- a/src/fvOptions/Make/files
+++ b/src/fvOptions/Make/files
@@ -49,11 +49,8 @@ $(interRegion)/interRegionExplicitPorositySource/interRegionExplicitPorositySour
 
 /* Constraints */
 
-generalConstraints=constraints/general
-$(generalConstraints)/explicitSetValue/explicitSetValue.C
-
-derivedConstraints=constraints/derived
-$(derivedConstraints)/fixedTemperatureConstraint/fixedTemperatureConstraint.C
+constraints/fixedValueConstraint/fixedValueConstraints.C
+constraints/fixedTemperatureConstraint/fixedTemperatureConstraint.C
 
 
 /* Corrections */
diff --git a/src/fvOptions/constraints/derived/fixedTemperatureConstraint/fixedTemperatureConstraint.C b/src/fvOptions/constraints/fixedTemperatureConstraint/fixedTemperatureConstraint.C
similarity index 100%
rename from src/fvOptions/constraints/derived/fixedTemperatureConstraint/fixedTemperatureConstraint.C
rename to src/fvOptions/constraints/fixedTemperatureConstraint/fixedTemperatureConstraint.C
diff --git a/src/fvOptions/constraints/derived/fixedTemperatureConstraint/fixedTemperatureConstraint.H b/src/fvOptions/constraints/fixedTemperatureConstraint/fixedTemperatureConstraint.H
similarity index 100%
rename from src/fvOptions/constraints/derived/fixedTemperatureConstraint/fixedTemperatureConstraint.H
rename to src/fvOptions/constraints/fixedTemperatureConstraint/fixedTemperatureConstraint.H
diff --git a/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValue.C b/src/fvOptions/constraints/fixedValueConstraint/FixedValueConstraint.C
similarity index 60%
rename from src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValue.C
rename to src/fvOptions/constraints/fixedValueConstraint/FixedValueConstraint.C
index 271dd70cb0c38798515a90ff3ae5fafe5f8edecb..68e2ea120308970d5d0c7dc87b45b883b76059c5 100644
--- a/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValue.C
+++ b/src/fvOptions/constraints/fixedValueConstraint/FixedValueConstraint.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2016 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -23,35 +23,15 @@ License
 
 \*---------------------------------------------------------------------------*/
 
-#include "ExplicitSetValue.H"
+#include "FixedValueConstraint.H"
 #include "fvMesh.H"
 #include "fvMatrices.H"
 #include "DimensionedField.H"
 
-// * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
-
-template<class Type>
-void Foam::fv::ExplicitSetValue<Type>::setFieldData(const dictionary& dict)
-{
-    fieldNames_.setSize(dict.toc().size());
-    injectionRate_.setSize(fieldNames_.size());
-
-    applied_.setSize(fieldNames_.size(), false);
-
-    label i = 0;
-    forAllConstIter(dictionary, dict, iter)
-    {
-        fieldNames_[i] = iter().keyword();
-        dict.lookup(iter().keyword()) >> injectionRate_[i];
-        i++;
-    }
-}
-
-
 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
 
 template<class Type>
-Foam::fv::ExplicitSetValue<Type>::ExplicitSetValue
+Foam::fv::FixedValueConstraint<Type>::FixedValueConstraint
 (
     const word& name,
     const word& modelType,
@@ -59,8 +39,7 @@ Foam::fv::ExplicitSetValue<Type>::ExplicitSetValue
     const fvMesh& mesh
 )
 :
-    cellSetOption(name, modelType, dict, mesh),
-    injectionRate_()
+    cellSetOption(name, modelType, dict, mesh)
 {
     read(dict);
 }
@@ -69,21 +48,47 @@ Foam::fv::ExplicitSetValue<Type>::ExplicitSetValue
 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
 
 template<class Type>
-void Foam::fv::ExplicitSetValue<Type>::constrain
+bool Foam::fv::FixedValueConstraint<Type>::read(const dictionary& dict)
+{
+    if (cellSetOption::read(dict))
+    {
+        const dictionary& fieldValuesDict = coeffs_.subDict("fieldValues");
+
+        fieldNames_.setSize(fieldValuesDict.size());
+        fieldValues_.setSize(fieldNames_.size());
+
+        label i = 0;
+        forAllConstIter(dictionary, fieldValuesDict, iter)
+        {
+            fieldNames_[i] = iter().keyword();
+            fieldValuesDict.lookup(iter().keyword()) >> fieldValues_[i];
+            i++;
+        }
+
+        applied_.setSize(fieldNames_.size(), false);
+
+        return true;
+    }
+    else
+    {
+        return false;
+    }
+}
+
+
+template<class Type>
+void Foam::fv::FixedValueConstraint<Type>::constrain
 (
     fvMatrix<Type>& eqn,
     const label fieldi
 )
 {
-    if (debug)
-    {
-        Info<< "ExplicitSetValue<"<< pTraits<Type>::typeName
-            << ">::constrain for source " << name_ << endl;
-    }
-
-    List<Type> values(cells_.size(), injectionRate_[fieldi]);
+    DebugInfo
+        << "FixedValueConstraint<"
+        << pTraits<Type>::typeName
+        << ">::constrain for source " << name_ << endl;
 
-    eqn.setValues(cells_, values);
+    eqn.setValues(cells_, List<Type>(cells_.size(), fieldValues_[fieldi]));
 }
 
 
diff --git a/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValue.H b/src/fvOptions/constraints/fixedValueConstraint/FixedValueConstraint.H
similarity index 66%
rename from src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValue.H
rename to src/fvOptions/constraints/fixedValueConstraint/FixedValueConstraint.H
index 5ad58ef3cf86fa18e1f25328be549cedba2e2f0b..4fc55ea5aa807b260be71afe20d0df95e4d5006b 100644
--- a/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValue.H
+++ b/src/fvOptions/constraints/fixedValueConstraint/FixedValueConstraint.H
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2016 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -22,20 +22,27 @@ License
     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
 
 Class
-    Foam::fv::explicitSetValue
+    Foam::fv::FixedValueConstraint
 
 Description
-    Set values field values explicity.
-
-    Sources described by:
+    Constrain the field values within a specified region.
 
+    For example to set the turbulence properties within a porous region:
     \verbatim
-    <Type>ExplicitSetValueCoeffs
+    porosityTurbulence
     {
-        injectionRate
+        type            scalarFixedValueConstraint;
+        active          yes;
+
+        scalarFixedValueConstraintCoeffs
         {
-            k           30.7;
-            epsilon     1.5;
+            selectionMode   cellZone;
+            cellZone        porosity;
+            fieldValues
+            {
+                k           1;
+                epsilon     150;
+            }
         }
     }
     \endverbatim
@@ -44,15 +51,15 @@ SeeAlso
     Foam::fvOption
 
 SourceFiles
-    explicitSetValue.C
+    FixedValueConstraint.C
+    fixedValueConstraints.C
 
 \*---------------------------------------------------------------------------*/
 
-#ifndef ExplicitSetValue_H
-#define ExplicitSetValue_H
+#ifndef FixedValueConstraint_H
+#define FixedValueConstraint_H
 
 #include "cellSetOption.H"
-#include "Tuple2.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
@@ -62,39 +69,30 @@ namespace fv
 {
 
 /*---------------------------------------------------------------------------*\
-                       Class explicitSetValue Declaration
+                       Class FixedValueConstraint Declaration
 \*---------------------------------------------------------------------------*/
 
 template<class Type>
-class ExplicitSetValue
+class FixedValueConstraint
 :
     public cellSetOption
 {
+    // Private member data
 
-protected:
-
-    // Protected data
-
-        //- Source value per field
-        List<Type> injectionRate_;
-
-
-    // Protected functions
-
-        //- Set the local field data
-        void setFieldData(const dictionary& dict);
+        //- Field values
+        List<Type> fieldValues_;
 
 
 public:
 
     //- Runtime type information
-    TypeName("explicitSetValue");
+    TypeName("FixedValueConstraint");
 
 
     // Constructors
 
         //- Construct from components
-        ExplicitSetValue
+        FixedValueConstraint
         (
             const word& name,
             const word& modelType,
@@ -105,16 +103,11 @@ public:
 
     // Member Functions
 
-        // Evaluation
-
-            //- Set value on field
-            virtual void constrain(fvMatrix<Type>& eqn, const label fieldi);
-
-
-        // IO
+        //- Read source dictionary
+        virtual bool read(const dictionary& dict);
 
-            //- Read source dictionary
-            virtual bool read(const dictionary& dict);
+        //- Set value on field
+        virtual void constrain(fvMatrix<Type>& eqn, const label fieldi);
 };
 
 
@@ -126,8 +119,7 @@ public:
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
 #ifdef NoRepository
-    #include "ExplicitSetValue.C"
-    #include "ExplicitSetValueIO.C"
+    #include "FixedValueConstraint.C"
 #endif
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
diff --git a/src/fvOptions/constraints/general/explicitSetValue/explicitSetValue.C b/src/fvOptions/constraints/fixedValueConstraint/fixedValueConstraints.C
similarity index 79%
rename from src/fvOptions/constraints/general/explicitSetValue/explicitSetValue.C
rename to src/fvOptions/constraints/fixedValueConstraint/fixedValueConstraints.C
index 2e4dea145faf003d583d146cc6f9c30c4d217445..bb593e49a4aad03ee326e920a36459ae4ef00f37 100644
--- a/src/fvOptions/constraints/general/explicitSetValue/explicitSetValue.C
+++ b/src/fvOptions/constraints/fixedValueConstraint/fixedValueConstraints.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011-2015 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2016 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -24,15 +24,14 @@ License
 \*---------------------------------------------------------------------------*/
 
 #include "makeFvOption.H"
-#include "ExplicitSetValue.H"
+#include "FixedValueConstraint.H"
 
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
-makeFvOption(ExplicitSetValue, scalar);
-makeFvOption(ExplicitSetValue, vector);
-makeFvOption(ExplicitSetValue, sphericalTensor);
-makeFvOption(ExplicitSetValue, symmTensor);
-makeFvOption(ExplicitSetValue, tensor);
-
+makeFvOption(FixedValueConstraint, scalar);
+makeFvOption(FixedValueConstraint, vector);
+makeFvOption(FixedValueConstraint, sphericalTensor);
+makeFvOption(FixedValueConstraint, symmTensor);
+makeFvOption(FixedValueConstraint, tensor);
 
 // ************************************************************************* //
diff --git a/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValueIO.C b/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValueIO.C
deleted file mode 100644
index 0f1c2f9d4939adffc8a605c6e9e7600bdeb7e864..0000000000000000000000000000000000000000
--- a/src/fvOptions/constraints/general/explicitSetValue/ExplicitSetValueIO.C
+++ /dev/null
@@ -1,45 +0,0 @@
-/*---------------------------------------------------------------------------*\
-  =========                 |
-  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
-   \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2011-2015 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 "ExplicitSetValue.H"
-
-// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
-
-template<class Type>
-bool Foam::fv::ExplicitSetValue<Type>::read(const dictionary& dict)
-{
-    if (cellSetOption::read(dict))
-    {
-        setFieldData(coeffs_.subDict("injectionRate"));
-        return true;
-    }
-    else
-    {
-        return false;
-    }
-}
-
-
-// ************************************************************************* //
diff --git a/tutorials/compressible/rhoSimpleFoam/angledDuctExplicitFixedCoeff/constant/fvOptions b/tutorials/compressible/rhoSimpleFoam/angledDuctExplicitFixedCoeff/constant/fvOptions
index f9bc8b076bafce50be361a2bfb2366efe4124302..bb431dc10a7f67d060db3e8f3bdbb0b7db1f9604 100644
--- a/tutorials/compressible/rhoSimpleFoam/angledDuctExplicitFixedCoeff/constant/fvOptions
+++ b/tutorials/compressible/rhoSimpleFoam/angledDuctExplicitFixedCoeff/constant/fvOptions
@@ -15,22 +15,7 @@ FoamFile
 }
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
-fixedTemperaure1
-{
-    type            fixedTemperatureConstraint;
-    active          yes;
-
-    fixedTemperatureConstraintCoeffs
-    {
-        selectionMode   cellZone;
-        cellZone        porosity;
-        mode            uniform;
-        temperature     350;
-    }
-}
-
-
-porosity1
+porosity
 {
     type            explicitPorositySource;
     active          yes;
@@ -65,4 +50,37 @@ porosity1
 }
 
 
+fixedTemperature
+{
+    type            fixedTemperatureConstraint;
+    active          yes;
+
+    fixedTemperatureConstraintCoeffs
+    {
+        selectionMode   cellZone;
+        cellZone        porosity;
+        mode            uniform;
+        temperature     350;
+    }
+}
+
+
+porosityTurbulence
+{
+    type            scalarFixedValueConstraint;
+    active          yes;
+
+    scalarFixedValueConstraintCoeffs
+    {
+        selectionMode   cellZone;
+        cellZone        porosity;
+        fieldValues
+        {
+            k           1;
+            epsilon     150;
+        }
+    }
+}
+
+
 // ************************************************************************* //