diff --git a/reco/eventbuilder/digis/CbmSeedFinderQa.cxx b/reco/eventbuilder/digis/CbmSeedFinderQa.cxx
index f90037f39f7b19575ec37ed7667105737c935239..2ccce441faaecd1b73deed4e0ebed1d1c1ea21c8 100644
--- a/reco/eventbuilder/digis/CbmSeedFinderQa.cxx
+++ b/reco/eventbuilder/digis/CbmSeedFinderQa.cxx
@@ -23,6 +23,12 @@ CbmSeedFinderQa::CbmSeedFinderQa() : fOutFolder("SeedFinderQA", "Seed finder QA"
     new TH1F("fhLinkedMCEventsPerTrigger", "Linked MC events per trigger (=0 for pure noise)", 5, -1, 4);
   fhLinkedMCEventsPerTrigger->SetCanExtend(TH1::kAllAxes);
 
+  fhLinkedTriggersPerMCEvent = new TH1F("fhLinkedTriggersPerMCEvent", "Linked triggers per MC event", 5, -1, 4);
+  fhLinkedTriggersPerMCEvent->SetCanExtend(TH1::kAllAxes);
+
+  fhMatchedTriggersPerMCEvent = new TH1F("fhMatchedTriggersPerMCEvent", "Matched triggers per MC event", 5, -1, 4);
+  fhMatchedTriggersPerMCEvent->SetCanExtend(TH1::kAllAxes);
+
   fhCorrectDigiRatio = new TH1F("fhCorrectDigiRatio", "Correct digis per seed [pct]", 416, -2, 102);
   fhCorrectDigiRatioNoNoise =
     new TH1F("fhCorrectDigiRatioNoNoise", "Correct digis per seed [pct], disregarding noise", 416, -2, 102);
@@ -44,9 +50,11 @@ CbmSeedFinderQa::CbmSeedFinderQa() : fOutFolder("SeedFinderQA", "Seed finder QA"
   histFolder->Add(fhCorrectVsFoundNoNoise);
   histFolder->Add(fhTimeOffset);
   histFolder->Add(fhLinkedMCEventsPerTrigger);
+  histFolder->Add(fhLinkedTriggersPerMCEvent);
+  histFolder->Add(fhMatchedTriggersPerMCEvent);
 
-  fCanv = new CbmQaCanvas("cSummary", "", 4 * 400, 2 * 400);
-  fCanv->Divide2D(8);
+  fCanv = new CbmQaCanvas("cSummary", "", 4 * 400, 3 * 400);
+  fCanv->Divide2D(10);
   fOutFolder.Add(fCanv);
 }
 
@@ -60,6 +68,8 @@ CbmSeedFinderQa::~CbmSeedFinderQa()
   delete fhCorrectVsFoundNoNoise;
   delete fhTimeOffset;
   delete fhLinkedMCEventsPerTrigger;
+  delete fhLinkedTriggersPerMCEvent;
+  delete fhMatchedTriggersPerMCEvent;
   delete fCanv;
 }
 
@@ -72,8 +82,10 @@ void CbmSeedFinderQa::Init()
   fEventList = (CbmMCEventList*) FairRootManager::Instance()->GetObject("MCEventList.");
 }
 
-void CbmSeedFinderQa::FillQaInfo(const int32_t WinStart, const int32_t WinEnd, const std::vector<CbmMatch>* vDigiMatch,
-                                 const double seedTime)
+void CbmSeedFinderQa::ResetPerTsStorage() { fvEventMatchesPerTs.clear(); }
+
+void CbmSeedFinderQa::FillQaSeedInfo(const int32_t WinStart, const int32_t WinEnd,
+                                     const std::vector<CbmMatch>* vDigiMatch, const double seedTime)
 {
   fvSeedTimesFull.push_back(seedTime);
 
@@ -117,6 +129,8 @@ void CbmSeedFinderQa::FillQaInfo(const int32_t WinStart, const int32_t WinEnd, c
   {
     seedMatch.AddLink(1.0, 0, -1, 0);
     fvEventMatches.push_back(seedMatch);
+    fvEventMatchesPerTs.push_back(seedMatch);
+
     //should never show up in histograms
     fvCorrectDigiRatio.push_back(std::numeric_limits<double>::quiet_NaN());
     fvCorrectDigiRatioNoNoise.push_back(std::numeric_limits<double>::quiet_NaN());
@@ -126,6 +140,7 @@ void CbmSeedFinderQa::FillQaInfo(const int32_t WinStart, const int32_t WinEnd, c
   }
   else {
     fvEventMatches.push_back(seedMatch);
+    fvEventMatchesPerTs.push_back(seedMatch);
   }
 
   //correct digis in seed window
@@ -160,6 +175,33 @@ void CbmSeedFinderQa::FillQaInfo(const int32_t WinStart, const int32_t WinEnd, c
   fvTimeOffset.push_back(timeDiff);
 }
 
+void CbmSeedFinderQa::FillQaMCInfo()
+{
+  const uint32_t nEvents = fEventList->GetNofEvents();
+  std::vector<uint32_t> vLinkedTriggersPerMCEvent;
+  std::vector<uint32_t> vMatchedTriggersPerMCEvent;
+  vLinkedTriggersPerMCEvent.resize(nEvents, 0);
+  vMatchedTriggersPerMCEvent.resize(nEvents, 0);
+
+  for (uint32_t iSeed = 0; iSeed < fvEventMatchesPerTs.size(); iSeed++) {
+
+    const CbmMatch eventMatch = fvEventMatchesPerTs.at(iSeed);
+    for (int32_t iLink = 0; iLink < eventMatch.GetNofLinks(); iLink++) {
+      const CbmLink eventLink = eventMatch.GetLink(iLink);
+      vLinkedTriggersPerMCEvent[fEventList->GetEventIndex(eventLink)]++;
+    }
+    const CbmLink matchedLink = eventMatch.GetMatchedLink();
+    vMatchedTriggersPerMCEvent[fEventList->GetEventIndex(matchedLink)]++;
+  }
+
+  for (const auto& value : vLinkedTriggersPerMCEvent) {
+    fhLinkedTriggersPerMCEvent->Fill(value);
+  }
+  for (const auto& value : vMatchedTriggersPerMCEvent) {
+    fhMatchedTriggersPerMCEvent->Fill(value);
+  }
+}
+
 void CbmSeedFinderQa::FillHistos()
 {
   for (uint32_t iEvent = 0; iEvent < fvEventMatches.size(); iEvent++) {
@@ -238,6 +280,12 @@ void CbmSeedFinderQa::WriteHistos()
   fCanv->cd(8);
   fhLinkedMCEventsPerTrigger->DrawCopy("colz", "");
 
+  fCanv->cd(9);
+  fhLinkedTriggersPerMCEvent->DrawCopy("colz", "");
+
+  fCanv->cd(10);
+  fhMatchedTriggersPerMCEvent->DrawCopy("colz", "");
+
   FairSink* sink = FairRootManager::Instance()->GetSink();
   sink->WriteObject(&fOutFolder, nullptr);
 }
diff --git a/reco/eventbuilder/digis/CbmSeedFinderQa.h b/reco/eventbuilder/digis/CbmSeedFinderQa.h
index 2722d32b6a1b99801bcca84c60852c746e1ce750..ed0d06c29032c76d8b347fc22af3384ff1019110 100644
--- a/reco/eventbuilder/digis/CbmSeedFinderQa.h
+++ b/reco/eventbuilder/digis/CbmSeedFinderQa.h
@@ -61,18 +61,27 @@ public:
   /** @brief Difference between true event time and seed time. */
   std::vector<double> fvTimeOffset;
 
+  /** @brief Matches that link constructed event seeds to MC events, current timeslice only. */
+  std::vector<CbmMatch> fvEventMatchesPerTs;
+
   /** @brief Gather QA Information. 
   * @params WinStart Starting position of seed window.
   * @params WinStart End position of seed window.
   * @params vDigiMatch Input vector of digi matches (should match input data to MC events). 
   * @params seedTime Current seed. 
   **/
-  void FillQaInfo(const int32_t WinStart, const int32_t WinEnd, const std::vector<CbmMatch>* vDigiMatch,
-                  const double seedTime);
+  void FillQaSeedInfo(const int32_t WinStart, const int32_t WinEnd, const std::vector<CbmMatch>* vDigiMatch,
+                      const double seedTime);
 
   /** @brief Output QA Information. */
   void OutputQa();
 
+  /** @brief Fill QA Information that uses the full list of MC events per TS. */
+  void FillQaMCInfo();
+
+  /** @brief Reset containers that are persistent for one TS. */
+  void ResetPerTsStorage();
+
   /** @brief Finalize histograms and canvases and write to file. */
   void WriteHistos();
 
@@ -84,14 +93,16 @@ private:
   TFolder fOutFolder;             /// output folder with histos and canvases
 
   /** Histograms **/
-  TH1F* fhLinkedMCEventsPerTrigger = nullptr;  /// linked MC events per trigger
-  TH1F* fhCorrectDigiRatio         = nullptr;  /// correct digis per event
-  TH1F* fhCorrectDigiRatioNoNoise  = nullptr;  /// correct digis per event, disregarding noise
-  TH1F* fhNoiseDigiRatio           = nullptr;  /// noise digis per event
-  TH1F* fhFoundDigiRatio           = nullptr;  /// digis found per event
-  TH2I* fhCorrectVsFound           = nullptr;  /// correct digis per event vs found digis per event
-  TH2I* fhCorrectVsFoundNoNoise    = nullptr;  /// correct digis per event vs found digis per event, disregarding noise
-  TH1F* fhTimeOffset               = nullptr;  /// difference between true event time and seed time
+  TH1F* fhLinkedTriggersPerMCEvent  = nullptr;  /// linked triggers per MC event
+  TH1F* fhMatchedTriggersPerMCEvent = nullptr;  /// matchted triggers per MC event
+  TH1F* fhLinkedMCEventsPerTrigger  = nullptr;  /// linked MC events per trigger
+  TH1F* fhCorrectDigiRatio          = nullptr;  /// correct digis per event
+  TH1F* fhCorrectDigiRatioNoNoise   = nullptr;  /// correct digis per event, disregarding noise
+  TH1F* fhNoiseDigiRatio            = nullptr;  /// noise digis per event
+  TH1F* fhFoundDigiRatio            = nullptr;  /// digis found per event
+  TH2I* fhCorrectVsFound            = nullptr;  /// correct digis per event vs found digis per event
+  TH2I* fhCorrectVsFoundNoNoise     = nullptr;  /// correct digis per event vs found digis per event, disregarding noise
+  TH1F* fhTimeOffset                = nullptr;  /// difference between true event time and seed time
 
   CbmQaCanvas* fCanv;  ///summary canvas
 
diff --git a/reco/eventbuilder/digis/CbmSeedFinderSlidingWindow.cxx b/reco/eventbuilder/digis/CbmSeedFinderSlidingWindow.cxx
index a9ed0943548a38dd076d8f5dd14b9a12691872b6..e2249d02ab55ffe60062590b8e71c45a17a33fdd 100644
--- a/reco/eventbuilder/digis/CbmSeedFinderSlidingWindow.cxx
+++ b/reco/eventbuilder/digis/CbmSeedFinderSlidingWindow.cxx
@@ -32,6 +32,8 @@ void CbmSeedFinderSlidingWindow::FillSeedTimes(const std::vector<inType>* vIn, c
   // Reset output array
   fvSeedTimes->clear();
 
+  if (fQa) { fQa->ResetPerTsStorage(); }
+
   // Ideal mode ignores digi input and copies MC event list
   if (fbIdealMode) {
 
@@ -53,7 +55,7 @@ void CbmSeedFinderSlidingWindow::FillSeedTimes(const std::vector<inType>* vIn, c
     }
     if (fQa) {
       for (uint32_t i = 0; i < fvSeedTimes->size(); i++) {
-        fQa->FillQaInfo(i, i, &eventMatches, fvSeedTimes->at(i));
+        fQa->FillQaSeedInfo(i, i, &eventMatches, fvSeedTimes->at(i));
       }
     }
     return;
@@ -98,7 +100,7 @@ void CbmSeedFinderSlidingWindow::FillSeedTimes(const std::vector<inType>* vIn, c
       fvSeedTimes->push_back(seedTime);
 
       if (vDigiMatch && fQa) {  // QA mode
-        fQa->FillQaInfo(winStartN, i, vDigiMatch, seedTime);
+        fQa->FillQaSeedInfo(winStartN, i, vDigiMatch, seedTime);
       }
       nDigisWin = 0;
 
@@ -116,6 +118,8 @@ void CbmSeedFinderSlidingWindow::FillSeedTimes(const std::vector<inType>* vIn, c
     }
   }
 
+  if (fQa) { fQa->FillQaMCInfo(); }
+
   //if (fQa && vDigiMatch) {  // QA mode
   //  std::cout << "CbmSeedFinderSlidingWindow::FillSeedTimes(): Found " << GetNofSeeds() << " seeds for this timeslice."
   //            << std::endl;