diff --git a/dumux/freeflow/nonisothermal/indices.hh b/dumux/freeflow/nonisothermal/indices.hh
index e434a97266961904a29d6f7bb6f2acf567a1044d..cb73441208d435f20ff0f1514c5692e1eae7085e 100644
--- a/dumux/freeflow/nonisothermal/indices.hh
+++ b/dumux/freeflow/nonisothermal/indices.hh
@@ -33,17 +33,19 @@ namespace Dumux
  * \ingroup NavierStokesNIModel
  * \brief Indices for the non-isothermal Navier-Stokes model.
  *
+ * \tparam dimension The dimension of the problem
+ * \tparam numEquations The number of model equations
+ * \tparam IsothermalIndices The isothermal indices class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
-class NavierStokesNonIsothermalIndices : public GET_PROP_TYPE(TypeTag, IsothermalIndices)
+template <int dimension, int numEquations, class IsothermalIndices, int PVOffset = 0>
+class NavierStokesNonIsothermalIndices : public IsothermalIndices
 {
 public:
-    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
-    static constexpr auto dim = GridView::dimension;
-    static constexpr auto numEq = GET_PROP_VALUE(TypeTag, NumEq);
+    static constexpr auto dim = dimension;
+    static constexpr auto numEq = numEquations;
 
-    static constexpr auto energyBalanceIdx = PVOffset + GET_PROP_VALUE(TypeTag, NumEq) - dim - 1;
+    static constexpr auto energyBalanceIdx = PVOffset + numEq - dim - 1;
     static constexpr int temperatureIdx = energyBalanceIdx;
 };
 } // end namespace
diff --git a/dumux/freeflow/nonisothermal/model.hh b/dumux/freeflow/nonisothermal/model.hh
index f58a6ffc03431f7a5c727fb9ba7fa0a06025813c..1fba624a3bcf115eec33f84f75cd0d164666c759 100644
--- a/dumux/freeflow/nonisothermal/model.hh
+++ b/dumux/freeflow/nonisothermal/model.hh
@@ -65,7 +65,15 @@ public:
 SET_BOOL_PROP(NavierStokesNonIsothermal, EnableEnergyBalance, true);
 
 //! The non-isothermal indices
-SET_TYPE_PROP(NavierStokesNonIsothermal, Indices, NavierStokesNonIsothermalIndices<TypeTag>);
+SET_PROP(NavierStokesNonIsothermal, Indices)
+{
+private:
+    static constexpr int dim = GET_PROP_TYPE(TypeTag, GridView)::dimension;
+    static constexpr int numEq = GET_PROP_VALUE(TypeTag, NumEq);
+    using IsothermalIndices = typename GET_PROP_TYPE(TypeTag, IsothermalIndices);
+public:
+    using type = NavierStokesNonIsothermalIndices<dim, numEq, IsothermalIndices>;
+};
 
 //! The non-isothermal vtk output fields
 SET_TYPE_PROP(NavierStokesNonIsothermal, VtkOutputFields, NavierStokesNonIsothermalVtkOutputFields<TypeTag>);
diff --git a/dumux/porousmediumflow/1pnc/indices.hh b/dumux/porousmediumflow/1pnc/indices.hh
index c354a534bf77fac4fd39895f6e4b33db6487f50b..925e9f52b6176baad7a8fa5b8ceaa97b923713c0 100644
--- a/dumux/porousmediumflow/1pnc/indices.hh
+++ b/dumux/porousmediumflow/1pnc/indices.hh
@@ -34,13 +34,14 @@ namespace Dumux
  * \ingroup OnePNCModel
  * \brief The indices for the isothermal one-phase n-component model.
  *
+ * \tparam phaseIndex The default phase index
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <int phaseIndex, int PVOffset = 0>
 struct OnePNCIndices
 {
     //! Set the default phase used by the fluid system to the first one
-    static const int phaseIdx = GET_PROP_VALUE(TypeTag, PhaseIdx);
+    static const int phaseIdx = phaseIndex;
 
     //! Component indices
     static const int phaseCompIdx = phaseIdx;//!< The index of the main component of the considered phase
diff --git a/dumux/porousmediumflow/1pnc/model.hh b/dumux/porousmediumflow/1pnc/model.hh
index 87d143f59f1898c30c4e0250914774bc9ca9bdd9..dc0c0959790b57c4a69e101cfa33b4551cf184a4 100644
--- a/dumux/porousmediumflow/1pnc/model.hh
+++ b/dumux/porousmediumflow/1pnc/model.hh
@@ -149,9 +149,17 @@ SET_TYPE_PROP(OnePNC, VolumeVariables, OnePNCVolumeVariables<TypeTag>);   //!< t
 SET_BOOL_PROP(OnePNC, EnableAdvection, true);                           //!< The one-phase model considers advection
 SET_BOOL_PROP(OnePNC, EnableMolecularDiffusion, true);                 //!< The one-phase model has no molecular diffusion
 SET_BOOL_PROP(OnePNC, EnableEnergyBalance, false);                      //!< Isothermal model by default
-SET_TYPE_PROP(OnePNC, Indices, OnePNCIndices <TypeTag, /*PVOffset=*/0>);                            //!< The indices required by the isothermal single-phase model
 SET_TYPE_PROP(OnePNC, VtkOutputFields, OnePNCVtkOutputFields<TypeTag>);   //!< Set the vtk output fields specific to this model
 
+//! The indices required by the isothermal single-phase model
+SET_PROP(OnePNC, Indices)
+{
+private:
+    static constexpr int phaseIdx = GET_PROP_VALUE(TypeTag, PhaseIdx);
+public:
+    using type = OnePNCIndices<phaseIdx>;
+};
+
 
 ///////////////////////////////////////////////////////////////////////////
 // properties for the non-isothermal single phase model
@@ -171,11 +179,19 @@ SET_BOOL_PROP(OnePNCNI, EnableEnergyBalance, true);
 SET_TYPE_PROP(OnePNCNI, IsothermalVtkOutputFields, OnePNCVtkOutputFields<TypeTag>);     //!< the isothermal vtk output fields
 SET_TYPE_PROP(OnePNCNI, IsothermalVolumeVariables, OnePNCVolumeVariables<TypeTag>);     //!< Vol vars of the isothermal model
 SET_TYPE_PROP(OnePNCNI, IsothermalLocalResidual, CompositionalLocalResidual<TypeTag>);   //!< Local residual of the isothermal model
-SET_TYPE_PROP(OnePNCNI, IsothermalIndices, OnePNCIndices <TypeTag, /*PVOffset=*/0>);                              //!< Indices of the isothermal model
 SET_TYPE_PROP(OnePNCNI,
               ThermalConductivityModel,
               ThermalConductivityAverage<typename GET_PROP_TYPE(TypeTag, Scalar)>); //!< Use the average for effective conductivities
 
+//! Indices of the isothermal model
+SET_PROP(OnePNCNI, IsothermalIndices)
+{
+private:
+    static constexpr int phaseIdx = GET_PROP_VALUE(TypeTag, PhaseIdx);
+public:
+    using type = OnePNCIndices<phaseIdx>;
+};
+
 } // end namespace Properties
 } // end namespace Dumux
 
diff --git a/dumux/porousmediumflow/1pncmin/model.hh b/dumux/porousmediumflow/1pncmin/model.hh
index 072e7b09fd09e0815c9ff1757f32868bb2de054a..5e34fbddb9083cbc40c4025bfa0bbba9e8fdfe7b 100644
--- a/dumux/porousmediumflow/1pncmin/model.hh
+++ b/dumux/porousmediumflow/1pncmin/model.hh
@@ -93,7 +93,6 @@ SET_TYPE_PROP(OnePNCMin, NonMineralizationVtkOutputFields, OnePNCVtkOutputFields
 SET_TYPE_PROP(OnePNCMinNI, IsothermalVolumeVariables, MineralizationVolumeVariables<TypeTag>);  //!< set isothermal VolumeVariables
 SET_TYPE_PROP(OnePNCMinNI, IsothermalVtkOutputFields, MineralizationVtkOutputFields<TypeTag>);  //!< set isothermal output fields
 SET_TYPE_PROP(OnePNCMinNI, IsothermalLocalResidual, MineralizationLocalResidual<TypeTag>);      //!< set isothermal output fields
-SET_TYPE_PROP(OnePNCMinNI, IsothermalIndices, OnePNCIndices <TypeTag, /*PVOffset=*/0>);         //!< use 1pnc indices for the isothermal indices
 
 SET_TYPE_PROP(OnePNCMinNI,
               ThermalConductivityModel,
@@ -106,6 +105,15 @@ public:
     static const int value = FluidSystem::numComponents + FluidSystem::numSPhases;
 };
 
+//! use 1pnc indices for the isothermal indices
+SET_PROP(OnePNCMinNI, IsothermalIndices)
+{
+private:
+    static constexpr int phaseIdx = GET_PROP_VALUE(TypeTag, PhaseIdx);
+public:
+    using type = OnePNCIndices<phaseIdx>;
+};
+
 } // end namespace Properties
 } // end namespace Dumux
 
diff --git a/dumux/porousmediumflow/2p/indices.hh b/dumux/porousmediumflow/2p/indices.hh
index 53b11b6bade8ad1a1018917a5f6820a3eb8825dc..08c97b4994a22cf40fb60168138145be5e39c3c9 100644
--- a/dumux/porousmediumflow/2p/indices.hh
+++ b/dumux/porousmediumflow/2p/indices.hh
@@ -45,14 +45,12 @@ struct TwoPFormulation
  * \ingroup TwoPModel
  * \brief Defines the indices required for the two-phase fully implicit model.
  *
- * \tparam TypeTag The problem type tag
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <class FluidSystem, int PVOffset = 0>
 struct TwoPCommonIndices
 {
-    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
-
     // Phase indices
     static const int wPhaseIdx = FluidSystem::wPhaseIdx; //!< index of the wetting phase
     static const int nPhaseIdx = FluidSystem::nPhaseIdx; //!< index of the non-wetting phase
@@ -73,15 +71,15 @@ struct TwoPCommonIndices
  * \brief The indices for the \f$p_w-S_n\f$ formulation of the
  *        isothermal two-phase model.
  *
- * \tparam TypeTag The problem type tag
+ * \tparam FluidSystem The fluid system class
  * \tparam formulation The formulation, either pwsn or pnsw
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag,
+template <class FluidSystem,
           int formulation = TwoPFormulation::pwsn,
           int PVOffset = 0>
 struct TwoPIndices
-: public TwoPCommonIndices<TypeTag, PVOffset>, TwoPFormulation
+: public TwoPCommonIndices<FluidSystem, PVOffset>, TwoPFormulation
 {
     // indices of the primary variables
     static const int pwIdx = PVOffset + 0; //!< index of the wetting phase pressure
@@ -93,12 +91,12 @@ struct TwoPIndices
  * \brief The indices for the \f$p_n-S_w\f$ formulation of the
  *        isothermal two-phase model.
  *
- * \tparam TypeTag The problem type tag
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset>
-struct TwoPIndices<TypeTag, TwoPFormulation::pnsw, PVOffset>
-: public TwoPCommonIndices<TypeTag, PVOffset>, TwoPFormulation
+template <class FluidSystem, int PVOffset>
+struct TwoPIndices<FluidSystem, TwoPFormulation::pnsw, PVOffset>
+: public TwoPCommonIndices<FluidSystem, PVOffset>, TwoPFormulation
 {
     // indices of the primary variables
     static const int pnIdx = PVOffset + 0; //!< index of the nonwetting phase pressure
diff --git a/dumux/porousmediumflow/2p/model.hh b/dumux/porousmediumflow/2p/model.hh
index cdba9ba5c0445bae4a82f267b6f51c791cc18dc8..5c04c81edd98ecb6f7ff9eb3d4e7a210c22d9284 100644
--- a/dumux/porousmediumflow/2p/model.hh
+++ b/dumux/porousmediumflow/2p/model.hh
@@ -107,9 +107,15 @@ SET_TYPE_PROP(TwoP, SpatialParams, FVSpatialParams<TypeTag>);                 //
 //! Set the vtk output fields specific to the twop model
 SET_TYPE_PROP(TwoP, VtkOutputFields, TwoPVtkOutputFields<typename GET_PROP_TYPE(TypeTag, Indices)>);
 
-SET_TYPE_PROP(TwoP,
-              Indices,
-              TwoPIndices<TypeTag, GET_PROP_VALUE(TypeTag, Formulation), 0>); //!< The indices required by the isothermal 2p model
+//! The indices required by the isothermal 2p model
+SET_PROP(TwoP, Indices)
+{
+private:
+    using Fluidsystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+    static constexpr int formulation = GET_PROP_VALUE(TypeTag, Formulation);
+public:
+    using type = TwoPIndices<Fluidsystem, formulation, 0>;
+};
 
 //! The two-phase model uses the immiscible fluid state
 SET_PROP(TwoP, FluidState)
@@ -135,9 +141,10 @@ SET_TYPE_PROP(TwoPNI, IsothermalVtkOutputFields, TwoPVtkOutputFields<typename GE
 SET_PROP(TwoPNI, IsothermalIndices)
 {
 private:
-    enum { Formulation = GET_PROP_VALUE(TypeTag, Formulation) };
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+    static constexpr int formulation = GET_PROP_VALUE(TypeTag, Formulation);
 public:
-    using type = TwoPIndices<TypeTag, Formulation, 0>;
+    using type = TwoPIndices<FluidSystem, formulation, 0>;
 };
 
 //! Somerton is used as default model to compute the effective thermal heat conductivity
diff --git a/dumux/porousmediumflow/2p1c/indices.hh b/dumux/porousmediumflow/2p1c/indices.hh
index 281c2a8c33c88912bd24f1685e684bfa0b88b122..74f0ddd5cb1a25dabb0397c7c9aa85193871a5c9 100644
--- a/dumux/porousmediumflow/2p1c/indices.hh
+++ b/dumux/porousmediumflow/2p1c/indices.hh
@@ -33,13 +33,12 @@ namespace Dumux {
  * \ingroup TwoPOneCModel
  * \brief The indices for the two-phase one-component model.
  *
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <class FluidSystem, int PVOffset = 0>
 class TwoPOneCIndices
 {
-    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
-
 public:
     // Phase indices
     static const int wPhaseIdx = FluidSystem::wPhaseIdx; //!< Index of the wetting phase.
diff --git a/dumux/porousmediumflow/2p1c/model.hh b/dumux/porousmediumflow/2p1c/model.hh
index bae135d9df9767c974eb389615f51ad4885518a0..aac8e7a92870988e3571c26b013f189fcb38b6e2 100644
--- a/dumux/porousmediumflow/2p1c/model.hh
+++ b/dumux/porousmediumflow/2p1c/model.hh
@@ -189,7 +189,13 @@ SET_TYPE_PROP(TwoPOneCNI, IsothermalVolumeVariables, TwoPOneCVolumeVariables<Typ
 SET_TYPE_PROP(TwoPOneCNI, IsothermalLocalResidual, ImmiscibleLocalResidual<TypeTag>);
 
 //! Set isothermal Indices.
-SET_TYPE_PROP(TwoPOneCNI, IsothermalIndices, TwoPOneCIndices<TypeTag, 0>);
+SET_PROP(TwoPOneCNI, IsothermalIndices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = TwoPOneCIndices<FluidSystem, 0>;
+};
 
 //! The isothermal vtk output fields.
 SET_TYPE_PROP(TwoPOneCNI, IsothermalVtkOutputFields, TwoPOneCVtkOutputFields<typename GET_PROP_TYPE(TypeTag, Indices)>);
diff --git a/dumux/porousmediumflow/2p2c/indices.hh b/dumux/porousmediumflow/2p2c/indices.hh
index 95f3aa74ef5b928538fd0d9cee960343c2d88b4e..f494d8521b144d97b2bdd8c7cea130286b499532 100644
--- a/dumux/porousmediumflow/2p2c/indices.hh
+++ b/dumux/porousmediumflow/2p2c/indices.hh
@@ -39,6 +39,7 @@ struct TwoPTwoCFormulation
  * \brief The indices for the isothermal two-phase two-component model.
  * \ingroup TwoPTwoCModel
  *
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
 template <class FluidSystem, int PVOffset = 0>
diff --git a/dumux/porousmediumflow/2pnc/indices.hh b/dumux/porousmediumflow/2pnc/indices.hh
index 9c8686bb2e205936c896a0dbfeb6a8829c42ad6d..a6b6d239547c6421e5e53e240b9bb15bebf147a3 100644
--- a/dumux/porousmediumflow/2pnc/indices.hh
+++ b/dumux/porousmediumflow/2pnc/indices.hh
@@ -45,14 +45,12 @@ struct TwoPNCFormulation
  * \ingroup TwoPNCModel
  * \brief The indices for the isothermal two-phase n-component model.
  *
- * \tparam TypeTag The problem Type Tag
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <class FluidSystem, int PVOffset = 0>
 class TwoPNCIndices
 {
-    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
-
 public:
     // Phase indices
     static const int wPhaseIdx = FluidSystem::wPhaseIdx;    //!< index of the wetting phase
diff --git a/dumux/porousmediumflow/2pnc/model.hh b/dumux/porousmediumflow/2pnc/model.hh
index c6c0efe7cf1a863a076b89674533d47154fed6cc..5cfff4e3457a49202ec15b8f2344f851f24e345d 100644
--- a/dumux/porousmediumflow/2pnc/model.hh
+++ b/dumux/porousmediumflow/2pnc/model.hh
@@ -129,7 +129,6 @@ public:
 
 SET_TYPE_PROP(TwoPNC, PrimaryVariableSwitch, TwoPNCPrimaryVariableSwitch<TypeTag>);         //!< The primary variable switch for the 2pnc model
 SET_TYPE_PROP(TwoPNC, VolumeVariables, TwoPNCVolumeVariables<TypeTag>);                     //!< the VolumeVariables property
-SET_TYPE_PROP(TwoPNC, Indices, TwoPNCIndices <TypeTag, /*PVOffset=*/0>);                    //!< The indices required by the isothermal 2pnc model
 SET_TYPE_PROP(TwoPNC, SpatialParams, FVSpatialParams<TypeTag>);                             //!< Use the FVSpatialParams by default
 SET_TYPE_PROP(TwoPNC, VtkOutputFields, TwoPNCVtkOutputFields<TypeTag>);                     //!< Set the vtk output fields specific to the TwoPNC model
 SET_TYPE_PROP(TwoPNC, LocalResidual, CompositionalLocalResidual<TypeTag>);                  //!< Use the compositional local residual
@@ -145,6 +144,14 @@ SET_BOOL_PROP(TwoPNC, EnableMolecularDiffusion, true);         //!< Enable molec
 SET_BOOL_PROP(TwoPNC, EnableEnergyBalance, false);             //!< This is the isothermal variant of the model
 SET_BOOL_PROP(TwoPNC, UseMoles, true);                         //!< Use mole fractions in the balance equations by default
 
+//! The indices required by the isothermal 2pnc model
+SET_PROP(TwoPNC, Indices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = TwoPNCIndices<FluidSystem, /*PVOffset=*/0>;
+};
 
 //! Use the model after Millington (1961) for the effective diffusivity
 SET_TYPE_PROP(TwoPNC, EffectiveDiffusivityModel, DiffusivityMillingtonQuirk<typename GET_PROP_TYPE(TypeTag, Scalar)>);
@@ -186,9 +193,18 @@ public:
 /////////////////////////////////////////////////
 SET_TYPE_PROP(TwoPNCNI, IsothermalVolumeVariables, TwoPNCVolumeVariables<TypeTag>);     //!< set isothermal VolumeVariables
 SET_TYPE_PROP(TwoPNCNI, IsothermalLocalResidual, CompositionalLocalResidual<TypeTag>);  //!< set isothermal LocalResidual
-SET_TYPE_PROP(TwoPNCNI, IsothermalIndices, TwoPNCIndices<TypeTag, /*PVOffset=*/0>);     //!< set isothermal Indices
 SET_TYPE_PROP(TwoPNCNI, IsothermalVtkOutputFields, TwoPNCVtkOutputFields<TypeTag>);     //!< set isothermal output fields
 
+//! set isothermal Indices
+SET_PROP(TwoPNCNI, IsothermalIndices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = TwoPNCIndices<FluidSystem, /*PVOffset=*/0>;
+};
+
+
 //! Somerton is used as default model to compute the effective thermal heat conductivity
 SET_PROP(TwoPNCNI, ThermalConductivityModel)
 {
diff --git a/dumux/porousmediumflow/3p/indices.hh b/dumux/porousmediumflow/3p/indices.hh
index 73281732b0df857f7c5c037e8efd3cbe1413fdbc..7c539c0e738719e313a2e2790bfa70f32e1bdbd0 100644
--- a/dumux/porousmediumflow/3p/indices.hh
+++ b/dumux/porousmediumflow/3p/indices.hh
@@ -33,13 +33,12 @@ namespace Dumux
  * \ingroup ThreePModel
  * \brief The common indices for the isothermal three-phase model.
  *
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <class FluidSystem, int PVOffset = 0>
 class ThreePIndices
 {
-    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
-
 public:
     // Phase indices
     static const int wPhaseIdx = FluidSystem::wPhaseIdx; //!< index of the wetting liquid phase
diff --git a/dumux/porousmediumflow/3p/model.hh b/dumux/porousmediumflow/3p/model.hh
index 734fbc6ab1e02a1b469c7fdca6390413eec8a18e..87e7da93cf7f0528e8f0ce0a9ad922c94502ae64 100644
--- a/dumux/porousmediumflow/3p/model.hh
+++ b/dumux/porousmediumflow/3p/model.hh
@@ -137,7 +137,13 @@ SET_BOOL_PROP(ThreeP, EnableEnergyBalance, false);
 SET_TYPE_PROP(ThreeP, VolumeVariables, ThreePVolumeVariables<TypeTag>);
 
 //! The indices required by the isothermal 3p model
-SET_TYPE_PROP(ThreeP, Indices, ThreePIndices<TypeTag,/*PVOffset=*/0>);
+SET_PROP(ThreeP, Indices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = ThreePIndices<FluidSystem,/*PVOffset=*/0>;
+};
 
 //! The spatial parameters to be employed.
 //! Use FVSpatialParams by default.
@@ -190,7 +196,14 @@ SET_TYPE_PROP(ThreePNI, IsothermalLocalResidual, ImmiscibleLocalResidual<TypeTag
 SET_TYPE_PROP(ThreePNI, IsothermalVtkOutputFields, ThreePVtkOutputFields<TypeTag>);
 
 //! Set isothermal Indices
-SET_TYPE_PROP(ThreePNI, IsothermalIndices, ThreePIndices<TypeTag,/*PVOffset=*/0>);
+SET_PROP(ThreePNI, IsothermalIndices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = ThreePIndices<FluidSystem,/*PVOffset=*/0>;
+};
+
 
 //! Set isothermal NumEq
 SET_INT_PROP(ThreePNI, IsothermalNumEq, 3);
diff --git a/dumux/porousmediumflow/3p3c/indices.hh b/dumux/porousmediumflow/3p3c/indices.hh
index 03ed1ca8884015548c980761dfe64c9eeac155b3..52c0273716bddaa19c808254208464ee00cdcf99 100644
--- a/dumux/porousmediumflow/3p3c/indices.hh
+++ b/dumux/porousmediumflow/3p3c/indices.hh
@@ -33,13 +33,12 @@ namespace Dumux
  * \ingroup ThreePThreeCModel
  * \brief The indices for the isothermal three-phase three-component model.
  *
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <class FluidSystem, int PVOffset = 0>
 class ThreePThreeCIndices
 {
-    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
-
 public:
     // Phase indices
     static const int wPhaseIdx = FluidSystem::wPhaseIdx; //!< index of the wetting liquid phase
diff --git a/dumux/porousmediumflow/3p3c/model.hh b/dumux/porousmediumflow/3p3c/model.hh
index 36d70f79635283453560abd0e8bd3043889407e5..5897c0957b9acdfd8865707ac41a8aa750bc2a05 100644
--- a/dumux/porousmediumflow/3p3c/model.hh
+++ b/dumux/porousmediumflow/3p3c/model.hh
@@ -183,7 +183,13 @@ SET_TYPE_PROP(ThreePThreeC, VolumeVariables, ThreePThreeCVolumeVariables<TypeTag
 SET_BOOL_PROP(ThreePThreeC, UseConstraintSolver, false);
 
 //! The indices required by the isothermal 3p3c model
-SET_TYPE_PROP(ThreePThreeC, Indices, ThreePThreeCIndices<TypeTag, /*PVOffset=*/0>);
+SET_PROP(ThreePThreeC, Indices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = ThreePThreeCIndices<FluidSystem, /*PVOffset=*/0>;
+};
 
 //! The spatial parameters to be employed.
 //! Use FVSpatialParams by default.
@@ -224,7 +230,13 @@ SET_TYPE_PROP(ThreePThreeCNI, IsothermalVolumeVariables, ThreePThreeCVolumeVaria
 SET_TYPE_PROP(ThreePThreeCNI, IsothermalLocalResidual, ThreePThreeCLocalResidual<TypeTag>);
 
 //! Set isothermal Indices
-SET_TYPE_PROP(ThreePThreeCNI, IsothermalIndices, ThreePThreeCIndices<TypeTag, /*PVOffset=*/0>);
+SET_PROP(ThreePThreeCNI, IsothermalIndices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+public:
+    using type = ThreePThreeCIndices<FluidSystem, /*PVOffset=*/0>;
+};
 
 //! Set isothermal NumEq
 SET_INT_PROP(ThreePThreeCNI, IsothermalNumEq, 3);
diff --git a/dumux/porousmediumflow/3pwateroil/indices.hh b/dumux/porousmediumflow/3pwateroil/indices.hh
index 52541a3205ba01e39c131a2ab268d5f5bb7dea36..5287096bcfd357ce4abdcd6a78df27029ba68853 100644
--- a/dumux/porousmediumflow/3pwateroil/indices.hh
+++ b/dumux/porousmediumflow/3pwateroil/indices.hh
@@ -32,14 +32,12 @@ namespace Dumux
  * \ingroup ThreePWaterOilModel
  * \brief The indices for the isothermal 3p2cni model.
  *
- * \tparam formulation The formulation, only pgSwSn
+ * \tparam FluidSystem The fluid system class
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
+template <class FluidSystem, int PVOffset = 0>
 class ThreePWaterOilIndices
 {
-    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
-
 public:
     // Phase indices
     static const int wPhaseIdx = FluidSystem::wPhaseIdx; //!< index of the wetting liquid phase
diff --git a/dumux/porousmediumflow/3pwateroil/model.hh b/dumux/porousmediumflow/3pwateroil/model.hh
index a985d08037617959562a7c49a91c41e39d0e7744..bb0519162184847b79f8dd6e44de294f438a04f7 100644
--- a/dumux/porousmediumflow/3pwateroil/model.hh
+++ b/dumux/porousmediumflow/3pwateroil/model.hh
@@ -220,9 +220,10 @@ SET_TYPE_PROP(ThreePWaterOilNI, IsothermalVtkOutputFields, ThreePWaterOilVtkOutp
 //set isothermal Indices
 SET_PROP(ThreePWaterOilNI, IsothermalIndices)
 {
-
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
 public:
-    using type = ThreePWaterOilIndices<TypeTag, /*PVOffset=*/0>;
+    using type = ThreePWaterOilIndices<FluidSystem, /*PVOffset=*/0>;
 };
 
 //set isothermal NumEq
diff --git a/dumux/porousmediumflow/mpnc/indices.hh b/dumux/porousmediumflow/mpnc/indices.hh
index 387dc0e58001fcfadcaef7a2ce4dfdd479f1b8ba..ff11f57995c49920c51acc6543dc254bc633c1f9 100644
--- a/dumux/porousmediumflow/mpnc/indices.hh
+++ b/dumux/porousmediumflow/mpnc/indices.hh
@@ -44,13 +44,16 @@ struct MpNcPressureFormulation
 /*!
  * \ingroup MPNCModel
  * \brief The primary variable and equation indices for the MpNc model.
+ *
+ * \tparam FluidSystem The fluid system class
+ * \tparam numEquationBalance Number of balance equations: all transport equations and the constraint equations
+ * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int BasePVOffset = 0>
+template <class FluidSystem, int numEquationBalance, int BasePVOffset = 0>
 class MPNCIndices
 {
-     using FluidSystem  = typename GET_PROP_TYPE(TypeTag, FluidSystem);
      enum { numPhases = FluidSystem::numPhases };
-     static const int numEqBalance = GET_PROP_VALUE(TypeTag, NumEqBalance);
+     static const int numEqBalance = numEquationBalance;
 public:
 
         // Phase indices
diff --git a/dumux/porousmediumflow/mpnc/model.hh b/dumux/porousmediumflow/mpnc/model.hh
index 7736df5cedf02bb80a49e00d181a5fae8afd3dca..a5fb39445f95dd967e0148b50483a3ee2cc05f08 100644
--- a/dumux/porousmediumflow/mpnc/model.hh
+++ b/dumux/porousmediumflow/mpnc/model.hh
@@ -136,7 +136,14 @@ NEW_TYPE_TAG(MPNCNonequil, INHERITS_FROM(MPNC, NonEquilibrium));
 //////////////////////////////////////////////////////////////////
 
 //! The indices required by the mpnc model
-SET_TYPE_PROP(MPNC, Indices, MPNCIndices<TypeTag, 0>);
+SET_PROP(MPNC, Indices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+    static constexpr int numEquationBalance = GET_PROP_VALUE(TypeTag, NumEqBalance);
+public:
+    using type = MPNCIndices<FluidSystem, numEquationBalance, /*PVOffset=*/0>;
+};
 //! Use ImplicitSpatialParams by default.
 SET_TYPE_PROP(MPNC, SpatialParams, FVSpatialParams<TypeTag>);
 
@@ -231,7 +238,15 @@ SET_INT_PROP(MPNC, NumEnergyEqSolid, 0);
 /////////////////////////////////////////////////
 SET_TYPE_PROP(MPNCNI, IsothermalVolumeVariables, MPNCVolumeVariables<TypeTag>);    //! set isothermal VolumeVariables
 SET_TYPE_PROP(MPNCNI, IsothermalLocalResidual, MPNCLocalResidual<TypeTag>); //! set isothermal LocalResidual
-SET_TYPE_PROP(MPNCNI, IsothermalIndices, MPNCIndices<TypeTag, /*PVOffset=*/0>);    //! set isothermal Indices
+//! set isothermal Indices
+SET_PROP(MPNCNI, IsothermalIndices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+    static constexpr int numEquationBalance = GET_PROP_VALUE(TypeTag, NumEqBalance);
+public:
+    using type = MPNCIndices<FluidSystem, numEquationBalance, /*PVOffset=*/0>;
+};
 //! set isothermal NumEq
 SET_INT_PROP(MPNCNI, IsothermalNumEq, GET_PROP_VALUE(TypeTag, NumEq));
 
@@ -241,7 +256,15 @@ SET_INT_PROP(MPNCNI, IsothermalNumEq, GET_PROP_VALUE(TypeTag, NumEq));
 
 SET_TYPE_PROP(MPNCNonequil, EquilibriumLocalResidual, MPNCLocalResidual<TypeTag>);
 SET_TYPE_PROP(MPNCNonequil, EquilibriumVtkOutputFields, MPNCVtkOutputFields<TypeTag>);
-SET_TYPE_PROP(MPNCNonequil, EquilibriumIndices, MPNCIndices<TypeTag, /*PVOffset=*/0>);
+
+SET_PROP(MPNCNonequil, EquilibriumIndices)
+{
+private:
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+    static constexpr int numEquationBalance = GET_PROP_VALUE(TypeTag, NumEqBalance);
+public:
+    using type = MPNCIndices<FluidSystem, numEquationBalance, /*PVOffset=*/0>;
+};
 
 //number of balance equations means all transport equations and the constraint equations, energy equations come on top of that
 SET_INT_PROP(MPNCNonequil, NumEqBalance, GET_PROP_VALUE(TypeTag, NumPhases)*GET_PROP_VALUE(TypeTag, NumComponents)+GET_PROP_VALUE(TypeTag, NumPhases));
diff --git a/dumux/porousmediumflow/nonequilibrium/indices.hh b/dumux/porousmediumflow/nonequilibrium/indices.hh
index 0bd3dca45cefb3e449d1a71a283d13bcae833864..3ff9b9d5906e5435504c07586744ba170e39ea65 100644
--- a/dumux/porousmediumflow/nonequilibrium/indices.hh
+++ b/dumux/porousmediumflow/nonequilibrium/indices.hh
@@ -33,22 +33,21 @@ namespace Dumux
  * \ingroup PorousmediumNonEquilibriumModel
  * \brief The primary variable and equation indices for the MpNc model.
  */
-template <class TypeTag, int BasePVOffset = 0>
-class NonEquilbriumIndices: public GET_PROP_TYPE(TypeTag, EquilibriumIndices)
+template <class EquilibriumIndices, class FluidSystem, int numEnergyEquationFluid, int numEnergyEquationSolid, int numEquationBalance, int BasePVOffset = 0>
+class NonEquilbriumIndices: public EquilibriumIndices
 {
 public:
-     using FluidSystem  = typename GET_PROP_TYPE(TypeTag, FluidSystem);
      enum { numPhases = FluidSystem::numPhases };
-     enum { numEnergyEqFluid = GET_PROP_VALUE(TypeTag, NumEnergyEqFluid) };
-     enum { numEnergyEqSolid = GET_PROP_VALUE(TypeTag, NumEnergyEqSolid) };
+     enum { numEnergyEqFluid = numEnergyEquationFluid };
+     enum { numEnergyEqSolid = numEnergyEquationSolid };
      /*! \todo Replacing the sum below with GET_PROP_VALUE(TypeTag, NumEq)
       *        yields a compilation error with clang, due to complex
       *        interdependencies of MPNC and NonEquilibrium type tags and
       *        Indices classes. This should be fixed.
       */
-     static const unsigned int numEq = GET_PROP_VALUE(TypeTag, NumEqBalance)
-                                     + GET_PROP_VALUE(TypeTag, NumEnergyEqFluid)
-                                     + GET_PROP_VALUE(TypeTag, NumEnergyEqSolid);
+     static const unsigned int numEq = numEquationBalance
+                                     + numEnergyEqFluid
+                                     + numEnergyEqSolid;
 
     /*!
      * \brief Index for the temperature of the wetting phase in a vector of primary
diff --git a/dumux/porousmediumflow/nonequilibrium/model.hh b/dumux/porousmediumflow/nonequilibrium/model.hh
index 21c3043e86ba0e66367deae9b99f3890dfee47b8..ce3385be411991ed1a2299a3e0238f7add8536e7 100644
--- a/dumux/porousmediumflow/nonequilibrium/model.hh
+++ b/dumux/porousmediumflow/nonequilibrium/model.hh
@@ -75,7 +75,17 @@ SET_TYPE_PROP(NonEquilibrium, HeatConductionType , FouriersLawNonEquilibrium<Typ
 SET_INT_PROP(NonEquilibrium, NumEq, GET_PROP_VALUE(TypeTag, NumEqBalance) +  GET_PROP_VALUE(TypeTag, NumEnergyEqFluid) + GET_PROP_VALUE(TypeTag, NumEnergyEqSolid));
 
 //! indices for non-isothermal models
-SET_TYPE_PROP(NonEquilibrium, Indices, NonEquilbriumIndices<TypeTag, 0>);
+SET_PROP(NonEquilibrium, Indices)
+{
+private:
+    using EquilibriumIndices = typename GET_PROP_TYPE(TypeTag, EquilibriumIndices);
+    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
+    static constexpr int numEnergyEquationFluid = GET_PROP_VALUE(TypeTag, NumEnergyEqFluid);
+    static constexpr int numEnergyEquationSolid = GET_PROP_VALUE(TypeTag, NumEnergyEqSolid);
+    static constexpr int numEquationBalance = GET_PROP_VALUE(TypeTag, NumEqBalance);
+public:
+    using type = NonEquilbriumIndices<EquilibriumIndices, FluidSystem, numEnergyEquationFluid, numEnergyEquationSolid, numEquationBalance, 0>;
+};
 
 SET_PROP(NonEquilibrium, FluidState)
 {
diff --git a/dumux/porousmediumflow/nonisothermal/indices.hh b/dumux/porousmediumflow/nonisothermal/indices.hh
index a53dfef4702382c79482d7ccf5c427e3cdc373c9..febb62ffdb89f8182910ff83d0cb34bd2f355aaf 100644
--- a/dumux/porousmediumflow/nonisothermal/indices.hh
+++ b/dumux/porousmediumflow/nonisothermal/indices.hh
@@ -37,11 +37,11 @@ namespace Dumux
  * \tparam formulation The formulation, either pwsn or pnsw.
  * \tparam PVOffset The first index in a primary variable vector.
  */
-template <class TypeTag, int PVOffset = 0>
-class EnergyIndices : public GET_PROP_TYPE(TypeTag, IsothermalIndices)
+template <class IsothermalIndices, int numEquation, int PVOffset = 0>
+class EnergyIndices : public IsothermalIndices
 {
 public:
-    static const int numEq = GET_PROP_VALUE(TypeTag, NumEq);
+    static const int numEq = numEquation;
     static const int temperatureIdx = PVOffset + numEq -1; //!< The index for temperature in primary variable vectors.
     static const int energyEqIdx = PVOffset + numEq -1; //!< The index for energy in equation vectors.
 
diff --git a/dumux/porousmediumflow/nonisothermal/model.hh b/dumux/porousmediumflow/nonisothermal/model.hh
index f0ecd5208aa2c3b58929578d73836a76afa24970..a7aec809d986d3e7ef763cb7c8f62fe4891d7ab7 100644
--- a/dumux/porousmediumflow/nonisothermal/model.hh
+++ b/dumux/porousmediumflow/nonisothermal/model.hh
@@ -68,7 +68,14 @@ SET_BOOL_PROP(NonIsothermal, EnableEnergyBalance, true);
 SET_INT_PROP(NonIsothermal, NumEq, GET_PROP_VALUE(TypeTag, IsothermalNumEq) + 1);
 
 //! indices for non-isothermal models
-SET_TYPE_PROP(NonIsothermal, Indices, EnergyIndices<TypeTag, 0>);
+SET_PROP(NonIsothermal, Indices)
+{
+private:
+    using IsothermalIndices = typename GET_PROP_TYPE(TypeTag, IsothermalIndices);
+    static constexpr int numEq = GET_PROP_VALUE(TypeTag, NumEq);
+public:
+    using type = EnergyIndices<IsothermalIndices, numEq, 0>;
+};
 
 //! indices for non-isothermal models
 SET_TYPE_PROP(NonIsothermal, VtkOutputFields, EnergyVtkOutputFields<TypeTag>);
diff --git a/dumux/porousmediumflow/richardsnc/indices.hh b/dumux/porousmediumflow/richardsnc/indices.hh
index fd9f255d74c52490c8f3db3fd60aaf46b7c63772..efaac711e0c945de167191dbd345a661784d67c2 100644
--- a/dumux/porousmediumflow/richardsnc/indices.hh
+++ b/dumux/porousmediumflow/richardsnc/indices.hh
@@ -34,7 +34,7 @@ namespace Dumux
  * \ingroup RichardsNCModel
  * \brief The indices for the isothermal Richards, n-component model.
  */
-template <class TypeTag, int PVOffset = 0>
+template <int PVOffset = 0>
 struct RichardsNCIndices
 {
 
diff --git a/dumux/porousmediumflow/richardsnc/model.hh b/dumux/porousmediumflow/richardsnc/model.hh
index 2f4035de370227e04a8d7bdfb4f89402e761c4de..ff747c5a41cb1c8baf7f73e89db3ac5369cc3c64 100644
--- a/dumux/porousmediumflow/richardsnc/model.hh
+++ b/dumux/porousmediumflow/richardsnc/model.hh
@@ -162,7 +162,7 @@ SET_PROP(RichardsNC, FluidState)
 SET_TYPE_PROP(RichardsNC, VtkOutputFields, RichardsNCVtkOutputFields<TypeTag>);           //!< Set the vtk output fields specific to the twop model
 
 //! Set the indices used
-SET_TYPE_PROP(RichardsNC, Indices, RichardsNCIndices<TypeTag>);
+SET_TYPE_PROP(RichardsNC, Indices, RichardsNCIndices<>);
 //! The spatial parameters to be employed.
 //! Use FVSpatialParamsOneP by default.
 SET_TYPE_PROP(RichardsNC, SpatialParams, FVSpatialParamsOneP<TypeTag>);
@@ -197,7 +197,7 @@ SET_TYPE_PROP(RichardsNCNI, IsothermalVolumeVariables, RichardsNCVolumeVariables
 SET_TYPE_PROP(RichardsNCNI, IsothermalLocalResidual, CompositionalLocalResidual<TypeTag>);
 
 //set isothermal Indices
-SET_TYPE_PROP(RichardsNCNI, IsothermalIndices, RichardsNCIndices<TypeTag>);
+SET_TYPE_PROP(RichardsNCNI, IsothermalIndices, RichardsNCIndices<>);
 
 //set isothermal NumEq
 SET_PROP(RichardsNCNI, IsothermalNumEq)
diff --git a/dumux/porousmediumflow/tracer/indices.hh b/dumux/porousmediumflow/tracer/indices.hh
index dbfe144747f41913e7c2910506c2356dcd880a0f..71e67488a2e92ffcc94e52cd0f307a2ad5af5b3b 100644
--- a/dumux/porousmediumflow/tracer/indices.hh
+++ b/dumux/porousmediumflow/tracer/indices.hh
@@ -33,7 +33,7 @@ namespace Dumux
  * \ingroup TracerModel
  * \brief Defines the primary variable and equation indices used by the isothermal tracer model.
  */
-template <class TypeTag, int PVOffset = 0>
+template <int PVOffset = 0>
 struct TracerIndices
 {
     /*!
diff --git a/dumux/porousmediumflow/tracer/model.hh b/dumux/porousmediumflow/tracer/model.hh
index 019fd7d17f52f4021da1fe13673b65dadb5c3faa..38749aa40f16a6c182103518f176981c4802b60f 100644
--- a/dumux/porousmediumflow/tracer/model.hh
+++ b/dumux/porousmediumflow/tracer/model.hh
@@ -107,7 +107,7 @@ SET_TYPE_PROP(Tracer, VolumeVariables, TracerVolumeVariables<TypeTag>);
 SET_TYPE_PROP(Tracer, AdvectionType, StationaryVelocityField<typename GET_PROP_TYPE(TypeTag, Scalar)>);
 
 //! Set the indices used by the tracer model
-SET_TYPE_PROP(Tracer, Indices, TracerIndices<TypeTag>);
+SET_TYPE_PROP(Tracer, Indices, TracerIndices<>);
 
 //! Use FVSpatialParamsOneP by default.
 SET_TYPE_PROP(Tracer, SpatialParams, FVSpatialParamsOneP<TypeTag>);