diff --git a/dune/gfe/riemannianpnsolver.cc b/dune/gfe/riemannianpnsolver.cc
index 76ad38a2c0229c7cddcf9661672f925a0066ec95..5707d3f8ba5ce204c2e39869fc25274355288d14 100644
--- a/dune/gfe/riemannianpnsolver.cc
+++ b/dune/gfe/riemannianpnsolver.cc
@@ -35,6 +35,13 @@ void RiemannianProximalNewtonSolver<Basis, TargetSpace, Assembler>::
         const Dune::BitSetVector<blocksize>& dirichletNodes,
         const Dune::ParameterTree& parameterSet)
 {
+    if(parameterSet.get("norm", "infinity") == "infinity")
+        normType_ = ErrorNormType::infinity;
+    else if(parameterSet.get("norm", "infinity") == "H1-Semi")
+        normType_ = ErrorNormType::H1semi;
+    else
+        DUNE_THROW(Dune::Exception, "Unknown norm type for stopping criterion!");
+
     setup(grid,
           assembler,
           x,
@@ -314,11 +321,18 @@ void RiemannianProximalNewtonSolver<Basis,TargetSpace,Assembler>::solve()
         corr = corr_global;
 #endif
 
-        // Make infinity norm of corr_global known on all processors
-        double corrNorm = corr.infinity_norm();
-        double corrGlobalInfinityNorm = grid_->comm().max(corrNorm);
+        double corrNorm;
+        if (normType_ == ErrorNormType::infinity)
+            corrNorm = corr.infinity_norm();
+        else if (normType_ == ErrorNormType::H1semi)
+            corrNorm = h1SemiNorm_->operator()(corr);
+        else
+            DUNE_THROW(Dune::Exception, "Unknown norm type for stopping criterion!");
+
+        // Make norm of corr_global known on all processors
+        double corrGlobalNorm = grid_->comm().max(corrNorm);
 
-        if (std::isnan(corrGlobalInfinityNorm))
+        if (std::isnan(corrGlobalNorm))
             solved = false;
 
         if (instrumented_) {
@@ -395,9 +409,14 @@ void RiemannianProximalNewtonSolver<Basis,TargetSpace,Assembler>::solve()
 
         if (solved) {
             if (this->verbosity_ == NumProc::FULL && rank==0)
-                std::cout << "Infinity norm of the correction: " << corrGlobalInfinityNorm << std::endl;
-
-            if (corrGlobalInfinityNorm < this->tolerance_ && corrGlobalInfinityNorm < 1/regularization) {
+                if (normType_ == ErrorNormType::infinity)
+                    std::cout << "infinity norm of the correction: " << corrGlobalNorm << std::endl;
+                else if (normType_ == ErrorNormType::H1semi)
+                    std::cout << "H1-semi norm of the correction: " << corrGlobalNorm << std::endl;
+                else
+                    DUNE_THROW(Dune::Exception, "Unknown norm type for stopping criterion!");
+
+            if (corrGlobalNorm < this->tolerance_ && corrGlobalNorm < 1/regularization) {
                 if (this->verbosity_ == NumProc::FULL and rank==0)
                     std::cout << "CORRECTION IS SMALL ENOUGH" << std::endl;
 
diff --git a/dune/gfe/riemannianpnsolver.hh b/dune/gfe/riemannianpnsolver.hh
index 915de346baec9100880cd1a1fe31c64904c86dbb..67aa29a0c62d686722b54ec929867c4897ebad0f 100644
--- a/dune/gfe/riemannianpnsolver.hh
+++ b/dune/gfe/riemannianpnsolver.hh
@@ -161,6 +161,9 @@ protected:
     /** \brief If set to true we log convergence speed and other stuff */
     bool instrumented_;
 
+    /** \brief Norm type used for stopping criterion (default is infinity norm) */
+    enum class ErrorNormType {infinity, H1semi} normType_ = ErrorNormType::infinity;
+
     /** \brief Store information about solver runs for unit testing */
     Statistics statistics_;
 
diff --git a/dune/gfe/riemanniantrsolver.cc b/dune/gfe/riemanniantrsolver.cc
index 7678a442dfc09472ce7c813e68b987e23b394bf1..5285d01c92e621bb0f7bde4c6085bd012c0b1309 100644
--- a/dune/gfe/riemanniantrsolver.cc
+++ b/dune/gfe/riemanniantrsolver.cc
@@ -36,6 +36,13 @@ void RiemannianTrustRegionSolver<Basis, TargetSpace, Assembler>::
         const Dune::BitSetVector<blocksize>& dirichletNodes,
         const Dune::ParameterTree& parameterSet)
 {
+    if(parameterSet.get("norm", "infinity") == "infinity")
+        normType_ = ErrorNormType::infinity;
+    else if(parameterSet.get("norm", "infinity") == "H1-Semi")
+        normType_ = ErrorNormType::H1semi;
+    else
+        DUNE_THROW(Dune::Exception, "Unknown norm type for stopping criterion!");
+
     setup(grid,
           assembler,
           x,
@@ -505,9 +512,16 @@ void RiemannianTrustRegionSolver<Basis,TargetSpace,Assembler>::solve()
         corr = corr_global;
 #endif
 
-        // Make infinity norm of corr_global known on all processors
-        double corrNorm = corr.infinity_norm();
-        double corrGlobalInfinityNorm = grid_->comm().max(corrNorm);
+        double corrNorm;
+        if (normType_ == ErrorNormType::infinity)
+            corrNorm = corr.infinity_norm();
+        else if (normType_ == ErrorNormType::H1semi)
+            corrNorm = h1SemiNorm_->operator()(corr);
+        else
+            DUNE_THROW(Dune::Exception, "Unknown norm type for stopping criterion!");
+
+        // Make norm of corr_global known on all processors
+        double corrGlobalNorm = grid_->comm().max(corrNorm);
 
         if (instrumented_) {
 
@@ -584,9 +598,14 @@ void RiemannianTrustRegionSolver<Basis,TargetSpace,Assembler>::solve()
 
         if (solved) {
             if (this->verbosity_ == NumProc::FULL and rank==0)
-                std::cout << "Infinity norm of the correction: " << corrGlobalInfinityNorm << std::endl;
-
-            if (corrGlobalInfinityNorm < this->tolerance_ && corrGlobalInfinityNorm < trustRegion.radius()*smallestScalingParameter) {
+                if (normType_ == ErrorNormType::infinity)
+                    std::cout << "infinity norm of the correction: " << corrGlobalNorm << std::endl;
+                else if (normType_ == ErrorNormType::H1semi)
+                    std::cout << "H1-semi norm of the correction: " << corrGlobalNorm << std::endl;
+                else
+                    DUNE_THROW(Dune::Exception, "Unknown norm type for stopping criterion!");
+
+            if (corrGlobalNorm < this->tolerance_ && corrGlobalNorm < trustRegion.radius()*smallestScalingParameter) {
                 if (this->verbosity_ == NumProc::FULL and rank==0)
                     std::cout << "CORRECTION IS SMALL ENOUGH" << std::endl;
 
diff --git a/dune/gfe/riemanniantrsolver.hh b/dune/gfe/riemanniantrsolver.hh
index 61df18e603de00f93cd808428ee37ab81c08910c..61017f956384efd1babde011c4e117895f451cc8 100644
--- a/dune/gfe/riemanniantrsolver.hh
+++ b/dune/gfe/riemanniantrsolver.hh
@@ -206,6 +206,9 @@ protected:
     /** \brief If set to true we log convergence speed and other stuff */
     bool instrumented_;
 
+    /** \brief Norm type used for stopping criterion (default is infinity norm) */
+    enum class ErrorNormType {infinity, H1semi} normType_ = ErrorNormType::infinity;
+
     /** \brief Store information about solver runs for unit testing */
     Statistics statistics_;