From b5925aaab29977c3f4703f154e1b0533b320c573 Mon Sep 17 00:00:00 2001
From: Timo Koch <timo.koch@iws.uni-stuttgart.de>
Date: Wed, 18 Jul 2018 00:05:17 +0200
Subject: [PATCH] [1pnc] Remove fluid system phase index

We can use one-phase adapter for reusing m-phase fluid systems now.
Reverts 7bb9a857 which introduced a bug for the diffusion coefficients
---
 dumux/porousmediumflow/1pnc/indices.hh        |   5 +-
 dumux/porousmediumflow/1pnc/model.hh          |  18 +--
 .../porousmediumflow/1pnc/volumevariables.hh  | 109 ++++++------------
 .../porousmediumflow/1pnc/vtkoutputfields.hh  |  17 ++-
 .../1pnc/implicit/1p2cniconductionproblem.hh  |  14 ++-
 .../1pnc/implicit/1p2cniconvectionproblem.hh  |  26 +++--
 .../1pnc/implicit/1p2ctestproblem.hh          |  36 +++---
 .../1pnc/implicit/test_1p2c_fv.cc             |   2 +-
 .../implicit/test_1p2cni_conduction_fv.cc     |   2 +-
 .../implicit/test_1p2cni_convection_fv.cc     |   2 +-
 10 files changed, 109 insertions(+), 122 deletions(-)

diff --git a/dumux/porousmediumflow/1pnc/indices.hh b/dumux/porousmediumflow/1pnc/indices.hh
index 8cc3114029..385e3c63ae 100644
--- a/dumux/porousmediumflow/1pnc/indices.hh
+++ b/dumux/porousmediumflow/1pnc/indices.hh
@@ -34,15 +34,12 @@ namespace Dumux {
  *
  * \tparam phaseIdx The index of the fluid phase in the fluid system
  */
-template<int phaseIdx>
 struct OnePNCIndices
 {
     //! Reference index for mass conservation equation.
     static constexpr int conti0EqIdx = 0;
     //! Index for wetting/non-wetting phase pressure (depending on formulation) in a solution vector
-    static constexpr int pressureIdx = phaseIdx;
-    //! The index of the fluid phase in the fluid system
-    static constexpr int fluidSystemPhaseIdx = phaseIdx;
+    static constexpr int pressureIdx = 0;
 };
 
 } // end namespace Dumux
diff --git a/dumux/porousmediumflow/1pnc/model.hh b/dumux/porousmediumflow/1pnc/model.hh
index 38ecb14625..80471fa797 100644
--- a/dumux/porousmediumflow/1pnc/model.hh
+++ b/dumux/porousmediumflow/1pnc/model.hh
@@ -82,12 +82,11 @@ namespace Dumux {
  *        consider a single-phase with multiple components.
  *
  * \tparam nComp the number of components to be considered.
- * \tparam fluidSystemPhaseIdx The index of the fluid phase in the fluid system
  */
-template<int nComp, int fluidSystemPhaseIdx, bool useM, int repCompEqIdx = nComp>
+template<int nComp, bool useM, int repCompEqIdx = nComp>
 struct OnePNCModelTraits
 {
-    using Indices = OnePNCIndices<fluidSystemPhaseIdx>;
+    using Indices = OnePNCIndices;
 
     static constexpr int numEq() { return nComp; }
     static constexpr int numPhases() { return 1; }
@@ -145,7 +144,7 @@ SET_PROP(OnePNC, ModelTraits)
 private:
     using FluidSystem = typename GET_PROP_TYPE(TypeTag, PTAG(FluidSystem));
 public:
-    using type = OnePNCModelTraits<FluidSystem::numComponents, GET_PROP_VALUE(TypeTag, PhaseIdx), GET_PROP_VALUE(TypeTag, UseMoles), GET_PROP_VALUE(TypeTag, ReplaceCompEqIdx)>;
+    using type = OnePNCModelTraits<FluidSystem::numComponents, GET_PROP_VALUE(TypeTag, UseMoles), GET_PROP_VALUE(TypeTag, ReplaceCompEqIdx)>;
 };
 
 
@@ -155,7 +154,8 @@ public:
  *        appropriately for the model ((non-)isothermal, equilibrium, ...).
  *        This can be done in the problem.
  */
-SET_PROP(OnePNC, FluidState){
+SET_PROP(OnePNC, FluidState)
+{
 private:
     using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
     using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
@@ -171,8 +171,8 @@ SET_TYPE_PROP(OnePNC, EffectiveDiffusivityModel,
 //! Use mole fractions in the balance equations by default
 SET_BOOL_PROP(OnePNC, UseMoles, true);
 
-SET_INT_PROP(OnePNC, PhaseIdx, 0); //!< The default phase index
-SET_TYPE_PROP(OnePNC, LocalResidual, CompositionalLocalResidual<TypeTag>);        //!< The local residual function
+//! The local residual function
+SET_TYPE_PROP(OnePNC, LocalResidual, CompositionalLocalResidual<TypeTag>);
 
 //! Set the volume variables property
 SET_PROP(OnePNC, VolumeVariables)
@@ -187,6 +187,8 @@ private:
     using PT = typename GET_PROP_TYPE(TypeTag, SpatialParams)::PermeabilityType;
     static_assert(FSY::numComponents == MT::numComponents(), "Number of components mismatch between model and fluid system");
     static_assert(FST::numComponents == MT::numComponents(), "Number of components mismatch between model and fluid state");
+    static_assert(FSY::numPhases == MT::numPhases(), "Number of phases mismatch between model and fluid system");
+    static_assert(FST::numPhases == MT::numPhases(), "Number of phases mismatch between model and fluid state");
 
     using Traits = OnePNCVolumeVariablesTraits<PV, FSY, FST, SSY, SST, PT, MT>;
 public:
@@ -213,7 +215,7 @@ SET_PROP(OnePNCNI, ModelTraits)
 {
 private:
     using FluidSystem = typename GET_PROP_TYPE(TypeTag, PTAG(FluidSystem));
-    using IsothermalTraits = OnePNCModelTraits<FluidSystem::numComponents, GET_PROP_VALUE(TypeTag, PhaseIdx), GET_PROP_VALUE(TypeTag, UseMoles), GET_PROP_VALUE(TypeTag, ReplaceCompEqIdx)>;
+    using IsothermalTraits = OnePNCModelTraits<FluidSystem::numComponents, GET_PROP_VALUE(TypeTag, UseMoles), GET_PROP_VALUE(TypeTag, ReplaceCompEqIdx)>;
 public:
     using type = PorousMediumFlowNIModelTraits<IsothermalTraits>;
 };
diff --git a/dumux/porousmediumflow/1pnc/volumevariables.hh b/dumux/porousmediumflow/1pnc/volumevariables.hh
index b5b9f1eca8..ca23fde92f 100644
--- a/dumux/porousmediumflow/1pnc/volumevariables.hh
+++ b/dumux/porousmediumflow/1pnc/volumevariables.hh
@@ -56,13 +56,8 @@ class OnePNCVolumeVariables
 
     enum
     {
-        fluidSystemPhaseIdx = Idx::fluidSystemPhaseIdx,
-
         // pressure primary variable index
-        pressureIdx = Idx::pressureIdx,
-
-        // main component index
-        mainCompMoleOrMassFracIdx = fluidSystemPhaseIdx
+        pressureIdx = Idx::pressureIdx
     };
 
 public:
@@ -104,24 +99,11 @@ public:
         // Could be avoided if diffusion coefficients also
         // became part of the fluid state.
         typename FluidSystem::ParameterCache paramCache;
-        paramCache.updatePhase(fluidState_, fluidSystemPhaseIdx);
+        paramCache.updatePhase(fluidState_, 0);
 
-        int compIIdx = mainCompMoleOrMassFracIdx;
-        for (unsigned int compJIdx = 0; compJIdx < numFluidComps; ++compJIdx)
-        {
-            diffCoeff_[compJIdx] = 0.0;
-            if(compIIdx != compJIdx)
-                diffCoeff_[compJIdx] = FluidSystem::binaryDiffusionCoefficient(fluidState_,
-                                                                               paramCache,
-                                                                               fluidSystemPhaseIdx,
-                                                                               compIIdx,
-                                                                               compJIdx);
-        }
-        // The diffusion coefficients are swapped, the general procedure for the flux
-        // calculation always goes form phaseIdx = 0 to numPhaseIdx.
-        // Swapping the coefficients ensures that the diffusion coefficient for compIdx == 0
-        // is always defined.
-        std::swap(diffCoeff_[0],diffCoeff_[mainCompMoleOrMassFracIdx]);
+        diffCoeff_[0] = 0.0; // the main component with itself doesn't have a binary diffusion coefficient
+        for (unsigned int compJIdx = 1; compJIdx < numFluidComps; ++compJIdx)
+            diffCoeff_[compJIdx] = FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, 0, 0, compJIdx);
     }
 
     /*!
@@ -144,45 +126,34 @@ public:
 
     {
         EnergyVolVars::updateTemperature(elemSol, problem, element, scv, fluidState, solidState);
-        fluidState.setSaturation(fluidSystemPhaseIdx, 1.);
+        fluidState.setSaturation(0, 1.0);
 
         const auto& priVars = elemSol[scv.localDofIndex()];
-        fluidState.setPressure(fluidSystemPhaseIdx, priVars[pressureIdx]);
-
-        // calculate the phase composition
-        Dune::FieldVector<Scalar, numFluidComps> moleFrac;
-
-        Scalar sumMoleFracNotMainComp = 0;
-        for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
-        {
-            if (compIdx != mainCompMoleOrMassFracIdx)
-            {
-                moleFrac[compIdx] = priVars[compIdx];
-                sumMoleFracNotMainComp += moleFrac[compIdx];
-            }
-        }
-        moleFrac[mainCompMoleOrMassFracIdx] = 1- sumMoleFracNotMainComp;
+        fluidState.setPressure(0, priVars[pressureIdx]);
 
         // Set fluid state mole fractions
-        for (int compIdx = 0; compIdx < numFluidComps; ++compIdx)
+        Scalar sumMoleFracNotMainComp = 0;
+        for (int compIdx = 1; compIdx < numFluidComps; ++compIdx)
         {
-            fluidState.setMoleFraction(fluidSystemPhaseIdx, compIdx, moleFrac[compIdx]);
+            fluidState.setMoleFraction(0, compIdx, priVars[compIdx]);
+            sumMoleFracNotMainComp += priVars[compIdx];
         }
+        fluidState.setMoleFraction(0, 0, 1.0 - sumMoleFracNotMainComp);
 
         typename FluidSystem::ParameterCache paramCache;
         paramCache.updateAll(fluidState);
 
-        Scalar rho = FluidSystem::density(fluidState, paramCache, fluidSystemPhaseIdx);
-        Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, fluidSystemPhaseIdx);
-        Scalar mu = FluidSystem::viscosity(fluidState, paramCache, fluidSystemPhaseIdx);
+        Scalar rho = FluidSystem::density(fluidState, paramCache, 0);
+        Scalar rhoMolar = FluidSystem::molarDensity(fluidState, paramCache, 0);
+        Scalar mu = FluidSystem::viscosity(fluidState, paramCache, 0);
 
-        fluidState.setDensity(fluidSystemPhaseIdx, rho);
-        fluidState.setMolarDensity(fluidSystemPhaseIdx, rhoMolar);
-        fluidState.setViscosity(fluidSystemPhaseIdx, mu);
+        fluidState.setDensity(0, rho);
+        fluidState.setMolarDensity(0, rhoMolar);
+        fluidState.setViscosity(0, mu);
 
         // compute and set the enthalpy
-        Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, fluidSystemPhaseIdx);
-        fluidState.setEnthalpy(fluidSystemPhaseIdx, h);
+        Scalar h = EnergyVolVars::enthalpy(fluidState, paramCache, 0);
+        fluidState.setEnthalpy(0, h);
     }
 
     /*!
@@ -204,9 +175,9 @@ public:
      * \note the phase index passed to this function is for compatibility reasons
      *       with multiphasic models.
      */
-    Scalar density(int phaseIdx = fluidSystemPhaseIdx) const
+    Scalar density(int phaseIdx = 0) const
     {
-        return fluidState_.density(fluidSystemPhaseIdx);
+        return fluidState_.density(0);
     }
 
     /*!
@@ -215,9 +186,9 @@ public:
      * \note the phase index passed to this function is for compatibility reasons
      *       with multiphasic models.
      */
-    Scalar molarDensity(int phaseIdx = fluidSystemPhaseIdx) const
+    Scalar molarDensity(int phaseIdx = 0) const
     {
-        return fluidState_.molarDensity(fluidSystemPhaseIdx);
+        return fluidState_.molarDensity(0);
     }
 
     /*!
@@ -226,7 +197,7 @@ public:
      * This method is here for compatibility reasons with other models. The saturation
      * is always 1.0 in a one-phasic context.
      */
-    Scalar saturation(int phaseIdx = fluidSystemPhaseIdx) const
+    Scalar saturation(int phaseIdx = 0) const
     { return 1.0; }
 
      /*!
@@ -242,7 +213,7 @@ public:
      {
          // make sure this is only called with admissible indices
          assert(compIdx < numFluidComps);
-         return fluidState_.moleFraction(fluidSystemPhaseIdx, compIdx);
+         return fluidState_.moleFraction(0, compIdx);
      }
 
      /*!
@@ -258,7 +229,7 @@ public:
      {
          // make sure this is only called with admissible indices
          assert(compIdx < numFluidComps);
-         return fluidState_.massFraction(fluidSystemPhaseIdx, compIdx);
+         return fluidState_.massFraction(0, compIdx);
      }
 
     /*!
@@ -270,9 +241,9 @@ public:
      * \note the phase index passed to this function is for compatibility reasons
      *       with multiphasic models.
      */
-    Scalar pressure(int phaseIdx = fluidSystemPhaseIdx) const
+    Scalar pressure(int phaseIdx = 0) const
     {
-        return fluidState_.pressure(fluidSystemPhaseIdx);
+        return fluidState_.pressure(0);
     }
 
     /*!
@@ -294,9 +265,9 @@ public:
      * \note the phase index passed to this function is for compatibility reasons
      *       with multiphasic models.
      */
-    Scalar mobility(int phaseIdx = fluidSystemPhaseIdx) const
+    Scalar mobility(int phaseIdx = 0) const
     {
-        return 1.0/fluidState_.viscosity(fluidSystemPhaseIdx);
+        return 1.0/fluidState_.viscosity(0);
     }
 
     /*!
@@ -306,9 +277,9 @@ public:
      * \note the phase index passed to this function is for compatibility reasons
      *       with multiphasic models.
      */
-    Scalar viscosity(int phaseIdx = fluidSystemPhaseIdx) const
+    Scalar viscosity(int phaseIdx = 0) const
     {
-        return fluidState_.viscosity(fluidSystemPhaseIdx);
+        return fluidState_.viscosity(0);
     }
 
     /*!
@@ -319,9 +290,6 @@ public:
 
     /*!
      * \brief Return the binary diffusion coefficient \f$\mathrm{[m^2/s]}\f$ in the fluid.
-     *
-     * \note For fluidSystemPhaseIdx > 0, the diffusion coefficients
-     *        diffCoeff_[0] and  diffCoeff_[mainCompMoleOrMassFracIdx] are swapped
      */
     Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
     {
@@ -337,7 +305,7 @@ public:
     Scalar molarity(int compIdx) const // [moles/m^3]
     {
         assert(compIdx < numFluidComps);
-        return fluidState_.molarity(fluidSystemPhaseIdx, compIdx);
+        return fluidState_.molarity(0, compIdx);
     }
 
      /*!
@@ -348,7 +316,7 @@ public:
      Scalar massFraction(int compIdx) const
      {
          assert(compIdx < numFluidComps);
-         return this->fluidState_.massFraction(fluidSystemPhaseIdx, compIdx);
+         return this->fluidState_.massFraction(0, compIdx);
      }
 
     /*!
@@ -362,10 +330,9 @@ protected:
     SolidState solidState_;
 
 private:
-    Scalar porosity_;        //!< Effective porosity within the control volume
-    PermeabilityType permeability_;
-    Scalar density_;
-    Dune::FieldVector<Scalar, numFluidComps> diffCoeff_;
+    Scalar porosity_; //!< Effective porosity within the control volume
+    PermeabilityType permeability_; //!< Effective permeability within the control volume
+    Dune::FieldVector<Scalar, numFluidComps> diffCoeff_; //!< Binary diffusion coefficients
 };
 
 } // end namespace Dumux
diff --git a/dumux/porousmediumflow/1pnc/vtkoutputfields.hh b/dumux/porousmediumflow/1pnc/vtkoutputfields.hh
index c2b3b3b1a3..7ac20b20b5 100644
--- a/dumux/porousmediumflow/1pnc/vtkoutputfields.hh
+++ b/dumux/porousmediumflow/1pnc/vtkoutputfields.hh
@@ -40,20 +40,19 @@ public:
     {
         using VolumeVariables = typename VtkOutputModule::VolumeVariables;
         using FluidSystem = typename VolumeVariables::FluidSystem;
-        using Indices = typename VolumeVariables::Indices;
 
-        vtk.addVolumeVariable([](const auto& volVars){ return volVars.pressure(Indices::fluidSystemPhaseIdx); }, "p");
-        vtk.addVolumeVariable([](const auto& volVars){ return volVars.density(Indices::fluidSystemPhaseIdx); }, "rho");
-        vtk.addVolumeVariable([](const auto& volVars){ return volVars.viscosity(Indices::fluidSystemPhaseIdx); }, "mu");
-        vtk.addVolumeVariable([](const auto& volVars){ return volVars.pressure(Indices::fluidSystemPhaseIdx) - 1e5; }, "delp");
+        vtk.addVolumeVariable([](const auto& volVars){ return volVars.pressure(0); }, "p");
+        vtk.addVolumeVariable([](const auto& volVars){ return volVars.density(0); }, "rho");
+        vtk.addVolumeVariable([](const auto& volVars){ return volVars.viscosity(0); }, "mu");
+        vtk.addVolumeVariable([](const auto& volVars){ return volVars.pressure(0) - 1e5; }, "delp");
 
         for (int i = 0; i < VolumeVariables::numComponents(); ++i)
-           vtk.addVolumeVariable([i](const auto& volVars){ return volVars.moleFraction(Indices::fluidSystemPhaseIdx, i); },
-                                     "x^" + std::string(FluidSystem::componentName(i)) + "_" + std::string(FluidSystem::phaseName(Indices::fluidSystemPhaseIdx)));
+           vtk.addVolumeVariable([i](const auto& volVars){ return volVars.moleFraction(0, i); },
+                                     "x^" + std::string(FluidSystem::componentName(i)) + "_" + std::string(FluidSystem::phaseName(0)));
 
         for (int i = 0; i < VolumeVariables::numComponents(); ++i)
-           vtk.addVolumeVariable([i](const auto& volVars){ return volVars.massFraction(Indices::fluidSystemPhaseIdx, i); },
-                                     "X^" + std::string(FluidSystem::componentName(i))+ "_" + std::string(FluidSystem::phaseName(Indices::fluidSystemPhaseIdx)));
+           vtk.addVolumeVariable([i](const auto& volVars){ return volVars.massFraction(0, i); },
+                                     "X^" + std::string(FluidSystem::componentName(i))+ "_" + std::string(FluidSystem::phaseName(0)));
     }
 };
 
diff --git a/test/porousmediumflow/1pnc/implicit/1p2cniconductionproblem.hh b/test/porousmediumflow/1pnc/implicit/1p2cniconductionproblem.hh
index aebaf03ff2..e720bbccf4 100644
--- a/test/porousmediumflow/1pnc/implicit/1p2cniconductionproblem.hh
+++ b/test/porousmediumflow/1pnc/implicit/1p2cniconductionproblem.hh
@@ -37,6 +37,7 @@
 #include <dumux/porousmediumflow/1pnc/model.hh>
 #include <dumux/porousmediumflow/problem.hh>
 
+#include <dumux/material/fluidsystems/1padapter.hh>
 #include <dumux/material/fluidsystems/h2on2.hh>
 #include <dumux/material/components/h2o.hh>
 #include "1pnctestspatialparams.hh"
@@ -65,9 +66,12 @@ SET_TYPE_PROP(OnePTwoCNIConductionTypeTag, Grid, Dune::YaspGrid<2>);
 SET_TYPE_PROP(OnePTwoCNIConductionTypeTag, Problem, OnePTwoCNIConductionProblem<TypeTag>);
 
 // Set fluid configuration
-SET_TYPE_PROP(OnePTwoCNIConductionTypeTag,
-              FluidSystem,
-              FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar)>);
+SET_PROP(OnePTwoCNIConductionTypeTag, FluidSystem)
+{
+    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
+    using H2ON2 = FluidSystems::H2ON2<Scalar, FluidSystems::H2ON2DefaultPolicy</*simplified=*/true>>;
+    using type = FluidSystems::OnePAdapter<H2ON2, H2ON2::liquidPhaseIdx>;
+};
 
 // Set the spatial parameters
 SET_TYPE_PROP(OnePTwoCNIConductionTypeTag, SpatialParams, OnePNCTestSpatialParams<TypeTag>);
@@ -126,6 +130,8 @@ class OnePTwoCNIConductionProblem : public PorousMediumFlowProblem<TypeTag>
         // indices of the primary variables
         pressureIdx = Indices::pressureIdx,
         temperatureIdx = Indices::temperatureIdx,
+
+        N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx)
     };
 
     //! property that defines whether mole or mass fractions are used
@@ -304,7 +310,7 @@ private:
     {
         PrimaryVariables priVars;
         priVars[pressureIdx] = 1e5; // initial condition for the pressure
-        priVars[FluidSystem::N2Idx] = 1e-5;  // initial condition for the N2 molefraction
+        priVars[N2Idx] = 1e-5;  // initial condition for the N2 molefraction
         priVars[temperatureIdx] = 290.;
         return priVars;
     }
diff --git a/test/porousmediumflow/1pnc/implicit/1p2cniconvectionproblem.hh b/test/porousmediumflow/1pnc/implicit/1p2cniconvectionproblem.hh
index a3a76a8afc..52b75eaf85 100644
--- a/test/porousmediumflow/1pnc/implicit/1p2cniconvectionproblem.hh
+++ b/test/porousmediumflow/1pnc/implicit/1p2cniconvectionproblem.hh
@@ -37,6 +37,7 @@
 #include <dumux/porousmediumflow/1pnc/model.hh>
 #include <dumux/porousmediumflow/problem.hh>
 
+#include <dumux/material/fluidsystems/1padapter.hh>
 #include <dumux/material/fluidsystems/h2on2.hh>
 #include <dumux/material/components/h2o.hh>
 #include "1pnctestspatialparams.hh"
@@ -65,9 +66,12 @@ SET_TYPE_PROP(OnePTwoCNIConvectionTypeTag, Grid, Dune::YaspGrid<2>);
 SET_TYPE_PROP(OnePTwoCNIConvectionTypeTag, Problem, OnePTwoCNIConvectionProblem<TypeTag>);
 
 // Set fluid configuration
-SET_TYPE_PROP(OnePTwoCNIConvectionTypeTag,
-              FluidSystem,
-              FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar)>);
+SET_PROP(OnePTwoCNIConvectionTypeTag, FluidSystem)
+{
+    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
+    using H2ON2 = FluidSystems::H2ON2<Scalar, FluidSystems::H2ON2DefaultPolicy</*simplified=*/true>>;
+    using type = FluidSystems::OnePAdapter<H2ON2, H2ON2::liquidPhaseIdx>;
+};
 
 // Set the spatial parameters
 SET_TYPE_PROP(OnePTwoCNIConvectionTypeTag, SpatialParams, OnePNCTestSpatialParams<TypeTag>);
@@ -129,9 +133,13 @@ class OnePTwoCNIConvectionProblem : public PorousMediumFlowProblem<TypeTag>
         pressureIdx = Indices::pressureIdx,
         temperatureIdx = Indices::temperatureIdx,
 
-        // equation indices
-        contiH2OIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx,
-        contiN2Idx = Indices::conti0EqIdx + FluidSystem::N2Idx,
+        // component indices
+        H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
+        N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
+
+        // indices of the equations
+        contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
+        contiN2EqIdx = Indices::conti0EqIdx + N2Idx,
         energyEqIdx = Indices::energyEqIdx
     };
 
@@ -295,8 +303,8 @@ public:
 
         if(globalPos[0] < eps_)
         {
-             flux[contiH2OIdx] = -darcyVelocity_*elemVolVars[scv].molarDensity();
-             flux[contiN2Idx] = -darcyVelocity_*elemVolVars[scv].molarDensity()*elemVolVars[scv].moleFraction(0, FluidSystem::N2Idx);
+             flux[contiH2OEqIdx] = -darcyVelocity_*elemVolVars[scv].molarDensity();
+             flux[contiN2EqIdx] = -darcyVelocity_*elemVolVars[scv].molarDensity()*elemVolVars[scv].moleFraction(0, N2Idx);
              flux[energyEqIdx] = -darcyVelocity_
                                  *elemVolVars[scv].density()
                                  *IapwsH2O::liquidEnthalpy(temperatureHigh_, elemVolVars[scv].pressure());
@@ -345,7 +353,7 @@ private:
     {
         PrimaryVariables priVars;
         priVars[pressureIdx] = pressureLow_; // initial condition for the pressure
-        priVars[FluidSystem::N2Idx] = 1e-10;  // initial condition for the N2 molefraction
+        priVars[N2Idx] = 1e-10;  // initial condition for the N2 molefraction
         priVars[temperatureIdx] = temperatureLow_;
         return priVars;
     }
diff --git a/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh b/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh
index cd01f340f7..9f76cef50c 100644
--- a/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh
+++ b/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh
@@ -39,16 +39,17 @@
 #include <dumux/porousmediumflow/problem.hh>
 
 #include <dumux/material/fluidsystems/h2on2.hh>
+#include <dumux/material/fluidsystems/1padapter.hh>
+
 #include "1pnctestspatialparams.hh"
 
-namespace Dumux
-{
+namespace Dumux {
 
 template <class TypeTag>
 class OnePTwoCTestProblem;
 
-namespace Properties
-{
+namespace Properties {
+
 NEW_TYPE_TAG(OnePTwoCTestTypeTag, INHERITS_FROM(OnePNC));
 NEW_TYPE_TAG(OnePTwoCTestBoxTypeTag, INHERITS_FROM(BoxModel, OnePTwoCTestTypeTag));
 NEW_TYPE_TAG(OnePTwoCTestCCTpfaTypeTag, INHERITS_FROM(CCTpfaModel, OnePTwoCTestTypeTag));
@@ -65,9 +66,12 @@ SET_TYPE_PROP(OnePTwoCTestTypeTag, Grid, Dune::YaspGrid<2>);
 SET_TYPE_PROP(OnePTwoCTestTypeTag, Problem, OnePTwoCTestProblem<TypeTag>);
 
 // Set fluid configuration
-SET_TYPE_PROP(OnePTwoCTestTypeTag,
-              FluidSystem,
-              FluidSystems::H2ON2<typename GET_PROP_TYPE(TypeTag, Scalar), FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/true>>);
+SET_PROP(OnePTwoCTestTypeTag, FluidSystem)
+{
+    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
+    using H2ON2 = FluidSystems::H2ON2<Scalar, FluidSystems::H2ON2DefaultPolicy</*simplified=*/true>>;
+    using type = FluidSystems::OnePAdapter<H2ON2, H2ON2::liquidPhaseIdx>;
+};
 
 // Set the spatial parameters
 SET_TYPE_PROP(OnePTwoCTestTypeTag, SpatialParams, OnePNCTestSpatialParams<TypeTag>);
@@ -125,9 +129,13 @@ class OnePTwoCTestProblem : public PorousMediumFlowProblem<TypeTag>
         // indices of the primary variables
         pressureIdx = Indices::pressureIdx,
 
+        // component indices
+        H2OIdx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::H2OIdx),
+        N2Idx = FluidSystem::compIdx(FluidSystem::MultiPhaseFluidSystem::N2Idx),
+
         // indices of the equations
-        contiH2OEqIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx,
-        contiN2EqIdx = Indices::conti0EqIdx + FluidSystem::N2Idx
+        contiH2OEqIdx = Indices::conti0EqIdx + H2OIdx,
+        contiN2EqIdx = Indices::conti0EqIdx + N2Idx
     };
 
     //! property that defines whether mole or mass fractions are used
@@ -202,7 +210,7 @@ public:
         // condition for the N2 molefraction at left boundary
         if (globalPos[0] < eps_ )
         {
-            values[FluidSystem::N2Idx] = 2.0e-5;
+            values[N2Idx] = 2.0e-5;
         }
 
         return values;
@@ -248,8 +256,8 @@ public:
         if(isBox && useNitscheTypeBc_)
         {
             flux[contiH2OEqIdx] = (volVars.pressure() - dirichletPressure) * 1e7;
-            flux[contiN2EqIdx] = flux[contiH2OEqIdx]  * (useMoles ? volVars.moleFraction(0, FluidSystem::N2Idx) :
-                                                                    volVars.massFraction(0, FluidSystem::N2Idx));
+            flux[contiN2EqIdx] = flux[contiH2OEqIdx]  * (useMoles ? volVars.moleFraction(0, N2Idx) :
+                                                                    volVars.massFraction(0, N2Idx));
             return flux;
         }
 
@@ -296,7 +304,7 @@ public:
         flux[contiH2OEqIdx] = tpfaFlux;
 
         // emulate an outflow condition for the component transport on the right side
-        flux[contiN2EqIdx] = tpfaFlux  * (useMoles ? volVars.moleFraction(0, FluidSystem::N2Idx) : volVars.massFraction(0, FluidSystem::N2Idx));
+        flux[contiN2EqIdx] = tpfaFlux  * (useMoles ? volVars.moleFraction(0, N2Idx) : volVars.massFraction(0, N2Idx));
 
         return flux;
     }
@@ -341,7 +349,7 @@ private:
     {
         PrimaryVariables priVars;
         priVars[pressureIdx] = 2e5 - 1e5*globalPos[0]; // initial condition for the pressure
-        priVars[FluidSystem::N2Idx] = 0.0;  // initial condition for the N2 molefraction
+        priVars[N2Idx] = 0.0;  // initial condition for the N2 molefraction
         return priVars;
     }
         static constexpr Scalar eps_ = 1e-6;
diff --git a/test/porousmediumflow/1pnc/implicit/test_1p2c_fv.cc b/test/porousmediumflow/1pnc/implicit/test_1p2c_fv.cc
index 9d795a1975..2259a2af78 100644
--- a/test/porousmediumflow/1pnc/implicit/test_1p2c_fv.cc
+++ b/test/porousmediumflow/1pnc/implicit/test_1p2c_fv.cc
@@ -111,7 +111,7 @@ int main(int argc, char** argv) try
     auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize");
 
     // intialize the vtk output module
-    VtkOutputModule<TypeTag, GET_PROP_VALUE(TypeTag, PhaseIdx)> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
+    VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
     using VtkOutputFields = typename GET_PROP_TYPE(TypeTag, VtkOutputFields);
     VtkOutputFields::init(vtkWriter); //!< Add model specific output fields
     vtkWriter.write(0.0);
diff --git a/test/porousmediumflow/1pnc/implicit/test_1p2cni_conduction_fv.cc b/test/porousmediumflow/1pnc/implicit/test_1p2cni_conduction_fv.cc
index f61f996cc7..52bf9985fe 100644
--- a/test/porousmediumflow/1pnc/implicit/test_1p2cni_conduction_fv.cc
+++ b/test/porousmediumflow/1pnc/implicit/test_1p2cni_conduction_fv.cc
@@ -110,7 +110,7 @@ int main(int argc, char** argv) try
     auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize");
 
     // intialize the vtk output module
-    VtkOutputModule<TypeTag, GET_PROP_VALUE(TypeTag, PhaseIdx)> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
+    VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
     using VtkOutputFields = typename GET_PROP_TYPE(TypeTag, VtkOutputFields);
     VtkOutputFields::init(vtkWriter); //!< Add model specific output fields
     vtkWriter.addField(problem->getExactTemperature(), "temperatureExact");
diff --git a/test/porousmediumflow/1pnc/implicit/test_1p2cni_convection_fv.cc b/test/porousmediumflow/1pnc/implicit/test_1p2cni_convection_fv.cc
index 1a555df69c..418edee9af 100644
--- a/test/porousmediumflow/1pnc/implicit/test_1p2cni_convection_fv.cc
+++ b/test/porousmediumflow/1pnc/implicit/test_1p2cni_convection_fv.cc
@@ -110,7 +110,7 @@ int main(int argc, char** argv) try
     auto maxDt = getParam<Scalar>("TimeLoop.MaxTimeStepSize");
 
     // intialize the vtk output module
-    VtkOutputModule<TypeTag, GET_PROP_VALUE(TypeTag, PhaseIdx)> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
+    VtkOutputModule<TypeTag> vtkWriter(*problem, *fvGridGeometry, *gridVariables, x, problem->name());
     using VtkOutputFields = typename GET_PROP_TYPE(TypeTag, VtkOutputFields);
     VtkOutputFields::init(vtkWriter); //!< Add model specific output fields
     vtkWriter.addField(problem->getExactTemperature(), "temperatureExact");
-- 
GitLab