diff --git a/src/thermophysicalModels/basic/Make/files b/src/thermophysicalModels/basic/Make/files
index 7d10c1c952aa1d0100331a8832804f0decbad48d..83ee06f8bf1fa7bd777a1e2f035467e2d7356641 100644
--- a/src/thermophysicalModels/basic/Make/files
+++ b/src/thermophysicalModels/basic/Make/files
@@ -6,6 +6,7 @@ psiThermo/psiThermos.C
 
 rhoThermo/rhoThermo.C
 rhoThermo/rhoThermos.C
+rhoThermo/liquidThermo.C
 
 derivedFvPatchFields/fixedEnergy/fixedEnergyFvPatchScalarField.C
 derivedFvPatchFields/gradientEnergy/gradientEnergyFvPatchScalarField.C
diff --git a/src/thermophysicalModels/basic/basicThermo/basicThermo.C b/src/thermophysicalModels/basic/basicThermo/basicThermo.C
index 9d12c6de9c1be56e52e2b99e23d503756787025f..a5acf77a3232f281f417e6e5f2f2a7090ffc6b5f 100644
--- a/src/thermophysicalModels/basic/basicThermo/basicThermo.C
+++ b/src/thermophysicalModels/basic/basicThermo/basicThermo.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) 2011-2017 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -449,14 +449,28 @@ Foam::wordList Foam::basicThermo::splitThermoName
         {
             cmpts[i] = thermoName.substr(beg, end-beg);
             cmpts[i++].replaceAll(">","");
+
+            // If the number of number of components in the name
+            // is greater than nCmpt return an empty list
+            if (i == nCmpt)
+            {
+                return wordList::null();
+            }
         }
         beg = end + 1;
     }
 
+    // If the number of number of components in the name is not equal to nCmpt
+    // return an empty list
+    if (i + 1 != nCmpt)
+    {
+        return wordList::null();
+    }
+
     if (beg < thermoName.size())
     {
         cmpts[i] = thermoName.substr(beg, string::npos);
-        cmpts[i++].replaceAll(">","");
+        cmpts[i].replaceAll(">","");
     }
 
     return cmpts;
diff --git a/src/thermophysicalModels/basic/basicThermo/basicThermo.H b/src/thermophysicalModels/basic/basicThermo/basicThermo.H
index c1d299dbb7449e15f9a0f9b1a17f865c74f6af52..5a5593886e10b5e4ece8b10e05c8bef9f87d2e03 100644
--- a/src/thermophysicalModels/basic/basicThermo/basicThermo.H
+++ b/src/thermophysicalModels/basic/basicThermo/basicThermo.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) 2011-2017 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -135,6 +135,17 @@ public:
 
     // Selectors
 
+        //- Generic lookup for thermodynamics package thermoTypeName
+        template<class Thermo, class Table>
+        static typename Table::iterator lookupThermo
+        (
+            const dictionary& thermoTypeDict,
+            Table* tablePtr,
+            const int nCmpt,
+            const char* cmptNames[],
+            const word& thermoTypeName
+        );
+
         //- Generic lookup for each of the related thermodynamics packages
         template<class Thermo, class Table>
         static typename Table::iterator lookupThermo
diff --git a/src/thermophysicalModels/basic/basicThermo/basicThermoTemplates.C b/src/thermophysicalModels/basic/basicThermo/basicThermoTemplates.C
index c9741c645b5a56b71b04eaf364c0f865dbb8579a..b589708163b1de3edc9d3776bbce54fd90f6ad24 100644
--- a/src/thermophysicalModels/basic/basicThermo/basicThermoTemplates.C
+++ b/src/thermophysicalModels/basic/basicThermo/basicThermoTemplates.C
@@ -2,7 +2,7 @@
   =========                 |
   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
    \\    /   O peration     |
-    \\  /    A nd           | Copyright (C) 2012-2015 OpenFOAM Foundation
+    \\  /    A nd           | Copyright (C) 2012-2017 OpenFOAM Foundation
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
 License
@@ -30,89 +30,153 @@ License
 template<class Thermo, class Table>
 typename Table::iterator Foam::basicThermo::lookupThermo
 (
-    const dictionary& thermoDict,
-    Table* tablePtr
+    const dictionary& thermoTypeDict,
+    Table* tablePtr,
+    const int nCmpt,
+    const char* cmptNames[],
+    const word& thermoTypeName
 )
 {
-    word thermoTypeName;
+    // Lookup the thermo package
+    typename Table::iterator cstrIter = tablePtr->find(thermoTypeName);
 
-    if (thermoDict.isDict("thermoType"))
+    // Print error message if package not found in the table
+    if (cstrIter == tablePtr->end())
     {
-        const dictionary& thermoTypeDict(thermoDict.subDict("thermoType"));
+        FatalErrorInFunction
+            << "Unknown " << Thermo::typeName << " type " << nl
+            << "thermoType" << thermoTypeDict << nl << nl
+            << "Valid " << Thermo::typeName << " types are:"
+            << nl << nl;
+
+        // Get the list of all the suitable thermo packages available
+        wordList validThermoTypeNames
+        (
+            tablePtr->sortedToc()
+        );
 
-        Info<< "Selecting thermodynamics package " << thermoTypeDict << endl;
+        // Build a table of the thermo packages constituent parts
+        // Note: row-0 contains the names of constituent parts
+        List<wordList> validThermoTypeNameCmpts
+        (
+            validThermoTypeNames.size() + 1
+        );
 
-        const int nCmpt = 7;
-        const char* cmptNames[nCmpt] =
+        validThermoTypeNameCmpts[0].setSize(nCmpt);
+        forAll(validThermoTypeNameCmpts[0], j)
         {
-            "type",
-            "mixture",
-            "transport",
-            "thermo",
-            "equationOfState",
-            "specie",
-            "energy"
-        };
-
-        // Construct the name of the thermo package from the components
-        thermoTypeName =
-            word(thermoTypeDict.lookup("type")) + '<'
-          + word(thermoTypeDict.lookup("mixture")) + '<'
-          + word(thermoTypeDict.lookup("transport")) + '<'
-          + word(thermoTypeDict.lookup("thermo")) + '<'
-          + word(thermoTypeDict.lookup("equationOfState")) + '<'
-          + word(thermoTypeDict.lookup("specie")) + ">>,"
-          + word(thermoTypeDict.lookup("energy")) + ">>>";
-
-        // Lookup the thermo package
-        typename Table::iterator cstrIter = tablePtr->find(thermoTypeName);
+            validThermoTypeNameCmpts[0][j] = cmptNames[j];
+        }
 
-        // Print error message if package not found in the table
-        if (cstrIter == tablePtr->end())
+        // Split the thermo package names into their constituent parts
+        // Removing incompatible entries from the list
+        label j = 0;
+        forAll(validThermoTypeNames, i)
         {
-            FatalErrorInFunction
-                << "Unknown " << Thermo::typeName << " type " << nl
-                << "thermoType" << thermoTypeDict << nl << nl
-                << "Valid " << Thermo::typeName << " types are:" << nl << nl;
-
-            // Get the list of all the suitable thermo packages available
-            wordList validThermoTypeNames
-            (
-                tablePtr->sortedToc()
-            );
-
-            // Build a table of the thermo packages constituent parts
-            // Note: row-0 contains the names of constituent parts
-            List<wordList> validThermoTypeNameCmpts
+            wordList names
             (
-                validThermoTypeNames.size() + 1
+                Thermo::splitThermoName(validThermoTypeNames[i], nCmpt)
             );
 
-            validThermoTypeNameCmpts[0].setSize(nCmpt);
-            forAll(validThermoTypeNameCmpts[0], j)
+            if (names.size())
             {
-                validThermoTypeNameCmpts[0][j] = cmptNames[j];
+                validThermoTypeNameCmpts[j++] = names;
             }
+        }
+        validThermoTypeNameCmpts.setSize(j);
 
-            // Split the thermo package names into their constituent parts
-            forAll(validThermoTypeNames, i)
-            {
-                validThermoTypeNameCmpts[i+1] =
-                    Thermo::splitThermoName(validThermoTypeNames[i], nCmpt);
-            }
+        // Print the table of available packages
+        // in terms of their constituent parts
+        printTable(validThermoTypeNameCmpts, FatalError);
+
+        FatalError<< exit(FatalError);
+    }
+
+    return cstrIter;
+}
+
+
+template<class Thermo, class Table>
+typename Table::iterator Foam::basicThermo::lookupThermo
+(
+    const dictionary& thermoDict,
+    Table* tablePtr
+)
+{
+    if (thermoDict.isDict("thermoType"))
+    {
+        const dictionary& thermoTypeDict(thermoDict.subDict("thermoType"));
 
-            // Print the table of available packages
-            // in terms of their constituent parts
-            printTable(validThermoTypeNameCmpts, FatalError);
+        Info<< "Selecting thermodynamics package " << thermoTypeDict << endl;
+
+        if (thermoTypeDict.found("properties"))
+        {
+            const int nCmpt = 4;
+            const char* cmptNames[nCmpt] =
+            {
+                "type",
+                "mixture",
+                "properties",
+                "energy"
+            };
+
+            // Construct the name of the thermo package from the components
+            const word thermoTypeName
+            (
+                word(thermoTypeDict.lookup("type")) + '<'
+              + word(thermoTypeDict.lookup("mixture")) + '<'
+              + word(thermoTypeDict.lookup("properties")) + ','
+              + word(thermoTypeDict.lookup("energy")) + ">>"
+            );
 
-            FatalError<< exit(FatalError);
+            return lookupThermo<Thermo, Table>
+            (
+                thermoTypeDict,
+                tablePtr,
+                nCmpt,
+                cmptNames,
+                thermoTypeName
+            );
         }
+        else
+        {
+            const int nCmpt = 7;
+            const char* cmptNames[nCmpt] =
+            {
+                "type",
+                "mixture",
+                "transport",
+                "thermo",
+                "equationOfState",
+                "specie",
+                "energy"
+            };
+
+            // Construct the name of the thermo package from the components
+            const word thermoTypeName
+            (
+                word(thermoTypeDict.lookup("type")) + '<'
+              + word(thermoTypeDict.lookup("mixture")) + '<'
+              + word(thermoTypeDict.lookup("transport")) + '<'
+              + word(thermoTypeDict.lookup("thermo")) + '<'
+              + word(thermoTypeDict.lookup("equationOfState")) + '<'
+              + word(thermoTypeDict.lookup("specie")) + ">>,"
+              + word(thermoTypeDict.lookup("energy")) + ">>>"
+            );
 
-        return cstrIter;
+            return lookupThermo<Thermo, Table>
+            (
+                thermoTypeDict,
+                tablePtr,
+                nCmpt,
+                cmptNames,
+                thermoTypeName
+            );
+        }
     }
     else
     {
-        thermoTypeName = word(thermoDict.lookup("thermoType"));
+        const word thermoTypeName(thermoDict.lookup("thermoType"));
 
         Info<< "Selecting thermodynamics package " << thermoTypeName << endl;
 
diff --git a/src/thermophysicalModels/basic/rhoThermo/liquidThermo.C b/src/thermophysicalModels/basic/rhoThermo/liquidThermo.C
new file mode 100644
index 0000000000000000000000000000000000000000..c023c506838c94851da381c0a38a6c28aea2449b
--- /dev/null
+++ b/src/thermophysicalModels/basic/rhoThermo/liquidThermo.C
@@ -0,0 +1,91 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2017 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 "rhoThermo.H"
+#include "heRhoThermo.H"
+#include "pureMixture.H"
+#include "thermo.H"
+#include "sensibleEnthalpy.H"
+#include "sensibleInternalEnergy.H"
+
+#include "thermophysicalPropertiesSelector.H"
+#include "liquidProperties.H"
+
+#include "addToRunTimeSelectionTable.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
+
+typedef heRhoThermo
+<
+    rhoThermo,
+    pureMixture
+    <
+        species::thermo
+        <
+            thermophysicalPropertiesSelector<liquidProperties>,
+            sensibleInternalEnergy
+        >
+    >
+> heRhoThermopureMixtureliquidProperties;
+
+defineTemplateTypeNameAndDebugWithName
+(
+    heRhoThermopureMixtureliquidProperties,
+    "heRhoThermo<pureMixture<liquid,sensibleInternalEnergy>>",
+    0
+);
+
+addToRunTimeSelectionTable
+(
+    basicThermo,
+    heRhoThermopureMixtureliquidProperties,
+    fvMesh
+);
+
+addToRunTimeSelectionTable
+(
+    fluidThermo,
+    heRhoThermopureMixtureliquidProperties,
+    fvMesh
+);
+
+addToRunTimeSelectionTable
+(
+    rhoThermo,
+    heRhoThermopureMixtureliquidProperties,
+    fvMesh
+);
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// ************************************************************************* //
diff --git a/src/thermophysicalModels/thermophysicalProperties/liquidProperties/H2O/H2O.C b/src/thermophysicalModels/thermophysicalProperties/liquidProperties/H2O/H2O.C
index 634007fc1db8827971921d2688cc0707668f7668..eb211a946f2e8730593ec61062577a326de9a30c 100644
--- a/src/thermophysicalModels/thermophysicalProperties/liquidProperties/H2O/H2O.C
+++ b/src/thermophysicalModels/thermophysicalProperties/liquidProperties/H2O/H2O.C
@@ -31,8 +31,6 @@ License
 namespace Foam
 {
     defineTypeNameAndDebug(H2O, 0);
-    addToRunTimeSelectionTable(thermophysicalProperties, H2O,);
-    addToRunTimeSelectionTable(thermophysicalProperties, H2O, dictionary);
     addToRunTimeSelectionTable(liquidProperties, H2O,);
     addToRunTimeSelectionTable(liquidProperties, H2O, dictionary);
 }
diff --git a/src/thermophysicalModels/thermophysicalProperties/liquidProperties/liquidProperties/liquidProperties.H b/src/thermophysicalModels/thermophysicalProperties/liquidProperties/liquidProperties/liquidProperties.H
index 939a80770a6a2a5137920bc531f80357fc767935..7df742adca56b698fde1fb1c479628e10273a548 100644
--- a/src/thermophysicalModels/thermophysicalProperties/liquidProperties/liquidProperties/liquidProperties.H
+++ b/src/thermophysicalModels/thermophysicalProperties/liquidProperties/liquidProperties/liquidProperties.H
@@ -85,7 +85,7 @@ class liquidProperties
 
 public:
 
-    TypeName("liquidProperties");
+    TypeName("liquid");
 
 
     // Declare run-time constructor selection tables
@@ -148,6 +148,15 @@ public:
     {}
 
 
+    // Static data
+
+        //- Is the equation of state is incompressible i.e. rho != f(p)
+        static const bool incompressible = true;
+
+        //- Is the equation of state is isochoric i.e. rho = const
+        static const bool isochoric = false;
+
+
     // Member Functions
 
         // Physical constants which define the specie
diff --git a/src/thermophysicalModels/thermophysicalProperties/solidProperties/solidProperties/solidProperties.H b/src/thermophysicalModels/thermophysicalProperties/solidProperties/solidProperties/solidProperties.H
index 652c81d8fd711fc94703b5c91ced63857d62a451..dce26f383e52a9b7acba5747231d6b8ed55094e2 100644
--- a/src/thermophysicalModels/thermophysicalProperties/solidProperties/solidProperties/solidProperties.H
+++ b/src/thermophysicalModels/thermophysicalProperties/solidProperties/solidProperties/solidProperties.H
@@ -72,7 +72,7 @@ class solidProperties
 public:
 
     //- Runtime type information
-    TypeName("solidProperties");
+    TypeName("solid");
 
 
     // Declare run-time constructor selection tables
diff --git a/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalProperties.H b/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalProperties.H
index 98ec171a547667d2fda7442e7e1bf32bcfba0dab..f44237ed70956df828f44b75e3da43b5e9726b62 100644
--- a/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalProperties.H
+++ b/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalProperties.H
@@ -115,12 +115,6 @@ public:
             //- Molecular weight [kg/kmol]
             inline scalar W() const;
 
-            //- Is the equation of state is incompressible i.e. rho != f(p)
-            static const bool incompressible = true;
-
-            //- Is the equation of state is isochoric i.e. rho = const
-            static const bool isochoric = false;
-
             //- Limit the temperature to be in the range Tlow_ to Thigh_
             inline scalar limit(const scalar T) const;
 
diff --git a/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalPropertiesI.H b/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalPropertiesI.H
index b6c04e00624849b5b645feae905d0a0ce1929904..87b99ce8e325976496f8741a39e1de1823527c0e 100644
--- a/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalPropertiesI.H
+++ b/src/thermophysicalModels/thermophysicalProperties/thermophysicalProperties/thermophysicalPropertiesI.H
@@ -23,6 +23,8 @@ License
 
 \*---------------------------------------------------------------------------*/
 
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
 inline Foam::scalar Foam::thermophysicalProperties::limit(const scalar T) const
 {
     return T;
diff --git a/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelector.C b/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelector.C
new file mode 100644
index 0000000000000000000000000000000000000000..59607a8160b0dfe3236929e6a54a0bb901580130
--- /dev/null
+++ b/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelector.C
@@ -0,0 +1,61 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2017 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 "thermophysicalPropertiesSelector.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class ThermophysicalProperties>
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::
+thermophysicalPropertiesSelector
+(
+    const word& name
+)
+:
+    propertiesPtr_(ThermophysicalProperties::New(name))
+{}
+
+
+template<class ThermophysicalProperties>
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::
+thermophysicalPropertiesSelector
+(
+    const dictionary& dict
+)
+{
+    const word name(dict.first()->keyword());
+
+    if (dict.isDict(name))
+    {
+        propertiesPtr_ = ThermophysicalProperties::New(dict.subDict(name));
+    }
+    else
+    {
+        propertiesPtr_ = ThermophysicalProperties::New(name);
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelector.H b/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelector.H
new file mode 100644
index 0000000000000000000000000000000000000000..690e4a950534d689aacf58899d9fc1c8fbfc5787
--- /dev/null
+++ b/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelector.H
@@ -0,0 +1,154 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2017 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::thermophysicalPropertiesSelector
+
+Description
+    Wrapper class providing run-time selection of thermophysicalProperties
+    for the templated thermodynamics packages.
+
+SourceFiles
+    thermophysicalPropertiesSelectorI.H
+    thermophysicalPropertiesSelector.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef thermophysicalPropertiesSelector_H
+#define thermophysicalPropertiesSelector_H
+
+#include "thermophysicalProperties.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                    Class thermophysicalPropertiesSelector Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class ThermophysicalProperties>
+class thermophysicalPropertiesSelector
+{
+    // Private member data
+
+        autoPtr<ThermophysicalProperties> propertiesPtr_;
+
+
+public:
+
+    // Constructors
+
+        //- Construct from name
+        thermophysicalPropertiesSelector(const word& name);
+
+        //- Construct from dictionary
+        thermophysicalPropertiesSelector(const dictionary& dict);
+
+
+    // Static data
+
+        //- Is the equation of state is incompressible i.e. rho != f(p)
+        static const bool incompressible =
+            ThermophysicalProperties::incompressible;
+
+        //- Is the equation of state is isochoric i.e. rho = const
+        static const bool isochoric =
+            ThermophysicalProperties::isochoric;
+
+
+    // Member Functions
+
+        // Physical constants which define the specie
+
+            //- Molecular weight [kg/kmol]
+            inline scalar W() const;
+
+            //- Limit the temperature to be in the range Tlow_ to Thigh_
+            inline scalar limit(const scalar T) const;
+
+
+        // Fundamental equation of state properties
+
+            //- Liquid density [kg/m^3]
+            inline scalar rho(scalar p, scalar T) const;
+
+            //- Liquid compressibility rho/p [s^2/m^2]
+            //  Note: currently it is assumed the liquid is incompressible
+            inline scalar psi(scalar p, scalar T) const;
+
+            //- Return (Cp - Cv) [J/(kg K]
+            //  Note: currently it is assumed the liquid is incompressible
+            //  so CpMCv 0
+            inline scalar CpMCv(scalar p, scalar T) const;
+
+
+        // Fundamental thermodynamic properties
+
+            //- Heat capacity at constant pressure [J/(kg K)]
+            inline scalar Cp(const scalar p, const scalar T) const;
+
+            //- Absolute Enthalpy [J/kg]
+            inline scalar Ha(const scalar p, const scalar T) const;
+
+            //- Sensible enthalpy [J/kg]
+            inline scalar Hs(const scalar p, const scalar T) const;
+
+            //- Chemical enthalpy [J/kg]
+            inline scalar Hc() const;
+
+            // Entropy [J/(kg K)]
+            inline scalar S(const scalar p, const scalar T) const;
+
+
+        // Physical properties
+
+            //- Liquid viscosity [Pa s]
+            inline scalar mu(scalar p, scalar T) const;
+
+            //- Liquid thermal conductivity  [W/(m K)]
+            inline scalar kappa(scalar p, scalar T) const;
+
+            //- Liquid thermal diffusivity of enthalpy [kg/ms]
+            inline scalar alphah(const scalar p, const scalar T) const;
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#include "thermophysicalPropertiesSelectorI.H"
+
+#ifdef NoRepository
+    #include "thermophysicalPropertiesSelector.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelectorI.H b/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelectorI.H
new file mode 100644
index 0000000000000000000000000000000000000000..870aa1dbed32356788193b0c45f498ade1d6383e
--- /dev/null
+++ b/src/thermophysicalModels/thermophysicalProperties/thermophysicalPropertiesSelector/thermophysicalPropertiesSelectorI.H
@@ -0,0 +1,175 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 2017 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/>.
+
+\*---------------------------------------------------------------------------*/
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::W() const
+{
+    return propertiesPtr_->W();
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::limit
+(
+    const scalar T
+) const
+{
+    return propertiesPtr_->limit(T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::rho
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->rho(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::psi
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->psi(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::CpMCv
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->CpMCv(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::Cp
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->Cp(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::Ha
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->Ha(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::Hs
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->Hs(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::Hc() const
+{
+    return propertiesPtr_->Hc();
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::S
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->S(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::mu
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->mu(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::kappa
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->kappa(p, T);
+}
+
+
+template<class ThermophysicalProperties>
+inline Foam::scalar
+Foam::thermophysicalPropertiesSelector<ThermophysicalProperties>::alphah
+(
+    scalar p,
+    scalar T
+) const
+{
+    return propertiesPtr_->alphah(p, T);
+}
+
+
+// ************************************************************************* //
diff --git a/tutorials/multiphase/compressibleInterFoam/laminar/depthCharge2D/constant/thermophysicalProperties.water b/tutorials/multiphase/compressibleInterFoam/laminar/depthCharge2D/constant/thermophysicalProperties.water
index 8c1db784a4e38b6df8b4c1af7b7fa434700bfe46..8e896f49c6cfc077e25f60fb4476edfa336806af 100644
--- a/tutorials/multiphase/compressibleInterFoam/laminar/depthCharge2D/constant/thermophysicalProperties.water
+++ b/tutorials/multiphase/compressibleInterFoam/laminar/depthCharge2D/constant/thermophysicalProperties.water
@@ -15,6 +15,21 @@ FoamFile
 }
 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
+thermoType
+{
+    type            heRhoThermo;
+    mixture         pureMixture;
+    properties      liquid;
+    energy          sensibleInternalEnergy;
+}
+
+mixture
+{
+    H2O;
+}
+
+
+/*
 thermoType
 {
     type            heRhoThermo;
@@ -48,6 +63,7 @@ mixture
         Pr          2.289;
     }
 }
+*/
 
 
 // ************************************************************************* //