diff --git a/dumux/material/components/brine.hh b/dumux/material/components/brine.hh
index 74d3dce37ee1710cb0238e7973f197d52be974fd..a6d4bcc60d61c5629a07374bc7e7269d3c748465 100644
--- a/dumux/material/components/brine.hh
+++ b/dumux/material/components/brine.hh
@@ -171,30 +171,11 @@ public:
         /* heat of dissolution for halite according to Michaelides 1971 */
         const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
 
-        /* Enthalpy of brine without air */
+        /* Enthalpy of brine without any dissolved gas */
         const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h; /* kJ/kg */
         return h_ls1*1E3; /*J/kg*/
     }
 
-    /*!
-     * \brief Specific isobaric heat capacity of liquid water \f$\mathrm{[J/kg]}\f$.
-     *
-     * \param temperature temperature of component in \f$\mathrm{[K]}\f$
-     * \param pressure pressure of component in \f$\mathrm{[Pa]}\f$
-     *
-     * See:
-     *
-     * IAPWS: "Revised Release on the IAPWS Industrial Formulation
-     * 1997 for the Thermodynamic Properties of Water and Steam",
-     * http://www.iapws.org/relguide/IF97-Rev.pdf  \cite IAPWS1997
-     */
-     DUNE_DEPRECATED_MSG("liquidHeatCapacity(Scalar temperature, Scalar pressure) is deprecated. Use liquidHeatCapacity(Scalar temperature, Scalar pressure, Scalar salinity) instead.")
-    static const Scalar liquidHeatCapacity(Scalar temperature,
-                                        Scalar pressure)
-    {
-        return liquidHeatCapacity(temperature, pressure, constantSalinity);
-    }
-
     /*!
      * \brief Specific isobaric heat capacity of liquid water \f$\mathrm{[J/kg]}\f$.
      *
@@ -313,7 +294,7 @@ public:
                     300*pMPa -
                     2400*pMPa*salinity +
                     TempC*(
-                        80.0 -
+                        80.0 +
                         3*TempC -
                         3300*salinity -
                         13*pMPa +
diff --git a/dumux/material/fluidsystems/brineco2.hh b/dumux/material/fluidsystems/brineco2.hh
index 01104058b5e5930d6233437478a6718fa454a5c4..8b448f56cf5b37c798e7cb860298d561d2183a01 100644
--- a/dumux/material/fluidsystems/brineco2.hh
+++ b/dumux/material/fluidsystems/brineco2.hh
@@ -82,8 +82,8 @@ public:
     static const int numComponents = 2;
     static const int numPhases = 2;
 
-    static const int lPhaseIdx = 0; // index of the liquid phase
-    static const int gPhaseIdx = 1; // index of the gas phase
+    static const int wPhaseIdx = 0; // index of the liquid phase
+    static const int nPhaseIdx = 1; // index of the gas phase
     static const int wCompIdx = 0;
     static const int nCompIdx = 1;
     static const int lCompIdx = wCompIdx;
@@ -117,7 +117,7 @@ public:
     {
         assert(0 <= phaseIdx && phaseIdx < numPhases);
 
-        return phaseIdx != gPhaseIdx;
+        return phaseIdx != nPhaseIdx;
     }
 
     /*!
@@ -210,7 +210,7 @@ public:
                                 startPressure, endPressure, pressureSteps);
         }
         // set the salinity of brine
-        BrineRawComponent::salinity = salinity;
+        BrineRawComponent::constantSalinity = salinity;
 
         if(Brine::isTabulated)
         {
@@ -238,12 +238,12 @@ public:
         Scalar temperature = fluidState.temperature(phaseIdx);
         Scalar pressure = fluidState.pressure(phaseIdx);
 
-        if (phaseIdx == lPhaseIdx) {
+        if (phaseIdx == wPhaseIdx) {
             // use normalized composition for to calculate the density
             // (the relations don't seem to take non-normalized
             // compositions too well...)
-            Scalar xlBrine = std::min(1.0, std::max(0.0, fluidState.moleFraction(lPhaseIdx, BrineIdx)));
-            Scalar xlCO2 = std::min(1.0, std::max(0.0, fluidState.moleFraction(lPhaseIdx, CO2Idx)));
+            Scalar xlBrine = std::min(1.0, std::max(0.0, fluidState.moleFraction(wPhaseIdx, BrineIdx)));
+            Scalar xlCO2 = std::min(1.0, std::max(0.0, fluidState.moleFraction(wPhaseIdx, CO2Idx)));
             Scalar sumx = xlBrine + xlCO2;
             xlBrine /= sumx;
             xlCO2 /= sumx;
@@ -257,13 +257,13 @@ public:
             return result;
         }
         else {
-            assert(phaseIdx == gPhaseIdx);
+            assert(phaseIdx == nPhaseIdx);
 
             // use normalized composition for to calculate the density
             // (the relations don't seem to take non-normalized
             // compositions too well...)
-            Scalar xgBrine = std::min(1.0, std::max(0.0, fluidState.moleFraction(gPhaseIdx, BrineIdx)));
-            Scalar xgCO2 = std::min(1.0, std::max(0.0, fluidState.moleFraction(gPhaseIdx, CO2Idx)));
+            Scalar xgBrine = std::min(1.0, std::max(0.0, fluidState.moleFraction(nPhaseIdx, BrineIdx)));
+            Scalar xgCO2 = std::min(1.0, std::max(0.0, fluidState.moleFraction(nPhaseIdx, CO2Idx)));
             Scalar sumx = xgBrine + xgCO2;
             xgBrine /= sumx;
             xgCO2 /= sumx;
@@ -294,7 +294,7 @@ public:
         Scalar pressure = fluidState.pressure(phaseIdx);
         Scalar result = 0;
 
-        if (phaseIdx == lPhaseIdx) {
+        if (phaseIdx == wPhaseIdx) {
             // assume pure brine for the liquid phase. TODO: viscosity
             // of mixture
             result = Brine::liquidViscosity(temperature, pressure);
@@ -341,7 +341,7 @@ public:
         assert(0 <= phaseIdx && phaseIdx < numPhases);
         assert(0 <= compIdx && compIdx < numComponents);
 
-        if (phaseIdx == gPhaseIdx)
+        if (phaseIdx == nPhaseIdx)
             // use the fugacity coefficients of an ideal gas. the
             // actual value of the fugacity is not relevant, as long
             // as the relative fluid compositions are observed,
@@ -408,12 +408,12 @@ public:
         // temperature and pressure.
         Brine_CO2::calculateMoleFractions(temperature,
                                                   pressure,
-                                                  BrineRawComponent::salinity,
+                                                  BrineRawComponent::constantSalinity,
                                                   /*knownPhaseIdx=*/-1,
                                                   xlCO2,
                                                   xgH2O);
 
-        if(phaseIdx == gPhaseIdx)
+        if(phaseIdx == nPhaseIdx)
         {
             return xgH2O;
         }
@@ -482,7 +482,7 @@ public:
 
         Scalar temperature = fluidState.temperature(phaseIdx);
         Scalar pressure = fluidState.pressure(phaseIdx);
-        if (phaseIdx == lPhaseIdx) {
+        if (phaseIdx == wPhaseIdx) {
             assert(compIIdx == BrineIdx);
             assert(compJIdx == CO2Idx);
 
@@ -491,7 +491,7 @@ public:
             return result;
         }
         else {
-            assert(phaseIdx == gPhaseIdx);
+            assert(phaseIdx == nPhaseIdx);
             assert(compIIdx == BrineIdx);
             assert(compJIdx == CO2Idx);
 
@@ -517,12 +517,12 @@ public:
         Scalar temperature = fluidState.temperature(phaseIdx);
         Scalar pressure = fluidState.pressure(phaseIdx);
 
-        if (phaseIdx == lPhaseIdx) {
+        if (phaseIdx == wPhaseIdx) {
             Scalar XlCO2 = fluidState.massFraction(phaseIdx, CO2Idx);
 
             Scalar result = liquidEnthalpyBrineCO2_(temperature,
                                                     pressure,
-                                                    BrineRawComponent::salinity,
+                                                    BrineRawComponent::constantSalinity,
                                                     XlCO2);
             Valgrind::CheckDefined(result);
             return result;
@@ -531,10 +531,10 @@ public:
             Scalar result = 0;
             result +=
                 Brine::gasEnthalpy(temperature, pressure) *
-                fluidState.massFraction(gPhaseIdx, BrineIdx);
+                fluidState.massFraction(nPhaseIdx, BrineIdx);
             result +=
                 CO2::gasEnthalpy(temperature, pressure) *
-                fluidState.massFraction(gPhaseIdx, CO2Idx);
+                fluidState.massFraction(nPhaseIdx, CO2Idx);
             Valgrind::CheckDefined(result);
             return result;
         }
@@ -551,7 +551,7 @@ public:
     static Scalar thermalConductivity(const FluidState &fluidState,
                                       int phaseIdx)
     {
-        if (phaseIdx == lPhaseIdx)
+        if (phaseIdx == wPhaseIdx)
         {
             return H2O::liquidThermalConductivity(fluidState.temperature(phaseIdx),
                                                   fluidState.pressure(phaseIdx));
@@ -576,7 +576,7 @@ public:
     static Scalar heatCapacity(const FluidState &fluidState,
                                int phaseIdx)
     {
-        if(phaseIdx == lPhaseIdx)
+        if(phaseIdx == wPhaseIdx)
             return H2O::liquidHeatCapacity(fluidState.temperature(phaseIdx),
                                            fluidState.pressure(phaseIdx));
         else
@@ -662,69 +662,22 @@ private:
     {
         /* X_CO2_w : mass fraction of CO2 in brine */
 
-        /* same function as enthalpy_brine, only extended by CO2 content */
-
-        /*Numerical coefficents from PALLISER*/
-        static const Scalar f[] = {
-            2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10
-        };
-
-        /*Numerical coefficents from MICHAELIDES for the enthalpy of brine*/
-        static const Scalar a[4][3] = {
-            { 9633.6, -4080.0, +286.49 },
-            { +166.58, +68.577, -4.6856 },
-            { -0.90963, -0.36524, +0.249667E-1 },
-            { +0.17965E-2, +0.71924E-3, -0.4900E-4 }
-        };
-
-        Scalar theta, h_NaCl;
-        Scalar m, h_ls, h_ls1, d_h;
-        Scalar S_lSAT, delta_h;
-        int i, j;
-        Scalar delta_hCO2, hg, hw;
-
-        theta = T - 273.15;
-
-        S_lSAT = f[0] + f[1]*theta + f[2]*theta*theta + f[3]*theta*theta*theta;
-        /*Regularization*/
-        if (S>S_lSAT) {
-            S = S_lSAT;
-        }
-
-        hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
-
-        /*DAUBERT and DANNER*/
-        /*U=*/h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
-                        +((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* kJ/kg */
-
-        m = (1E3/58.44)*(S/(1-S));
-        i = 0;
-        j = 0;
-        d_h = 0;
-
-        for (i = 0; i<=3; i++) {
-            for (j=0; j<=2; j++) {
-                d_h = d_h + a[i][j] * pow(theta, i) * pow(m, j);
-            }
-        }
-        /* heat of dissolution for halite according to Michaelides 1971 */
-        delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
-
-        /* Enthalpy of brine without CO2 */
-        h_ls1 =(1-S)*hw + S*h_NaCl + S*delta_h; /* kJ/kg */
+        const Scalar h_ls1 = BrineRawComponent::liquidEnthalpy(T, p, S)/1E3; /* J/kg */
 
         /* heat of dissolution for CO2 according to Fig. 6 in Duan and Sun 2003. (kJ/kg)
            In the relevant temperature ranges CO2 dissolution is
            exothermal */
-        delta_hCO2 = (-57.4375 + T * 0.1325) * 1000/44;
+        const Scalar delta_hCO2 = (-57.4375 + T * 0.1325) * 1000/44;
+
+        const Scalar hw = H2O::liquidEnthalpy(T, p) /1E3; /* kJ/kg */
 
         /* enthalpy contribution of CO2 (kJ/kg) */
-        hg = CO2::liquidEnthalpy(T, p)/1E3 + delta_hCO2;
+        const Scalar hg = CO2::liquidEnthalpy(T, p)/1E3 + delta_hCO2;
 
         /* Enthalpy of brine with dissolved CO2 */
-        h_ls = (h_ls1 - X_CO2_w*hw + hg*X_CO2_w)*1E3; /*J/kg*/
+        const Scalar h_ls = (h_ls1 - X_CO2_w*hw + hg*X_CO2_w)*1E3; /*J/kg*/
 
-        return (h_ls);
+        return h_ls;
     }
 };
 } // end namespace FluidSystems