CbmRoot
Loading...
Searching...
No Matches
DigiEventQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese [committer], P.-A. Loizeau */
4
5#include "DigiEventQa.h"
6
7#include "Histogram.h"
8
9#include <functional>
10#include <iomanip>
11#include <iostream>
12#include <sstream>
13#include <string>
14
16using std::string;
17using std::vector;
18
19namespace cbm::algo::evbuild
20{
21
22 // --- Execution --------------------------------------------------------
23 DigiEventQaData DigiEventQa::operator()(const vector<DigiEvent>& events) const
24 {
25 // --- Instantiate return object
26 DigiEventQaData result;
27 for (const auto& entry : fConfig.fData) {
28 ECbmModuleId subsystem = entry.first;
29 auto& detConfig = fConfig.fData.at(subsystem);
30 result.fDigiTimeHistos[subsystem] =
31 &(result.fHistContainer.fvH1.emplace_front(DigiEventQaConfig::GetDigiTimeHistoName(subsystem), "",
32 detConfig.fNumBins, detConfig.fMinValue, detConfig.fMaxValue));
33 }
34
35 // --- Event loop. Fill histograms.
36 for (const auto& event : events) { //
37 for (auto& subsystem : result.fDigiTimeHistos) { //
38 QaDigiTimeInEvent(event, subsystem.first, subsystem.second);
39 }
40 }
41 result.fNumEvents = events.size();
42
43 return result;
44 }
45 // --------------------------------------------------------------------------
46
47
48 // --- QA: digi time within event ----------------------------------------
49 void DigiEventQa::QaDigiTimeInEvent(const DigiEvent& event, ECbmModuleId system, H1D* histo) const
50 {
51 switch (system) {
52
53 case ECbmModuleId::kBmon: FillDeltaT<CbmBmonDigi>(event.fBmon, event.fTime, histo); break;
54 case ECbmModuleId::kSts: FillDeltaT<CbmStsDigi>(event.fSts, event.fTime, histo); break;
55 case ECbmModuleId::kMuch: FillDeltaT<CbmMuchDigi>(event.fMuch, event.fTime, histo); break;
56 case ECbmModuleId::kRich: FillDeltaT<CbmRichDigi>(event.fRich, event.fTime, histo); break;
57 case ECbmModuleId::kTrd: FillDeltaT<CbmTrdDigi>(event.fTrd, event.fTime, histo); break;
58 case ECbmModuleId::kTrd2d: FillDeltaT<CbmTrdDigi>(event.fTrd2d, event.fTime, histo); break;
59 case ECbmModuleId::kTof: FillDeltaT<CbmTofDigi>(event.fTof, event.fTime, histo); break;
60 case ECbmModuleId::kFsd: FillDeltaT<CbmFsdDigi>(event.fFsd, event.fTime, histo); break;
61 case ECbmModuleId::kPsd: FillDeltaT<CbmPsdDigi>(event.fPsd, event.fTime, histo); break;
62 default: throw std::runtime_error("DigiEventQa: Invalid system Id " + ::ToString(system));
63 }
64 }
65 // --------------------------------------------------------------------------
66
67
68 // ----- Info to string -------------------------------------------------
69 std::string DigiEventQa::ToString() const
70 {
71 std::stringstream out;
72 out << "--- Using DigiEventQa with parameters:";
73 for (const auto& entry : fConfig.fData) {
74 out << "\n " << std::left << std::setw(5) << ::ToString(entry.first) << ": ";
75 out << "bins" << std::right << std::setw(5) << entry.second.fNumBins;
76 out << " range [" << std::right << std::setw(5) << entry.second.fMinValue;
77 out << ", " << std::right << std::setw(5) << entry.second.fMaxValue << "] ns";
78 }
79 return out.str();
80 }
81 // --------------------------------------------------------------------------
82
83
84} // namespace cbm::algo::evbuild
ECbmModuleId
Definition CbmDefs.h:39
@ kTrd
Transition Radiation Detector.
@ kTof
Time-of-flight Detector.
@ kPsd
Projectile spectator detector.
@ kSts
Silicon Tracking System.
@ kTrd2d
TRD-FASP Detector (FIXME)
@ kMuch
Muon detection system.
@ kFsd
Forward spectator detector.
@ kRich
Ring-Imaging Cherenkov Detector.
ROOT-free implementation of a histogram.
DigiEventQaData operator()(const std::vector< DigiEvent > &events) const
Execution.
std::string ToString() const
Info to string.
void QaDigiTimeInEvent(const DigiEvent &event, ECbmModuleId system, qa::H1D *histo) const
Fill histogram with digi time within event.
void FillDeltaT(gsl::span< const Digi > digis, double eventTime, qa::H1D *histo) const
Fill histogram with digi time within event.
1D-histogram
PODVector< CbmRichDigi > fRich
Unpacked RICH digis.
Definition DigiData.h:38
PODVector< CbmTrdDigi > fTrd
Unpacked TRD digis.
Definition DigiData.h:36
PODVector< CbmStsDigi > fSts
Unpacked STS digis.
Definition DigiData.h:32
PODVector< CbmTrdDigi > fTrd2d
Unpacked TRD2D digis.
Definition DigiData.h:37
PODVector< CbmFsdDigi > fFsd
Unpacked FSD digis.
Definition DigiData.h:40
PODVector< CbmTofDigi > fTof
Unpacked TOF digis.
Definition DigiData.h:34
PODVector< CbmPsdDigi > fPsd
Unpacked PSD digis.
Definition DigiData.h:39
PODVector< CbmMuchDigi > fMuch
Unpacked MUCH digis.
Definition DigiData.h:33
PODVector< CbmBmonDigi > fBmon
Unpacked Bmon digis.
Definition DigiData.h:35
Event data with event number and trigger time.
Definition DigiData.h:79
double fTime
Event trigger time [ns].
Definition DigiData.h:82
std::map< ECbmModuleId, DigiEventQaDetConfig > fData
Definition DigiEventQa.h:57
static std::string GetDigiTimeHistoName(const ECbmModuleId &subsystem)
Definition DigiEventQa.h:75
QA results for CbmDigiEvent objects.
Definition DigiEventQa.h:26
qa::HistogramContainer fHistContainer
Definition DigiEventQa.h:27
std::unordered_map< ECbmModuleId, qa::H1D * > fDigiTimeHistos
Definition DigiEventQa.h:28
std::forward_list< qa::H1D > fvH1
List of 1D-histograms.