diff --git a/CHANGELOG.md b/CHANGELOG.md
index 7e34a7e9d26ca74f5b9fffe3aef5513a98a92367..328c900e9944aba53dba918ef81655b2210191dc 100644
--- a/CHANGELOG.md
+++ b/CHANGELOG.md
@@ -21,7 +21,10 @@ Differences Between DuMuX 3.1 and DuMuX 3.0
   using NewtonConvergenceWriter = Dumux::NewtonConvergenceWriter<FVGridGeometry, SolutionVector>;
   auto convergenceWriter = std::make_shared<NewtonConvergenceWriter>(*fvGridGeometry);
   ```
-
+- The interface of the abstract `TimeLoopBase` class has been extended by `virtual void advanceTimeStep()`, `virtual void setTimeStepSize(Scalar dt)`,
+  `virtual Scalar maxTimeStepSize()`, and `virtual bool finished()`, thus forcing the inheriting classes to implement those functions.
+  `TimeLoop` is no longer a template argument in `NewtonSolver`'s `solve()` method, thereby increasing
+  type safety and providing better compiler error messages.
 
 ### Deprecated classes/files, to be removed after 3.1:
 
diff --git a/dumux/common/timeloop.hh b/dumux/common/timeloop.hh
index d934a45d417fd09c54569d6832e7a126802bd900..7d1529d7d2f3d23be0a34a7c14e2c32e0239409a 100644
--- a/dumux/common/timeloop.hh
+++ b/dumux/common/timeloop.hh
@@ -75,11 +75,30 @@ public:
     virtual Scalar time() const = 0;
 
     /*!
-     * \brief Returns the suggested time step length \f$\mathrm{[s]}\f$ so that we
-     *        don't miss the beginning of the next episode or cross
-     *        the end of the simulation.
+     * \brief Returns the suggested time step length \f$\mathrm{[s]}\f$
      */
     virtual Scalar timeStepSize() const = 0;
+
+    /*!
+     * \brief Get the maximum possible time step size \f$\mathrm{[s]}\f$
+     */
+    virtual Scalar maxTimeStepSize() const = 0;
+
+    /*!
+     * \brief Advance to the next time step.
+     */
+    virtual void advanceTimeStep() = 0;
+
+    /*!
+     * \brief Set the current time step size to a given value.
+     * \param dt The new value for the time step size \f$\mathrm{[s]}\f$
+     */
+    virtual void setTimeStepSize(Scalar dt) = 0;
+
+    /*!
+     * \brief Returns true if the simulation is finished.
+     */
+    virtual bool finished() const = 0;
 };
 
 /*!
@@ -155,7 +174,7 @@ public:
     /*!
      * \brief Advance time step.
      */
-    void advanceTimeStep()
+    void advanceTimeStep() override
     {
         timeStepIdx_++;
         time_ += timeStepSize_;
@@ -193,7 +212,7 @@ public:
      * To get the time after the time integration you have to add timeStepSize() to
      * time().
      */
-    Scalar time() const override
+    Scalar time() const final
     { return time_; }
 
     /*!
@@ -227,15 +246,13 @@ public:
      * episode, the timeStep() method will take care that the step
      * size won't exceed the episode or the end of the simulation,
      * though.
-     * \note Always call this after TimeLoop::advanceTimeStep() and TimeLoop::reportTimeStep()
      *
      * \param dt The new value for the time step size \f$\mathrm{[s]}\f$
      */
-    void setTimeStepSize(Scalar dt)
+    void setTimeStepSize(Scalar dt) final
     {
         using std::min;
-        computeMaxTimeStepSize_();
-        timeStepSize_ = min(dt, maxTimeStepSize_);
+        timeStepSize_ = min(dt, maxTimeStepSize());
     }
 
     /*!
@@ -246,10 +263,8 @@ public:
      */
     void setMaxTimeStepSize(Scalar maxDt)
     {
-        using std::min;
         userSetMaxTimeStepSize_ = maxDt;
-        computeMaxTimeStepSize_();
-        timeStepSize_ = min(timeStepSize_, maxTimeStepSize_);
+        setTimeStepSize(timeStepSize_);
     }
 
     /*!
@@ -257,7 +272,7 @@ public:
      *        don't miss the beginning of the next episode or cross
      *        the end of the simulation.
      */
-    Scalar timeStepSize() const override
+    Scalar timeStepSize() const final
     { return timeStepSize_; }
 
     /*!
@@ -283,7 +298,7 @@ public:
      * This is the case if either setFinished(true) has been called or
      * if the end time is reached.
      */
-    bool finished() const
+    bool finished() const override
     {
         return finished_ || endTime_-time_ < 1e-10*time_;
     }
@@ -302,10 +317,13 @@ public:
      * \note This gets aligned on every setTimeStepSize call to end time
      *       and other possible check points
      */
-    Scalar maxTimeStepSize() const
+    Scalar maxTimeStepSize() const override
     {
-        using std::min;
-        return min(userSetMaxTimeStepSize_, maxTimeStepSize_);
+        if (finished())
+            return 0.0;
+
+        using std::min; using std::max;
+        return min(userSetMaxTimeStepSize_, max<Scalar>(0.0, endTime_ - time_));
     }
 
     /*!
@@ -360,28 +378,14 @@ public:
      * @}
      */
 
-private:
-    //! Computes the maximum timestep size respecting end time
-    void computeMaxTimeStepSize_()
-    {
-        if (finished())
-        {
-            maxTimeStepSize_ = 0.0;
-            return;
-        }
-
-        using std::max;
-        using std::min;
-        maxTimeStepSize_ = min(userSetMaxTimeStepSize_, max<Scalar>(0.0, endTime_ - time_));
-    }
-
+protected:
     Dune::Timer timer_;
     Scalar time_;
     Scalar endTime_;
 
     Scalar timeStepSize_;
     Scalar lastTimeStepSize_;
-    Scalar maxTimeStepSize_, userSetMaxTimeStepSize_;
+    Scalar userSetMaxTimeStepSize_;
     Scalar timeAfterLastTimeStep_, timeStepWallClockTime_;
     int timeStepIdx_;
     bool finished_;
@@ -408,21 +412,21 @@ public:
     /*!
      * \brief Advance time step.
      */
-    void advanceTimeStep()
+    void advanceTimeStep() override
     {
-        // advance time index and time
-        TimeLoop<Scalar>::advanceTimeStep();
+        const auto dt = this->timeStepSize();
+        const auto nextTime = this->time()+dt;
 
         //! Check point management, TimeLoop::isCheckPoint() has to be called after this!
         // if we reached a periodic check point
-        if (periodicCheckPoints_ && Dune::FloatCmp::eq(this->time(), lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_, 1e-8*this->timeStepSize()))
+        if (periodicCheckPoints_ && Dune::FloatCmp::eq(nextTime, lastPeriodicCheckPoint_ + deltaPeriodicCheckPoint_, 1e-8*dt))
         {
             lastPeriodicCheckPoint_ += deltaPeriodicCheckPoint_;
             isCheckPoint_ = true;
         }
 
         // or a manually set check point
-        else if (!checkPoints_.empty() && Dune::FloatCmp::eq(this->time(), checkPoints_.front(), 1e-8*this->timeStepSize()))
+        else if (!checkPoints_.empty() && Dune::FloatCmp::eq(nextTime, checkPoints_.front(), 1e-8*dt))
         {
             checkPoints_.pop();
             isCheckPoint_ = true;
@@ -434,25 +438,21 @@ public:
             isCheckPoint_ = false;
         }
 
-        // make sure to respect future check check points
-        this->setTimeStepSize(this->timeStepSize());
+        // advance the time step like in the parent class
+        TimeLoop<Scalar>::advanceTimeStep();
     }
 
     /*!
-     * \brief Set the current time step size to a given value.
-     *
-     * If the step size would exceed the length of the current
-     * episode, the timeStep() method will take care that the step
-     * size won't exceed the episode or the end of the simulation,
-     * though.
-     * \note Always call this after TimeLoop::advanceTimeStep()
-     *
-     * \param dt The new value for the time step size \f$\mathrm{[s]}\f$
+     * \brief The current maximum time step size
+     * \note This gets aligned on every setTimeStepSize call to end time
+     *       and other possible check points
      */
-    void setTimeStepSize(Scalar dt)
+    Scalar maxTimeStepSize() const override
     {
         using std::min;
-        TimeLoop<Scalar>::setTimeStepSize(min(dt, computeStepSizeRespectingCheckPoints_()));
+        const auto maxCheckPointDt = computeStepSizeRespectingCheckPoints_();
+        const auto maxDtParent = TimeLoop<Scalar>::maxTimeStepSize();
+        return min(maxDtParent, maxCheckPointDt);
     }
 
     /*!
@@ -481,12 +481,12 @@ public:
             std::cout << "Enabled periodic check points every " << interval
                       << " seconds with the next check point at " << lastPeriodicCheckPoint_ + interval << " seconds." << std::endl;
 
-        // make sure we respect this check point on the next time step
-        setTimeStepSize(this->timeStepSize());
-
         // check if the current time point is a check point
         if (Dune::FloatCmp::eq(this->time(), lastPeriodicCheckPoint_, 1e-8*interval))
             isCheckPoint_ = true;
+
+        // make sure we respect this check point on the next time step
+        this->setTimeStepSize(this->timeStepSize());
     }
 
     //! disable periodic check points
@@ -520,7 +520,7 @@ public:
         setCheckPoint_(t);
 
         // make sure we respect this check point on the next time step
-        setTimeStepSize(this->timeStepSize());
+        this->setTimeStepSize(this->timeStepSize());
     }
 
     /*!
@@ -547,7 +547,7 @@ public:
             setCheckPoint_(*first);
 
         // make sure we respect this check point on the next time step
-        setTimeStepSize(this->timeStepSize());
+        this->setTimeStepSize(this->timeStepSize());
     }
 
 private:
@@ -578,7 +578,7 @@ private:
             std::cout << "Set check point at t = " << t << " seconds." << std::endl;
     }
 
-     /*!
+    /*!
      * \brief Aligns dt to the next check point
      */
     Scalar computeStepSizeRespectingCheckPoints_() const
diff --git a/dumux/nonlinear/newtonsolver.hh b/dumux/nonlinear/newtonsolver.hh
index 1838235f998e540e8b780dbc6659eec9632f14a0..9e43590dc94d69e3b02449f1f1500d8e26a3c66d 100644
--- a/dumux/nonlinear/newtonsolver.hh
+++ b/dumux/nonlinear/newtonsolver.hh
@@ -43,6 +43,7 @@
 #include <dumux/common/exceptions.hh>
 #include <dumux/common/typetraits/vector.hh>
 #include <dumux/common/typetraits/isvalid.hh>
+#include <dumux/common/timeloop.hh>
 #include <dumux/linear/linearsolveracceptsmultitypematrix.hh>
 #include <dumux/linear/matrixconverter.hh>
 #include <dumux/assembly/partialreassembler.hh>
@@ -92,6 +93,7 @@ class NewtonSolver
     using JacobianMatrix = typename Assembler::JacobianMatrix;
     using SolutionVector = typename Assembler::ResidualType;
     using ConvergenceWriter = ConvergenceWriterInterface<SolutionVector>;
+    using TimeLoop = TimeLoopBase<Scalar>;
 
     using PrimaryVariableSwitch = typename Detail::GetPVSwitch<Assembler>::type;
     using HasPriVarsSwitch = typename Detail::GetPVSwitch<Assembler>::value_t; // std::true_type or std::false_type
@@ -194,7 +196,6 @@ public:
      * \brief Run the Newton method to solve a non-linear system.
      *        Does time step control when the Newton fails to converge
      */
-    template<class TimeLoop>
     DUNE_DEPRECATED_MSG("Use attachConvergenceWriter(convWriter) and solve(x, *timeLoop) instead")
     void solve(SolutionVector& uCurrentIter, TimeLoop& timeLoop,
                std::shared_ptr<ConvergenceWriter> convWriter)
@@ -204,11 +205,11 @@ public:
 
         solve(uCurrentIter, timeLoop);
     }
+
     /*!
      * \brief Run the Newton method to solve a non-linear system.
      *        Does time step control when the Newton fails to converge
      */
-    template<class TimeLoop>
     void solve(SolutionVector& uCurrentIter, TimeLoop& timeLoop)
     {
         if (assembler_->isStationaryProblem())