diff --git a/dune/vtk/pvdwriter.hh b/dune/vtk/pvdwriter.hh
index c5345b44ff3c08264efff1b9f0c1275dbda2cf24..919cab9b54d417a392cf197dcef053f336a5bb69 100644
--- a/dune/vtk/pvdwriter.hh
+++ b/dune/vtk/pvdwriter.hh
@@ -6,13 +6,14 @@
 #include <tuple>
 
 #include <dune/vtk/vtktypes.hh>
-// #include <dune/vtk/vtkwriterinterface.hh>
+#include <dune/vtk/filewriter.hh>
 
 namespace Dune
 {
   /// File-Writer for ParaView .pvd files
   template <class VtkWriter>
   class PvdWriter
+      : public FileWriter
   {
     using Self = PvdWriter;
 
@@ -29,7 +30,18 @@ namespace Dune
     }
 
     /// Write the attached data to the file
-    void writeTimestep (double time, std::string const& fn) const;
+    /**
+     * Create timestep files for the data associated to the current timestep `time`.
+     *
+     * \param time  The time value of the written data
+     * \param fn  Filename of the file to write to. Is stored in \ref timesteps_.
+     * \param writeCollection  Create a collection .pvd file directly
+     **/
+    void writeTimestep (double time, std::string const& fn, bool writeCollection = true) const;
+
+    /// Writes collection of timesteps to .pvd file.
+    // NOTE: requires an aforging call to \ref writeTimestep
+    virtual void write (std::string const& fn) const override;
 
     /// Attach point data to the writer, \see VtkFunction for possible arguments
     template <class Function, class... Args>
@@ -49,7 +61,7 @@ namespace Dune
 
   protected:
     /// Write a series of vtk files in a .pvd ParaView Data file
-    void writeFile (double time, std::ofstream& out) const;
+    void writeFile (std::ofstream& out) const;
 
   protected:
     VtkWriter vtkWriter_;
diff --git a/dune/vtk/pvdwriter.impl.hh b/dune/vtk/pvdwriter.impl.hh
index 669a894da8c0c6a6c11a46276a495db01c6b95b7..ddafacc260e2757b4c022f8c17374767ae1da2af 100644
--- a/dune/vtk/pvdwriter.impl.hh
+++ b/dune/vtk/pvdwriter.impl.hh
@@ -9,7 +9,7 @@ namespace Dune {
 
 template <class W>
 void PvdWriter<W>
-  ::writeTimestep (double time, std::string const& fn) const
+  ::writeTimestep (double time, std::string const& fn, bool writeCollection) const
 {
   auto p = filesystem::path(fn);
   auto name = p.stem();
@@ -31,6 +31,36 @@ void PvdWriter<W>
   timesteps_.emplace_back(time, filename + ext);
   vtkWriter_.write(filename + ext);
 
+  if (rank == 0 && writeCollection) {
+    std::ofstream out(p.string() + ".pvd", std::ios_base::ate | std::ios::binary);
+    assert(out.is_open());
+
+    out.imbue(std::locale::classic());
+    out << std::setprecision(datatype_ == Vtk::FLOAT32
+      ? std::numeric_limits<float>::digits10+2
+      : std::numeric_limits<double>::digits10+2);
+
+    writeFile(out);
+  }
+}
+
+
+template <class W>
+void PvdWriter<W>
+  ::write (std::string const& fn) const
+{
+  auto p = filesystem::path(fn);
+  auto name = p.stem();
+  p.remove_filename();
+  p /= name.string();
+
+  int rank = 0;
+  int num_ranks = 1;
+#ifdef HAVE_MPI
+  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
+  MPI_Comm_size(MPI_COMM_WORLD, &num_ranks);
+#endif
+
   if (rank == 0) {
     std::ofstream out(p.string() + ".pvd", std::ios_base::ate | std::ios::binary);
     assert(out.is_open());
@@ -40,14 +70,14 @@ void PvdWriter<W>
       ? std::numeric_limits<float>::digits10+2
       : std::numeric_limits<double>::digits10+2);
 
-    writeFile(time, out);
+    writeFile(out);
   }
 }
 
 
 template <class W>
 void PvdWriter<W>
-  ::writeFile (double time, std::ofstream& out) const
+  ::writeFile (std::ofstream& out) const
 {
   out << "<?xml version=\"1.0\"?>\n";
   out << "<VTKFile"
diff --git a/dune/vtk/vtktimeserieswriter.hh b/dune/vtk/vtktimeserieswriter.hh
index 13e8507c0dccdda3a50ce064876afdf253f89cb2..f337ffd97b0c4b862bc095b79bad5b2123e5111e 100644
--- a/dune/vtk/vtktimeserieswriter.hh
+++ b/dune/vtk/vtktimeserieswriter.hh
@@ -47,16 +47,26 @@ namespace Dune
       filesystem::create_directories(tmpDir_);
     }
 
-    ~VtkTimeseriesWriter ()
-    {
-      std::remove(tmpDir_.string().c_str());
-    }
-
-    /// Write the attached data to the file with \ref Vtk::FormatTypes and \ref Vtk::DataTypes
-    void writeTimestep (double time, std::string const& fn) const;
+    /// Remove all written intermediate files and remove temporary directory
+    ~VtkTimeseriesWriter ();
+
+    /// Write the attached data to the file
+    /**
+     * Create intermediate files for the data associated to the current timestep `time`.
+     *
+     * \param time  The time value of the written data
+     * \param fn  Filename of the file to write to. Only the base part
+     *            (without dir and extentsion) is used to write the intermediate
+     *            file into a tmp directory.
+     * \param writeCollection  Create a timeseries file directly
+     **/
+    void writeTimestep (double time, std::string const& fn, bool writeCollection = true) const;
 
     /// Writes all timesteps to single timeseries file.
     // NOTE: requires an aforging call to \ref writeTimestep
+    /**
+     * Create a timeseries file with all timesteps written by \ref writeTimestep.
+     **/
     virtual void write (std::string const& fn) const override;
 
     /// Attach point data to the writer, \see VtkFunction for possible arguments
diff --git a/dune/vtk/vtktimeserieswriter.impl.hh b/dune/vtk/vtktimeserieswriter.impl.hh
index 4e53490c9d495ab167556dd032b19008e6b1d0e4..b51ad6368c026b7db6f1f49bb451a2c69324750f 100644
--- a/dune/vtk/vtktimeserieswriter.impl.hh
+++ b/dune/vtk/vtktimeserieswriter.impl.hh
@@ -22,15 +22,30 @@
 
 namespace Dune {
 
+template <class W>
+VtkTimeseriesWriter<W>::~VtkTimeseriesWriter ()
+{
+  if (initialized_) {
+    int ec = std::remove(filenameMesh_.c_str());
+    assert(ec == 0);
+    for (auto const& timestep : timesteps_) {
+      ec = std::remove(timestep.second.c_str());
+      assert(ec == 0);
+    }
+  }
+  std::remove(tmpDir_.string().c_str());
+}
+
+
 template <class W>
 void VtkTimeseriesWriter<W>
-  ::writeTimestep (double time, std::string const& fn)
+  ::writeTimestep (double time, std::string const& fn, bool writeCollection) const
 {
   auto name = filesystem::path(fn).stem();
   auto tmp = tmpDir_;
   tmp /= name.string();
 
-  vtkWriter_.update();
+  vtkWriter_.dataCollector_.update();
 
   std::string filenameBase = tmp.string();
 
@@ -55,6 +70,9 @@ void VtkTimeseriesWriter<W>
   std::ofstream out(filenameData, std::ios_base::ate | std::ios::binary);
   vtkWriter_.writeDataAppended(out, blocks_);
   timesteps_.emplace_back(time, filenameData);
+
+  if (writeCollection)
+    write(fn);
 }
 
 
@@ -63,7 +81,6 @@ void VtkTimeseriesWriter<W>
   ::write (std::string const& fn) const
 {
   assert( initialized_ );
-  assert( timesteps_.size() < 1000 ); // maximal number of allowed timesteps in timeseries file
 
   auto p = filesystem::path(fn);
   auto name = p.stem();
@@ -109,16 +126,6 @@ void VtkTimeseriesWriter<W>
     vtkWriter_.writeTimeseriesParallelFile(parallel_out, p.string() + "_ts", num_ranks, timesteps_);
   }
 #endif
-
-  // remove all temporary data files
-  int ec = std::remove(filenameMesh_.c_str());
-  assert(ec == 0);
-  for (auto const& timestep : timesteps_) {
-    ec = std::remove(timestep.second.c_str());
-    assert(ec == 0);
-  }
-  timesteps_.clear();
-  initialized_ = false;
 }
 
 } // end namespace Dune
diff --git a/dune/vtk/vtkwriterinterface.hh b/dune/vtk/vtkwriterinterface.hh
index ffc1a5ba0e1e96d5ffb120e8b84704bee597ffba..c32a0625efed470354e3b46d92a1fd4c826b43f0 100644
--- a/dune/vtk/vtkwriterinterface.hh
+++ b/dune/vtk/vtkwriterinterface.hh
@@ -147,11 +147,6 @@ namespace Dune
       return datatype_;
     }
 
-    void update ()
-    {
-      dataCollector_.update();
-    }
-
   protected:
     mutable DataCollector dataCollector_;
 
diff --git a/src/pvdwriter.cc b/src/pvdwriter.cc
index 2c466ed138bf1578077139e5e0c37a51c26337ae..242b6cd159d5ca3c94d1057e4e2011c603c99faf 100644
--- a/src/pvdwriter.cc
+++ b/src/pvdwriter.cc
@@ -41,10 +41,12 @@ void write (std::string prefix, GridView const& gridView)
   using Writer = VtkUnstructuredGridWriter<GridView>;
   PvdWriter<Writer> pvdWriter(gridView, Vtk::ASCII, Vtk::FLOAT32);
 
+  std::string filename = prefix + "_" + std::to_string(GridView::dimensionworld) + "d_ascii.vtu";
+
   pvdWriter.addPointData(p1Analytic, "p1");
   pvdWriter.addCellData(p1Analytic, "p0");
   for (double t = 0.0; t < 10.0; t += 1.0) {
-    pvdWriter.writeTimestep(t, prefix + "_" + std::to_string(GridView::dimensionworld) + "d_ascii.vtu");
+    pvdWriter.writeTimestep(t, filename);
   }
 }
 
diff --git a/src/timeserieswriter.cc b/src/timeserieswriter.cc
index e73314442bb7d13aecb0dafcf5350148d4e9b551..7b54dd41e6753dd7b9b22faf099f9104936f0f06 100644
--- a/src/timeserieswriter.cc
+++ b/src/timeserieswriter.cc
@@ -25,7 +25,8 @@ template <class GridView>
 void write (std::string prefix, GridView const& gridView)
 {
   FieldVector<double,GridView::dimensionworld> c{11.0, 7.0, 3.0};
-  auto p1Analytic = makeAnalyticGridViewFunction([&c](auto const& x) -> float { return c.dot(x); }, gridView);
+  double shift = 0.0;
+  auto p1Analytic = makeAnalyticGridViewFunction([&c,&shift](auto const& x) -> float { return c.dot(x) + shift; }, gridView);
 
   using Writer = VtkUnstructuredGridWriter<GridView>;
   VtkTimeseriesWriter<Writer> seriesWriter(gridView, Vtk::BINARY, Vtk::FLOAT32);
@@ -33,7 +34,8 @@ void write (std::string prefix, GridView const& gridView)
   seriesWriter.addCellData(p1Analytic, "q0");
   std::string filename = prefix + "_" + std::to_string(GridView::dimensionworld) + "d_binary32.vtu";
   for (double t = 0.0; t < 5; t += 0.5) {
-    seriesWriter.writeTimestep(t, filename);
+    seriesWriter.writeTimestep(t, filename, false);
+    shift += 0.25;
   }
   seriesWriter.write(filename);