diff --git a/reco/tasks/CbmReco.cxx b/reco/tasks/CbmReco.cxx
index ccfd00a3167782bea1354549e17e8d5f2de8b7db..80a3e7a0bef0ad9c14e01d4d45752ca3585f5e1b 100644
--- a/reco/tasks/CbmReco.cxx
+++ b/reco/tasks/CbmReco.cxx
@@ -6,6 +6,7 @@
 
 #include "CbmSourceTs.h"
 #include "CbmTaskBuildEvents.h"
+#include "CbmTaskDigiEventQa.h"
 #include "CbmTaskTriggerDigi.h"
 #include "CbmTaskUnpack.h"
 #include "CbmTsEventHeader.h"
@@ -14,6 +15,9 @@
 #include <FairRootFileSink.h>
 #include <FairRunOnline.h>
 
+#include <THttpServer.h>
+#include <TRootSniffer.h>
+
 #include <iostream>
 #include <memory>
 #include <string>
@@ -102,6 +106,10 @@ int32_t CbmReco::Run()
     evtBuild->SetEventWindow(entry.first, entry.second.first, entry.second.second);
   evtBuild->SetOutputBranchPersistent("DigiEvent", fConfig.fStoreEvents);
 
+  // --- Event QA
+  auto evtQa = make_unique<CbmTaskDigiEventQa>();
+  evtQa->Config(fConfig);
+
   // --- Run configuration
   FairRunOnline run(source.release());
   run.SetSink(sink.release());
@@ -109,6 +117,13 @@ int32_t CbmReco::Run()
   run.AddTask(unpack.release());
   run.AddTask(trigger.release());
   run.AddTask(evtBuild.release());
+  run.AddTask(evtQa.release());
+
+  // ----- HttpServer for online monitoring
+  Int_t serverHttpPort    = 8080;
+  Int_t serverRefreshRate = 100;  // timeslices
+  run.ActivateHttpServer(serverRefreshRate, serverHttpPort);
+  run.GetHttpServer()->GetSniffer()->SetScanGlobalDir(kFALSE);
 
   // --- Initialise and start run
   timer.Stop();
diff --git a/reco/tasks/CbmTaskDigiEventQa.cxx b/reco/tasks/CbmTaskDigiEventQa.cxx
index c775574f50cf303637df45dfbe138bb06e387604..fa73bb7fe21f2339af78e276615e56c35a352e1a 100644
--- a/reco/tasks/CbmTaskDigiEventQa.cxx
+++ b/reco/tasks/CbmTaskDigiEventQa.cxx
@@ -8,6 +8,7 @@
 #include "CbmDigiManager.h"
 #include "CbmDigiTimeslice.h"
 #include "CbmModuleList.h"
+#include "CbmReco.h"  // for CbmRecoConfig
 
 #include <FairLogger.h>
 #include <FairRunOnline.h>
@@ -32,6 +33,19 @@ CbmTaskDigiEventQa::~CbmTaskDigiEventQa() {}
 // ---------------------------------------------------------------------------
 
 
+// -----   Configuration   ---------------------------------------------------
+void CbmTaskDigiEventQa::Config(const CbmRecoConfig& config)
+{
+  auto limitIt = config.fEvtbuildWindows.find(ECbmModuleId::kSts);
+  if (limitIt != config.fEvtbuildWindows.end()) {
+    double lower     = limitIt->second.first - 10.;
+    double upper     = limitIt->second.second + 10.;
+    fHistDigiTimeSts = new TH1F("hDigiTimeSts", "STS digi time in event", 100, lower, upper);
+  }
+}
+// ---------------------------------------------------------------------------
+
+
 // -----   Execution   -------------------------------------------------------
 void CbmTaskDigiEventQa::Exec(Option_t*)
 {
@@ -55,7 +69,7 @@ void CbmTaskDigiEventQa::Exec(Option_t*)
     for (const auto& digi : event.fData.fSts.fDigis) {
       numDigis++;
       const double tDigi = digi.GetTime() - event.fTime;
-      fHistDigiTime->Fill(tDigi);
+      if (fHistDigiTimeSts) fHistDigiTimeSts->Fill(tDigi);
       sumT += tDigi;
       sumTsq += tDigi * tDigi;
     }
@@ -102,10 +116,11 @@ void CbmTaskDigiEventQa::Finish()
   LOG(info) << "Events     : " << fNumEvents;
   LOG(info) << "Digis      : " << fNumDigis;
   LOG(info) << "Exec time  : " << fixed << setprecision(2) << 1000. * fExecTime / double(fNumTs) << " ms / TS";
-  LOG(info) << "Digi times : entries " << fHistDigiTime->GetEntries() << ", mean " << fHistDigiTime->GetMean()
-            << ", stddev " << fHistDigiTime->GetStdDev();
+  if (fHistDigiTimeSts)
+    LOG(info) << "STS digi times : entries " << fHistDigiTimeSts->GetEntries() << ", mean "
+              << fHistDigiTimeSts->GetMean() << ", stddev " << fHistDigiTimeSts->GetStdDev();
   LOG(info) << "=====================================";
-  fHistDigiTime->Write();
+  if (fHistDigiTimeSts) fHistDigiTimeSts->Write();
 }
 // ----------------------------------------------------------------------------
 
@@ -129,10 +144,12 @@ InitStatus CbmTaskDigiEventQa::Init()
   }
   LOG(info) << "--- Found branch DigiEvent";
 
-  // --- Histograms
-  fHistDigiTime       = new TH1F("hDigiTime", "Digi time in event", 100, -50., 50.);
+  // --- Register histograms
   THttpServer* server = FairRunOnline::Instance()->GetHttpServer();
-  if (server) server->Register("DigiTime", fHistDigiTime);
+  if (server) {
+    LOG(info) << "--- Http server present; registering histograms";
+    if (fHistDigiTimeSts) server->Register("DigiEvent", fHistDigiTimeSts);
+  }
 
   LOG(info) << "==================================================";
   std::cout << std::endl;
@@ -141,4 +158,5 @@ InitStatus CbmTaskDigiEventQa::Init()
 }
 // ----------------------------------------------------------------------------
 
+
 ClassImp(CbmTaskDigiEventQa)
diff --git a/reco/tasks/CbmTaskDigiEventQa.h b/reco/tasks/CbmTaskDigiEventQa.h
index 6a5d60ebf8128e3a334625391a96f7c85cb748d8..32452a9133ed177375ee9d1a39584bf67ebc9277 100644
--- a/reco/tasks/CbmTaskDigiEventQa.h
+++ b/reco/tasks/CbmTaskDigiEventQa.h
@@ -12,15 +12,20 @@
 
 #include <vector>
 
-#include "EventBuilder.h"
+//#include "EventBuilder.h"
 
 class CbmDigiManager;
+class CbmRecoConfig;
 class TH1F;
 
 /** @class CbmTaskDigiEventQa
  ** @brief QA task class for digi events produced by the event builder
  ** @author Volker Friese <v.friese@gsi.de>
  ** @since 15.03.2022
+ **
+ ** Currently implemented functionality: histogram the STS digi time within each event.
+ ** To be expanded for other detector systems and probably more QA figures.
+ ** The histograms are published to the THttpServer.
  **/
 class CbmTaskDigiEventQa : public FairTask {
 
@@ -37,6 +42,14 @@ public:
   virtual ~CbmTaskDigiEventQa();
 
 
+  /** @brief Configuration
+   ** @param config  Reconstruction configuration
+   **
+   ** Histograms are created with limits adjusted to the windows use by the event builder.
+   **/
+  void Config(const CbmRecoConfig& config);
+
+
   /** @brief Task execution **/
   virtual void Exec(Option_t* opt);
 
@@ -60,7 +73,7 @@ private:                                               // members
   size_t fNumEvents                        = 0;        ///< Number of analysed events
   size_t fNumDigis                         = 0;        ///< Number of analysed digis
   double fExecTime                         = 0.;       ///< Execution time [s]
-  TH1F* fHistDigiTime                      = nullptr;  //!
+  TH1F* fHistDigiTimeSts                   = nullptr;  //!
 
 
   ClassDef(CbmTaskDigiEventQa, 1);