diff --git a/dumux/adaptive/initializationindicator.hh b/dumux/adaptive/initializationindicator.hh
index b1e1461af9cfaf6e39442cdbb46db0af7a12e466..04580f67f8d8373e176991336d5cf0a882de3f4e 100644
--- a/dumux/adaptive/initializationindicator.hh
+++ b/dumux/adaptive/initializationindicator.hh
@@ -226,8 +226,8 @@ public:
                     // Get bcTypes and maybe mark for refinement on Dirichlet boundaries
                     for (const auto& scv : scvs(fvGeometry))
                     {
-                        bcTypes[scv.indexInElement()] = problem_->boundaryTypes(element, scv);
-                        if (refineAtDirichletBC_ && bcTypes[scv.indexInElement()].hasDirichlet())
+                        bcTypes[scv.localDofIndex()] = problem_->boundaryTypes(element, scv);
+                        if (refineAtDirichletBC_ && bcTypes[scv.localDofIndex()].hasDirichlet())
                         {
                             indicatorVector_[eIdx] = true;
                             break; // element is marked, escape scv loop
diff --git a/dumux/assembly/boxlocalassembler.hh b/dumux/assembly/boxlocalassembler.hh
index adab4567da70eaee8c55ebdbeb2196ad728446e6..38cb6e678d6f61077121780f166147b2086ec266 100644
--- a/dumux/assembly/boxlocalassembler.hh
+++ b/dumux/assembly/boxlocalassembler.hh
@@ -86,13 +86,13 @@ public:
         {
             const auto residual = this->asImp_().evalLocalResidual(); // forward to the internal implementation
             for (const auto& scv : scvs(this->fvGeometry()))
-                res[scv.dofIndex()] += residual[scv.indexInElement()];
+                res[scv.dofIndex()] += residual[scv.localDofIndex()];
         }
         else
         {
             const auto residual = this->asImp_().assembleJacobianAndResidualImpl(jac, gridVariables, partialReassembler); // forward to the internal implementation
             for (const auto& scv : scvs(this->fvGeometry()))
-                res[scv.dofIndex()] += residual[scv.indexInElement()];
+                res[scv.dofIndex()] += residual[scv.localDofIndex()];
         }
 
         auto applyDirichlet = [&] (const auto& scvI,
@@ -104,7 +104,7 @@ public:
             for (const auto& scvJ : scvs(this->fvGeometry()))
             {
                 jac[scvI.dofIndex()][scvJ.dofIndex()][eqIdx] = 0.0;
-                if (scvI.indexInElement() == scvJ.indexInElement())
+                if (scvI.localDofIndex() == scvJ.localDofIndex())
                     jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
             }
         };
@@ -129,7 +129,7 @@ public:
             for (const auto& scvJ : scvs(this->fvGeometry()))
             {
                 jac[scvI.dofIndex()][scvJ.dofIndex()][eqIdx] = 0.0;
-                if (scvI.indexInElement() == scvJ.indexInElement())
+                if (scvI.localDofIndex() == scvJ.localDofIndex())
                     jac[scvI.dofIndex()][scvI.dofIndex()][eqIdx][pvIdx] = 1.0;
             }
         };
@@ -146,7 +146,7 @@ public:
         const auto residual = this->evalLocalResidual();
 
         for (const auto& scv : scvs(this->fvGeometry()))
-            res[scv.dofIndex()] += residual[scv.indexInElement()];
+            res[scv.dofIndex()] += residual[scv.localDofIndex()];
 
         auto applyDirichlet = [&] (const auto& scvI,
                                    const auto& dirichletValues,
@@ -171,7 +171,7 @@ public:
         {
             for (const auto& scvI : scvs(this->fvGeometry()))
             {
-                const auto bcTypes = this->elemBcTypes()[scvI.indexInElement()];
+                const auto bcTypes = this->elemBcTypes()[scvI.localDofIndex()];
                 if (bcTypes.hasDirichlet())
                 {
                     const auto dirichletValues = this->problem().dirichlet(this->element(), scvI);
@@ -281,7 +281,7 @@ public:
                 auto evalResiduals = [&](Scalar priVar)
                 {
                     // update the volume variables and compute element residual
-                    elemSol[scv.indexInElement()][pvIdx] = priVar;
+                    elemSol[scv.localDofIndex()][pvIdx] = priVar;
                     curVolVars.update(elemSol, this->problem(), element, scv);
                     return this->evalLocalResidual();
                 };
@@ -289,8 +289,8 @@ public:
                 // derive the residuals numerically
                 static const NumericEpsilon<Scalar, numEq> eps_{GET_PROP_VALUE(TypeTag, ModelParameterGroup)};
                 static const int numDiffMethod = getParamFromGroup<int>(GET_PROP_VALUE(TypeTag, ModelParameterGroup), "Assembly.NumericDifferenceMethod");
-                NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.indexInElement()][pvIdx], partialDerivs, origResiduals,
-                                                          eps_(elemSol[scv.indexInElement()][pvIdx], pvIdx), numDiffMethod);
+                NumericDifferentiation::partialDerivative(evalResiduals, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
+                                                          eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
 
                 // update the global stiffness matrix with the current partial derivatives
                 for (auto&& scvJ : scvs(fvGeometry))
@@ -305,7 +305,7 @@ public:
                             // the residual of equation 'eqIdx' at dof 'i'
                             // depending on the primary variable 'pvIdx' at dof
                             // 'col'.
-                            A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.indexInElement()][eqIdx];
+                            A[scvJ.dofIndex()][dofIdx][eqIdx][pvIdx] += partialDerivs[scvJ.localDofIndex()][eqIdx];
                         }
                     }
                 }
@@ -314,7 +314,7 @@ public:
                 curVolVars = origVolVars;
 
                 // restore the original element solution
-                elemSol[scv.indexInElement()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
+                elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
                 // TODO additional dof dependencies
             }
         }
@@ -401,7 +401,7 @@ public:
                 auto evalStorage = [&](Scalar priVar)
                 {
                     // auto partialDerivsTmp = partialDerivs;
-                    elemSol[scv.indexInElement()][pvIdx] = priVar;
+                    elemSol[scv.localDofIndex()][pvIdx] = priVar;
                     curVolVars.update(elemSol, this->problem(), element, scv);
                     return this->evalLocalStorageResidual();
                 };
@@ -409,8 +409,8 @@ public:
                 // derive the residuals numerically
                 static const NumericEpsilon<Scalar, numEq> eps_{GET_PROP_VALUE(TypeTag, ModelParameterGroup)};
                 static const int numDiffMethod = getParamFromGroup<int>(GET_PROP_VALUE(TypeTag, ModelParameterGroup), "Assembly.NumericDifferenceMethod");
-                NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.indexInElement()][pvIdx], partialDerivs, origResiduals,
-                                                          eps_(elemSol[scv.indexInElement()][pvIdx], pvIdx), numDiffMethod);
+                NumericDifferentiation::partialDerivative(evalStorage, elemSol[scv.localDofIndex()][pvIdx], partialDerivs, origResiduals,
+                                                          eps_(elemSol[scv.localDofIndex()][pvIdx], pvIdx), numDiffMethod);
 
                 // update the global stiffness matrix with the current partial derivatives
                 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
@@ -419,14 +419,14 @@ public:
                     // the residual of equation 'eqIdx' at dof 'i'
                     // depending on the primary variable 'pvIdx' at dof
                     // 'col'.
-                    A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.indexInElement()][eqIdx];
+                    A[dofIdx][dofIdx][eqIdx][pvIdx] += partialDerivs[scv.localDofIndex()][eqIdx];
                 }
 
                 // restore the original state of the scv's volume variables
                 curVolVars = origVolVars;
 
                 // restore the original element solution
-                elemSol[scv.indexInElement()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
+                elemSol[scv.localDofIndex()][pvIdx] = curSol[scv.dofIndex()][pvIdx];
                 // TODO additional dof dependencies
             }
         }
@@ -534,7 +534,7 @@ public:
             else
             {
                 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
-                if (this->elemBcTypes()[insideScv.indexInElement()].hasNeumann())
+                if (this->elemBcTypes()[insideScv.localDofIndex()].hasNeumann())
                 {
                     // add flux term derivatives
                     this->localResidual().addRobinFluxDerivatives(A[insideScv.dofIndex()],
diff --git a/dumux/assembly/boxlocalresidual.hh b/dumux/assembly/boxlocalresidual.hh
index 12fd39ca72b18093dbd8f23d7fdb92ddf66f9e81..eaadcedddf6ae9fd64388e42bb5e31f84f77775b 100644
--- a/dumux/assembly/boxlocalresidual.hh
+++ b/dumux/assembly/boxlocalresidual.hh
@@ -73,13 +73,13 @@ public:
         {
             const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
             const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
-            residual[insideScv.indexInElement()] += flux;
-            residual[outsideScv.indexInElement()] -= flux;
+            residual[insideScv.localDofIndex()] += flux;
+            residual[outsideScv.localDofIndex()] -= flux;
         }
         else
         {
             const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
-            residual[insideScv.indexInElement()] += flux;
+            residual[insideScv.localDofIndex()] += flux;
         }
     }
 
@@ -104,7 +104,7 @@ public:
         else
         {
             const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
-            const auto& bcTypes = elemBcTypes[scv.indexInElement()];
+            const auto& bcTypes = elemBcTypes[scv.localDofIndex()];
 
             // Neumann and Robin ("solution dependent Neumann") boundary conditions
             if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
diff --git a/dumux/assembly/cclocalresidual.hh b/dumux/assembly/cclocalresidual.hh
index 03314a9ca57eb7b4d050e82b67a21adaebae2a03..f78a759d21d5a38960dabe6424e7efe4ba2c6782 100644
--- a/dumux/assembly/cclocalresidual.hh
+++ b/dumux/assembly/cclocalresidual.hh
@@ -64,7 +64,7 @@ public:
                   const SubControlVolumeFace& scvf) const
     {
         const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
-        const auto localScvIdx = scv.indexInElement();
+        const auto localScvIdx = scv.localDofIndex();
         residual[localScvIdx] += this->asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
     }
 
diff --git a/dumux/assembly/fvlocalresidual.hh b/dumux/assembly/fvlocalresidual.hh
index 73a5bbc0abaa2e393deff679449231aa46f4f1ba..cb1a30e54111df61840953b9bcfe6a923360bfe4 100644
--- a/dumux/assembly/fvlocalresidual.hh
+++ b/dumux/assembly/fvlocalresidual.hh
@@ -112,7 +112,7 @@ public:
         // all sub control volumes
         for (auto&& scv : scvs(fvGeometry))
         {
-            auto localScvIdx = scv.indexInElement();
+            auto localScvIdx = scv.localDofIndex();
             const auto& volVars = elemVolVars[scv];
             storage[localScvIdx] = asImp().computeStorage(scv, volVars);
             storage[localScvIdx] *= scv.volume() * volVars.extrusionFactor();
@@ -382,7 +382,7 @@ public:
         storage *= scv.volume();
         storage /= timeLoop_->timeStepSize();
 
-        residual[scv.indexInElement()] += storage;
+        residual[scv.localDofIndex()] += storage;
     }
 
     /*!
@@ -411,7 +411,7 @@ public:
         source *= scv.volume()*curVolVars.extrusionFactor();
 
         //! subtract source from local rate (sign convention in user interface)
-        residual[scv.indexInElement()] -= source;
+        residual[scv.localDofIndex()] -= source;
     }
 
     /*!
diff --git a/dumux/common/fvproblem.hh b/dumux/common/fvproblem.hh
index 3e9d4c316e8e4fb5275dfab8fc42ddf99b042611..fb7aa0e46c878a8c82beb1adfce0dbda164431e8 100644
--- a/dumux/common/fvproblem.hh
+++ b/dumux/common/fvproblem.hh
@@ -407,7 +407,7 @@ public:
                                 const SubControlVolume &scv) const
     {
         NumEqVector source(0);
-        auto scvIdx = scv.indexInElement();
+        auto scvIdx = scv.localDofIndex();
         auto key = std::make_pair(fvGridGeometry_->elementMapper().index(element), scvIdx);
         if (pointSourceMap_.count(key))
         {
diff --git a/dumux/common/pointsource.hh b/dumux/common/pointsource.hh
index b4bdac1fad551b773e40267a3f1520f7d93be48c..97961938cf54afa8806fe1d51f524d99afecb9cd 100644
--- a/dumux/common/pointsource.hh
+++ b/dumux/common/pointsource.hh
@@ -303,7 +303,7 @@ public:
                     for (auto&& scv : scvs(fvGeometry))
                     {
                         if (intersectsPointGeometry(globalPos, scv.geometry()))
-                            scvIndices.push_back(scv.indexInElement());
+                            scvIndices.push_back(scv.localDofIndex());
                     }
                     // for all scvs that where tested positiv add the point sources
                     // to the element/scv to point source map
diff --git a/dumux/discretization/box/darcyslaw.hh b/dumux/discretization/box/darcyslaw.hh
index bf7d53b48bcc94327ec8906f2be205d6e7f5ebf1..c8136ad8af4a679a596425b9e9b496c6db9dfd69 100644
--- a/dumux/discretization/box/darcyslaw.hh
+++ b/dumux/discretization/box/darcyslaw.hh
@@ -98,10 +98,10 @@ public:
             const auto& volVars = elemVolVars[scv];
 
             if (enableGravity)
-                rho += volVars.density(phaseIdx)*shapeValues[scv.indexInElement()][0];
+                rho += volVars.density(phaseIdx)*shapeValues[scv.localDofIndex()][0];
 
             // the global shape function gradient
-            gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.indexInElement()));
+            gradP.axpy(volVars.pressure(phaseIdx), fluxVarCache.gradN(scv.localDofIndex()));
         }
 
         if (enableGravity)
@@ -135,8 +135,8 @@ public:
 
         std::vector<Scalar> ti(fvGeometry.numScv());
         for (const auto& scv : scvs(fvGeometry))
-            ti[scv.indexInElement()] =
-                -1.0*scvf.area()*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.indexInElement()));
+            ti[scv.localDofIndex()] =
+                -1.0*scvf.area()*vtmv(scvf.unitOuterNormal(), K, fluxVarCache.gradN(scv.localDofIndex()));
 
         return ti;
     }
diff --git a/dumux/discretization/box/elementboundarytypes.hh b/dumux/discretization/box/elementboundarytypes.hh
index e7d8eb45d814e4d609f8c852337b18617a5eac7e..fdcd924e939239f34caaf09522c2bed031259b82 100644
--- a/dumux/discretization/box/elementboundarytypes.hh
+++ b/dumux/discretization/box/elementboundarytypes.hh
@@ -62,7 +62,7 @@ public:
 
         for (const auto& scv : scvs(fvGeometry))
         {
-            int scvIdxLocal = scv.indexInElement();
+            int scvIdxLocal = scv.localDofIndex();
             vertexBCTypes_[scvIdxLocal].reset();
 
             if (fvGeometry.fvGridGeometry().dofOnBoundary(scv.dofIndex()))
diff --git a/dumux/discretization/box/elementsolution.hh b/dumux/discretization/box/elementsolution.hh
index d77d714ca2b9c9ede53a9d1b13ce3065c5d3b722..0bfd6e9136eecf52a8a6073555a2c95b9604260a 100644
--- a/dumux/discretization/box/elementsolution.hh
+++ b/dumux/discretization/box/elementsolution.hh
@@ -63,7 +63,7 @@ public:
         const auto numVert = element.subEntities(GridView::dimension);
         priVars_.resize(numVert);
         for (const auto& scv : scvs(fvGeometry))
-            priVars_[scv.indexInElement()] = elemVolVars[scv].priVars();
+            priVars_[scv.localDofIndex()] = elemVolVars[scv].priVars();
     }
 
     //! extract the element solution from the solution vector using a mapper
@@ -85,7 +85,7 @@ public:
         const auto numVert = element.subEntities(GridView::dimension);
         priVars_.resize(numVert);
         for (const auto& scv : scvs(fvGeometry))
-            priVars_[scv.indexInElement()] = sol[scv.dofIndex()];
+            priVars_[scv.localDofIndex()] = sol[scv.dofIndex()];
     }
 
     //! bracket operator const access
diff --git a/dumux/discretization/box/elementvolumevariables.hh b/dumux/discretization/box/elementvolumevariables.hh
index 1c9cc0ca05d873658f8c4921ff4ac305c03db0c6..6c67ef0bc3e6c8f384d55a4c5ce0a30a94157776 100644
--- a/dumux/discretization/box/elementvolumevariables.hh
+++ b/dumux/discretization/box/elementvolumevariables.hh
@@ -63,7 +63,7 @@ public:
 
     template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
     const VolumeVariables& operator [](const SubControlVolume& scv) const
-    { return gridVolVars().volVars(eIdx_, scv.indexInElement()); }
+    { return gridVolVars().volVars(eIdx_, scv.localDofIndex()); }
 
     // For compatibility reasons with the case of not storing the vol vars.
     // function to be called before assembling an element, preparing the vol vars within the stencil
@@ -133,7 +133,7 @@ public:
         // resize volume variables to the required size
         volumeVariables_.resize(fvGeometry.numScv());
         for (auto&& scv : scvs(fvGeometry))
-            volumeVariables_[scv.indexInElement()].update(elemSol, gridVolVars().problem(), element, scv);
+            volumeVariables_[scv.localDofIndex()].update(elemSol, gridVolVars().problem(), element, scv);
     }
 
     const VolumeVariables& operator [](std::size_t scvIdx) const
@@ -144,11 +144,11 @@ public:
 
     template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
     const VolumeVariables& operator [](const SubControlVolume& scv) const
-    { return volumeVariables_[scv.indexInElement()]; }
+    { return volumeVariables_[scv.localDofIndex()]; }
 
     template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
     VolumeVariables& operator [](const SubControlVolume& scv)
-    { return volumeVariables_[scv.indexInElement()]; }
+    { return volumeVariables_[scv.localDofIndex()]; }
 
     //! The global volume variables object we are a restriction of
     const GridVolumeVariables& gridVolVars() const
diff --git a/dumux/discretization/box/fickslaw.hh b/dumux/discretization/box/fickslaw.hh
index 34b2e7add035f472b2521b4e68fb030dc2ae921e..c78cec06311655bedec09d69b943df7d3d54be00 100644
--- a/dumux/discretization/box/fickslaw.hh
+++ b/dumux/discretization/box/fickslaw.hh
@@ -93,7 +93,7 @@ public:
         // density interpolation
         Scalar rho(0.0);
         for (auto&& scv : scvs(fvGeometry))
-            rho += elemVolVars[scv].molarDensity(phaseIdx)*shapeValues[scv.indexInElement()][0];
+            rho += elemVolVars[scv].molarDensity(phaseIdx)*shapeValues[scv.localDofIndex()][0];
 
         for (int compIdx = 0; compIdx < numComponents; compIdx++)
         {
@@ -119,7 +119,7 @@ public:
             // the mole/mass fraction gradient
             GlobalPosition gradX(0.0);
             for (auto&& scv : scvs(fvGeometry))
-                gradX.axpy(elemVolVars[scv].moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.indexInElement()));
+                gradX.axpy(elemVolVars[scv].moleFraction(phaseIdx, compIdx), fluxVarsCache.gradN(scv.localDofIndex()));
 
             // compute the diffusive flux
             componentFlux[compIdx] = -1.0*rho*vtmv(scvf.unitOuterNormal(), D, gradX)*scvf.area();
@@ -142,7 +142,7 @@ public:
         Scalar rho(0.0);
         const auto& shapeValues = fluxVarCache.shapeValues();
         for (auto&& scv : scvs(fvGeometry))
-            rho += elemVolVars[scv].molarDensity(phaseIdx)*shapeValues[scv.indexInElement()][0];
+            rho += elemVolVars[scv].molarDensity(phaseIdx)*shapeValues[scv.localDofIndex()][0];
 
         const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
         const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
@@ -171,7 +171,7 @@ public:
 
             ti[compIdx].resize(fvGeometry.numScv());
             for (auto&& scv : scvs(fvGeometry))
-                ti[compIdx][scv.indexInElement()] = -rho*vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.indexInElement()))*scvf.area();
+                ti[compIdx][scv.localDofIndex()] = -rho*vtmv(scvf.unitOuterNormal(), D, fluxVarCache.gradN(scv.localDofIndex()))*scvf.area();
         }
 
         return ti;
diff --git a/dumux/discretization/box/fourierslaw.hh b/dumux/discretization/box/fourierslaw.hh
index 475b0adb370edfecd34787e523148a3ea7c08535..df78a006a093ddb2f0b32c1736c52550513fa2d5 100644
--- a/dumux/discretization/box/fourierslaw.hh
+++ b/dumux/discretization/box/fourierslaw.hh
@@ -96,7 +96,7 @@ public:
         // compute the temperature gradient with the shape functions
         GlobalPosition gradTemp(0.0);
         for (auto&& scv : scvs(fvGeometry))
-            gradTemp.axpy(elemVolVars[scv].temperature(), fluxVarsCache.gradN(scv.indexInElement()));
+            gradTemp.axpy(elemVolVars[scv].temperature(), fluxVarsCache.gradN(scv.localDofIndex()));
 
         // comute the heat conduction flux
         return -1.0*vtmv(scvf.unitOuterNormal(), lambda, gradTemp)*scvf.area();
diff --git a/dumux/discretization/box/fourierslawnonequilibrium.hh b/dumux/discretization/box/fourierslawnonequilibrium.hh
index 9422b57c4c681e01136329004af3fbd780c0fe6b..76e2ba58c1a8e434570ce56f29d9dd90dae4cb93 100644
--- a/dumux/discretization/box/fourierslawnonequilibrium.hh
+++ b/dumux/discretization/box/fourierslawnonequilibrium.hh
@@ -122,9 +122,9 @@ public:
         {
             // compute the temperature gradient with the shape functions
             if (phaseIdx < numEnergyEqFluid)
-                gradTemp.axpy(elemVolVars[scv].temperature(phaseIdx), fluxVarsCache.gradN(scv.indexInElement()));
+                gradTemp.axpy(elemVolVars[scv].temperature(phaseIdx), fluxVarsCache.gradN(scv.localDofIndex()));
             else
-               gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.indexInElement()));
+               gradTemp.axpy(elemVolVars[scv].temperatureSolid(), fluxVarsCache.gradN(scv.localDofIndex()));
         }
 
         // comute the heat conduction flux
diff --git a/dumux/discretization/box/gridvolumevariables.hh b/dumux/discretization/box/gridvolumevariables.hh
index a5fe851d5eeb3b65cc9d5cfe6ed5e332ec3c4774..5cf94cf4fda826da3956a71d0b51826475ad8598 100644
--- a/dumux/discretization/box/gridvolumevariables.hh
+++ b/dumux/discretization/box/gridvolumevariables.hh
@@ -88,17 +88,17 @@ public:
             // update the volvars of the element
             volumeVariables_[eIdx].resize(fvGeometry.numScv());
             for (auto&& scv : scvs(fvGeometry))
-                volumeVariables_[eIdx][scv.indexInElement()].update(elemSol, problem(), element, scv);
+                volumeVariables_[eIdx][scv.localDofIndex()].update(elemSol, problem(), element, scv);
         }
     }
 
     template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
     const VolumeVariables& volVars(const SubControlVolume& scv) const
-    { return volumeVariables_[scv.elementIndex()][scv.indexInElement()]; }
+    { return volumeVariables_[scv.elementIndex()][scv.localDofIndex()]; }
 
     template<class SubControlVolume, typename std::enable_if_t<!std::is_integral<SubControlVolume>::value, int> = 0>
     VolumeVariables& volVars(const SubControlVolume& scv)
-    { return volumeVariables_[scv.elementIndex()][scv.indexInElement()]; }
+    { return volumeVariables_[scv.elementIndex()][scv.localDofIndex()]; }
 
     const VolumeVariables& volVars(const std::size_t eIdx, const std::size_t scvIdx) const
     { return volumeVariables_[eIdx][scvIdx]; }
diff --git a/dumux/discretization/box/maxwellstefanslaw.hh b/dumux/discretization/box/maxwellstefanslaw.hh
index b0f61b7a7d5f2109f567ce13b4aed8f06b99d740..571aa5bb2b7e87d85e0df0f83979d00a7f42aacf 100644
--- a/dumux/discretization/box/maxwellstefanslaw.hh
+++ b/dumux/discretization/box/maxwellstefanslaw.hh
@@ -103,14 +103,14 @@ public:
             const auto& volVars = elemVolVars[scv];
 
             // density interpolation
-            rho +=  volVars.molarDensity(phaseIdx)*shapeValues[scv.indexInElement()][0];
+            rho +=  volVars.molarDensity(phaseIdx)*shapeValues[scv.localDofIndex()][0];
             // the ansatz function gradient
-            jacInvT.mv(shapeJacobian[scv.indexInElement()][0], gradN[scv.indexInElement()]);
+            jacInvT.mv(shapeJacobian[scv.localDofIndex()][0], gradN[scv.localDofIndex()]);
 
             //interpolate the mole fraction for the diffusion matrix
             for (int compIdx = 0; compIdx < numComponents; compIdx++)
             {
-              moleFrac[compIdx] += volVars.moleFraction(phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0];
+              moleFrac[compIdx] += volVars.moleFraction(phaseIdx, compIdx)*shapeValues[scv.localDofIndex()][0];
             }
         }
 
@@ -124,7 +124,7 @@ public:
                 const auto& volVars = elemVolVars[scv];
 
                 // the mole/mass fraction gradient
-                gradX.axpy(volVars.moleFraction(phaseIdx, compIdx), gradN[scv.indexInElement()]);
+                gradX.axpy(volVars.moleFraction(phaseIdx, compIdx), gradN[scv.localDofIndex()]);
             }
 
            normalX[compIdx] = gradX *scvf.unitOuterNormal();
diff --git a/dumux/discretization/box/subcontrolvolume.hh b/dumux/discretization/box/subcontrolvolume.hh
index cfcfe504a4c8e2afbe5fcab2f1874b707bdc3bf3..6040e904cb63e05c79649c50e3c4e016d45d608e 100644
--- a/dumux/discretization/box/subcontrolvolume.hh
+++ b/dumux/discretization/box/subcontrolvolume.hh
@@ -146,7 +146,7 @@ public:
     }
 
     //! The global index of this scv
-    LocalIndexType indexInElement() const
+    LocalIndexType localDofIndex() const
     {
         return scvIdx_;
     }
diff --git a/dumux/discretization/cellcentered/subcontrolvolume.hh b/dumux/discretization/cellcentered/subcontrolvolume.hh
index 8ec5d6c7d03c60bcc710b3fa9e4be063f420f57e..22dbc4dccace654c0f78e3a48117f4b898053514 100644
--- a/dumux/discretization/cellcentered/subcontrolvolume.hh
+++ b/dumux/discretization/cellcentered/subcontrolvolume.hh
@@ -139,7 +139,7 @@ public:
     }
 
     //! The global index of this scv
-    LocalIndexType indexInElement() const
+    LocalIndexType localDofIndex() const
     {
         return 0;
     }
diff --git a/dumux/discretization/subcontrolvolumebase.hh b/dumux/discretization/subcontrolvolumebase.hh
index 6479200dd3092c1199e10ca667411800ab824d1d..7d993f9d797394a1db629bc4725a2c00a0d1ba04 100644
--- a/dumux/discretization/subcontrolvolumebase.hh
+++ b/dumux/discretization/subcontrolvolumebase.hh
@@ -65,9 +65,9 @@ public:
     }
 
     //! The index of the dof this scv is embedded in (box)
-    LocalIndexType indexInElement() const
+    LocalIndexType localDofIndex() const
     {
-        return asImp_().indexInElement();
+        return asImp_().localDofIndex();
     }
 
     // The position of the dof this scv is embedded in
diff --git a/dumux/porousmediumflow/1p/incompressiblelocalresidual.hh b/dumux/porousmediumflow/1p/incompressiblelocalresidual.hh
index 72237fdb8ab11f0b86faafb829343d24b5488446..001549dd23b6110dedcc08a887788156b62916b8 100644
--- a/dumux/porousmediumflow/1p/incompressiblelocalresidual.hh
+++ b/dumux/porousmediumflow/1p/incompressiblelocalresidual.hh
@@ -195,7 +195,7 @@ public:
                                  / curElemVolVars[scvf.insideScvIdx()].viscosity();
         for (const auto& scv : scvs(fvGeometry))
         {
-            auto d = up*ti[scv.indexInElement()];
+            auto d = up*ti[scv.localDofIndex()];
             A[insideScv.dofIndex()][scv.dofIndex()][conti0EqIdx][pressureIdx] += d;
             A[outsideScv.dofIndex()][scv.dofIndex()][conti0EqIdx][pressureIdx] -= d;
         }
diff --git a/dumux/porousmediumflow/2p/griddatatransfer.hh b/dumux/porousmediumflow/2p/griddatatransfer.hh
index 748ffedeb9808c56e6523f13a45ba17bdb27f356..df4bf68def921b26ccd3415f7968b3a0bed34e5c 100644
--- a/dumux/porousmediumflow/2p/griddatatransfer.hh
+++ b/dumux/porousmediumflow/2p/griddatatransfer.hh
@@ -234,7 +234,7 @@ public:
                     volVars.update(elemSol, *problem_, element, scv);
 
                     // write solution at dof in current solution vector
-                    sol_[scv.dofIndex()] = elemSol[scv.indexInElement()];
+                    sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
 
                     const auto dofIdxGlobal = scv.dofIndex();
                     // For cc schemes, overwrite the saturation by a mass conservative one here
@@ -329,7 +329,7 @@ public:
                     ElementSolution elemSolSon(element, sol_, *fvGridGeometry_);
                     const auto fatherGeometry = fatherElement.geometry();
                     for (const auto& scv : scvs(fvGeometry))
-                        elemSolSon[scv.indexInElement()] = evalSolution(fatherElement,
+                        elemSolSon[scv.localDofIndex()] = evalSolution(fatherElement,
                                                                         fatherGeometry,
                                                                         adaptedValuesFather.u,
                                                                         scv.dofPosition());
@@ -355,7 +355,7 @@ public:
                         }
 
                         // store constructed (pressure) values of son in the current solution (saturation comes later)
-                        sol_[dofIdxGlobal] = elemSolSon[scv.indexInElement()];
+                        sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
                     }
                 }
             }
diff --git a/dumux/porousmediumflow/2p/incompressiblelocalresidual.hh b/dumux/porousmediumflow/2p/incompressiblelocalresidual.hh
index 421b5ce7af25fdc5018db0e6b7f4730f930d9b31..c873858c120254401464cb0e938696ed0a314c1e 100644
--- a/dumux/porousmediumflow/2p/incompressiblelocalresidual.hh
+++ b/dumux/porousmediumflow/2p/incompressiblelocalresidual.hh
@@ -345,7 +345,7 @@ public:
         for (const auto& scvJ : scvs(fvGeometry))
         {
             const auto globalJ = scvJ.dofIndex();
-            const auto localJ = scvJ.indexInElement();
+            const auto localJ = scvJ.localDofIndex();
 
             // the transmissibily associated with the scvJ
             const auto tj = ti[localJ];
diff --git a/dumux/porousmediumflow/fluxvariablescache.hh b/dumux/porousmediumflow/fluxvariablescache.hh
index 32699d57d963f18931851a77f24bd16bb24fb9e5..f92143eaa6eeab12fb0852847b24651d94c71286 100644
--- a/dumux/porousmediumflow/fluxvariablescache.hh
+++ b/dumux/porousmediumflow/fluxvariablescache.hh
@@ -99,7 +99,7 @@ public:
         // compute the gradN at for every scv/dof
         gradN_.resize(fvGeometry.numScv());
         for (const auto& scv: scvs(fvGeometry))
-            jacInvT_.mv(shapeJacobian_[scv.indexInElement()][0], gradN_[scv.indexInElement()]);
+            jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.localDofIndex()]);
 
     }
 
diff --git a/dumux/porousmediumflow/mineralization/volumevariables.hh b/dumux/porousmediumflow/mineralization/volumevariables.hh
index 1855d728961dedc0f39e7e9c25dc6c12522adc5a..87a55ba74e0affb58b5c879e499725fd32fe2b4c 100644
--- a/dumux/porousmediumflow/mineralization/volumevariables.hh
+++ b/dumux/porousmediumflow/mineralization/volumevariables.hh
@@ -56,7 +56,7 @@ public:
         ParentType::update(elemSol, problem, element, scv);
 
         // calculate the remaining quantities
-        auto&& priVars = elemSol[scv.indexInElement()];
+        auto&& priVars = elemSol[scv.localDofIndex()];
 
         sumPrecipitates_ = 0.0;
         for(int sPhaseIdx = 0; sPhaseIdx < ModelTraits::numSPhases(); ++sPhaseIdx)
diff --git a/dumux/porousmediumflow/mpnc/localresidual.hh b/dumux/porousmediumflow/mpnc/localresidual.hh
index 14cb90a15f6a11ef8ab25551a1dc06f4c3af3a59..31e1d1db28c190635ce1717140f286159884de89 100644
--- a/dumux/porousmediumflow/mpnc/localresidual.hh
+++ b/dumux/porousmediumflow/mpnc/localresidual.hh
@@ -72,7 +72,7 @@ public:
         for (auto&& scv : scvs(fvGeometry))
         {
             //here we need to set the constraints of the mpnc model into the residual
-            const auto localScvIdx = scv.indexInElement();
+            const auto localScvIdx = scv.localDofIndex();
             for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
                 residual[localScvIdx][phase0NcpIdx + phaseIdx] = elemVolVars[scv].phaseNcp(phaseIdx);
         }
diff --git a/dumux/porousmediumflow/nonequilibrium/localresidual.hh b/dumux/porousmediumflow/nonequilibrium/localresidual.hh
index 5ac069eeec690471443d6922376762e10d351599..089815a3b89bfe8f79626ed654c3d99c2b798bdc 100644
--- a/dumux/porousmediumflow/nonequilibrium/localresidual.hh
+++ b/dumux/porousmediumflow/nonequilibrium/localresidual.hh
@@ -281,7 +281,7 @@ public:
         // In the case of a kinetic consideration, mass transfer
         // between phases is realized via source terms there is a
         // balance equation for each component in each phase
-        const auto& localScvIdx = scv.indexInElement();
+        const auto& localScvIdx = scv.localDofIndex();
         const auto& volVars = elemVolVars[localScvIdx];
         ComponentVector componentIntoPhaseMassTransfer[numPhases];
 #define FUNKYMASSTRANSFER 0
diff --git a/dumux/porousmediumflow/nonequilibrium/thermal/localresidual.hh b/dumux/porousmediumflow/nonequilibrium/thermal/localresidual.hh
index 6a4fff355eb393ae46aef8c8bfbd204db3479a40..cc81d628a259149c7447f5f042c1dae43bf14902 100644
--- a/dumux/porousmediumflow/nonequilibrium/thermal/localresidual.hh
+++ b/dumux/porousmediumflow/nonequilibrium/thermal/localresidual.hh
@@ -164,7 +164,7 @@ public:
                                     const SubControlVolume &scv)
     {
         //specialization for 2 fluid phases
-        const auto& localScvIdx = scv.indexInElement();
+        const auto& localScvIdx = scv.localDofIndex();
         const auto& volVars = elemVolVars[localScvIdx];
         const auto& fs = volVars.fluidState() ;
         const Scalar characteristicLength = volVars.characteristicLength()  ;
@@ -415,7 +415,7 @@ public:
                                     const SubControlVolume &scv)
     {
         //specialization for 2 fluid phases
-        const auto &localScvIdx = scv.indexInElement();
+        const auto &localScvIdx = scv.localDofIndex();
         const auto &volVars = elemVolVars[localScvIdx];
 
         const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
diff --git a/dumux/porousmediumflow/nonequilibrium/volumevariables.hh b/dumux/porousmediumflow/nonequilibrium/volumevariables.hh
index b1c113862f4c235d72316bcae0adcceebf1b3d4c..801de358020c9a420939bfabb3a67b6f1d4140ef 100644
--- a/dumux/porousmediumflow/nonequilibrium/volumevariables.hh
+++ b/dumux/porousmediumflow/nonequilibrium/volumevariables.hh
@@ -174,20 +174,20 @@ public:
             // retrieve temperature from solution vector
             for(int phaseIdx=0; phaseIdx < numEnergyEqFluid; ++phaseIdx)
             {
-                const auto T = elemSol[scv.indexInElement()][Indices::temperature0Idx + phaseIdx];
+                const auto T = elemSol[scv.localDofIndex()][Indices::temperature0Idx + phaseIdx];
                 fluidState.setTemperature(phaseIdx, T);
             }
         }
         else
         {
-            const auto T = elemSol[scv.indexInElement()][Indices::temperature0Idx];
+            const auto T = elemSol[scv.localDofIndex()][Indices::temperature0Idx];
             fluidState.setTemperature(T);
         }
 
         // set solid temperatures
         static_assert(numEnergyEqSolid == 1, "Solid system not yet implemented");
         for (int solidPhaseIdx = numEnergyEqFluid; solidPhaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++solidPhaseIdx)
-            temperatureSolid_ = elemSol[scv.indexInElement()][Indices::temperature0Idx + solidPhaseIdx];
+            temperatureSolid_ = elemSol[scv.localDofIndex()][Indices::temperature0Idx + solidPhaseIdx];
     }
 
     /*!
@@ -606,12 +606,12 @@ public:
         for(int phaseIdx=0; phaseIdx < numEnergyEqFluid; ++phaseIdx)
         {
             // retrieve temperature from solution vector
-            const Scalar T = elemSol[scv.indexInElement()][Indices::temperature0Idx + phaseIdx];
+            const Scalar T = elemSol[scv.localDofIndex()][Indices::temperature0Idx + phaseIdx];
             fluidState.setTemperature(phaseIdx, T);
         }
 
         for(int solidPhaseIdx = numEnergyEqFluid; solidPhaseIdx < numEnergyEqFluid+numEnergyEqSolid; ++solidPhaseIdx)
-            temperatureSolid_ = elemSol[scv.indexInElement()][Indices::temperature0Idx + solidPhaseIdx];
+            temperatureSolid_ = elemSol[scv.localDofIndex()][Indices::temperature0Idx + solidPhaseIdx];
     }
 
     /*!
diff --git a/dumux/porousmediumflow/tracer/localresidual.hh b/dumux/porousmediumflow/tracer/localresidual.hh
index 5783b45f07e530d52cc3dfc160fda869e04b0f2a..0efdd882a15474f8a0653a68f8dba3263fa3ccaf 100644
--- a/dumux/porousmediumflow/tracer/localresidual.hh
+++ b/dumux/porousmediumflow/tracer/localresidual.hh
@@ -297,8 +297,8 @@ public:
             for (const auto& scv : scvs(fvGeometry))
             {
                 // diffusive term
-                const auto diffDeriv = useMoles ? ti[compIdx][scv.indexInElement()]
-                                                : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
+                const auto diffDeriv = useMoles ? ti[compIdx][scv.localDofIndex()]
+                                                : ti[compIdx][scv.localDofIndex()]*FluidSystem::molarMass(compIdx);
                 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
                 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
             }
diff --git a/dumux/porousmediumflow/velocityoutput.hh b/dumux/porousmediumflow/velocityoutput.hh
index cb10fe222a47853c0903f559014dbc97aef3aaef..ded01c7d86749178f1ecda1bb7a852df5ca5c0e8 100644
--- a/dumux/porousmediumflow/velocityoutput.hh
+++ b/dumux/porousmediumflow/velocityoutput.hh
@@ -239,7 +239,7 @@ public:
                 // calculate the subcontrolvolume velocity by the Piola transformation
                 Dune::FieldVector<CoordScalar, dimWorld> scvVelocity(0);
 
-                jacobianT2.mtv(scvVelocities[scv.indexInElement()], scvVelocity);
+                jacobianT2.mtv(scvVelocities[scv.localDofIndex()], scvVelocity);
                 scvVelocity /= geometry.integrationElement(localPos)*cellNum_[vIdxGlobal];
                 // add up the wetting phase subcontrolvolume velocities for each vertex
                 velocity[vIdxGlobal] += scvVelocity;
diff --git a/dumux/porousmediumflow/volumevariables.hh b/dumux/porousmediumflow/volumevariables.hh
index 589b590c118c32054a94b6613430deaf21837560..eb7961aa8231518a5a2e6f001b75cfb36a3fb010 100644
--- a/dumux/porousmediumflow/volumevariables.hh
+++ b/dumux/porousmediumflow/volumevariables.hh
@@ -96,7 +96,7 @@ public:
     template<class ElemSol, class Scv>
     static const PrimaryVariables& extractDofPriVars(const ElemSol& elemSol,
                                                      const Scv& scv)
-    { return elemSol[scv.indexInElement()]; }
+    { return elemSol[scv.localDofIndex()]; }
 
     /*!
      * \brief Return the vector of primary variables
diff --git a/test/discretization/box/test_boxfvgeometry.cc b/test/discretization/box/test_boxfvgeometry.cc
index 570b78a8f2c17dadf8775755f3b85d2f1161dac4..dca65f1e850b3e6ee36833b8106ee6407e79693e 100644
--- a/test/discretization/box/test_boxfvgeometry.cc
+++ b/test/discretization/box/test_boxfvgeometry.cc
@@ -91,7 +91,7 @@ int main (int argc, char *argv[]) try
 
         for (auto&& scv : scvs(fvGeometry))
         {
-            std::cout << "-- scv " << scv.indexInElement() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
+            std::cout << "-- scv " << scv.localDofIndex() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
         }
 
         auto range2 = scvfs(fvGeometry);
diff --git a/test/discretization/staggered/test_staggeredfvgeometry.cc b/test/discretization/staggered/test_staggeredfvgeometry.cc
index b23013d12abe9f97ef1f68399bfd703b50788c04..8140c0f1069531306fecbd543988475443805676 100644
--- a/test/discretization/staggered/test_staggeredfvgeometry.cc
+++ b/test/discretization/staggered/test_staggeredfvgeometry.cc
@@ -125,7 +125,7 @@ int main (int argc, char *argv[]) try
 
         for (auto&& scv : scvs(fvGeometry))
         {
-            std::cout << "-- scv " << scv.indexInElement() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
+            std::cout << "-- scv " << scv.localDofIndex() << " center at: " << scv.center() << " , volume: " << scv.volume()  << std::endl;
         }
 
         auto range2 = scvfs(fvGeometry);
diff --git a/test/porousmediumflow/1p/implicit/1pniconvectionproblem.hh b/test/porousmediumflow/1p/implicit/1pniconvectionproblem.hh
index 9dcce949d475a0dc3e45dae836e6c94d79f1d6ec..f8237a8b0720fe791734377cf38027fae0b1ae5d 100644
--- a/test/porousmediumflow/1p/implicit/1pniconvectionproblem.hh
+++ b/test/porousmediumflow/1p/implicit/1pniconvectionproblem.hh
@@ -270,9 +270,9 @@ public:
      * Negative values indicate an inflow.
      */
     NumEqVector neumann(const Element& element,
-                             const FVElementGeometry& fvGeometry,
-                             const ElementVolumeVariables& elemVolvars,
-                             const SubControlVolumeFace& scvf) const
+                        const FVElementGeometry& fvGeometry,
+                        const ElementVolumeVariables& elemVolvars,
+                        const SubControlVolumeFace& scvf) const
     {
         NumEqVector values(0.0);
         const auto globalPos = scvf.ipGlobal();
diff --git a/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh b/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh
index addb94784e82ec2ce444228f63c26dbe013ad079..0392a229c6063d87013776b2574bc747435e691c 100644
--- a/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh
+++ b/test/porousmediumflow/1pnc/implicit/1p2ctestproblem.hh
@@ -256,7 +256,7 @@ public:
             if(isBox)
                 for(auto&& scvf : scvfs(fvGeometry))
                     if(scvf.center()[0] > this->fvGridGeometry().bBoxMax()[0] - eps_)
-                        sol[fvGeometry.scv(scvf.insideScvIdx()).indexInElement()][pressureIdx] = dirichletPressure;
+                        sol[fvGeometry.scv(scvf.insideScvIdx()).localDofIndex()][pressureIdx] = dirichletPressure;
 
             return sol;
         }();