CbmRoot
Loading...
Searching...
No Matches
CbmRichMCbmDenoiseQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2024 UGiessen, Giessen
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Martin Beyer [committer] */
4
6
7#include "CbmEvent.h"
8#include "CbmHistManager.h"
9#include "CbmRichHit.h"
10#include "CbmRichRing.h"
11
12#include <FairRootManager.h>
13#include <Logger.h>
14
15#include <TCanvas.h>
16#include <TClonesArray.h>
17#include <TDirectory.h>
18#include <TEllipse.h>
19#include <TFile.h>
20#include <TH2D.h>
21
23{
24 FairRootManager* manager = FairRootManager::Instance();
25 if (!manager) LOG(fatal) << GetName() << "::Init: No FairRootManager";
26
27 fCbmEvents = static_cast<TClonesArray*>(manager->GetObject("CbmEvent"));
28 LOG(info) << GetName() << "::Init: CbmEvent array " << (fCbmEvents ? "found" : "not found");
29
30 fRichHits = static_cast<TClonesArray*>(manager->GetObject("RichHit"));
31 if (!fRichHits) LOG(fatal) << GetName() << "::Init: No RichHit array!";
32
33 fRichRings = static_cast<TClonesArray*>(manager->GetObject("RichRing"));
34 if (!fRichRings) LOG(fatal) << GetName() << "::Init: No RichRing array!";
35
37
38 return kSUCCESS;
39}
40
42{
43 fHM = std::make_unique<CbmHistManager>();
44 fHMSed = std::make_unique<CbmHistManager>();
45
46 // 1D Histograms
47 fHM->Create1<TH1D>("fhRingsPerEvent", "number of ring per event;# rings per event;Entries", 6, -0.5, 5.5);
48
49 fHM->Create1<TH1D>("fhRingRadius", "ring radius;ring radius [cm];Entries", 100, 0., 7.);
50
51 fHM->Create1<TH1D>("fhRingRadiusUp", "ring radius upper half;ring radius [cm]; count", 100, 0., 7.);
52
53 fHM->Create1<TH1D>("fhRingRadiusDown", "ring radius lower half;ring radius [cm]; count", 100, 0., 7.);
54
55 fHM->Create1<TH1D>("fhRingHits", "number of ring hits;# ring hits;Entries", 50, -0.5, 49.5);
56
57 fHM->Create1<TH1D>("fhRingHitsUp", "number of ring hits upper half;# ring hits;Entries", 50, -0.5, 49.5);
58
59 fHM->Create1<TH1D>("fhRingHitsDown", "number of ring hits lower half;# ring hits;Entries", 50, -0.5, 49.5);
60
61 fHM->Create1<TH1D>("fhRingHitsTimeDiff", "ringhit to ring time difference;hittime-ringtime [ns]; count", 100, -10.,
62 10.);
63
64 // 2D Histograms
65 fHM->Create2<TH2D>("fhRichHitsXY", "position of rich hits;X [cm];Y [cm];Entries", 37, -11.2, 11.2, 82, -24.115,
66 23.885);
67
68 fHM->Create2<TH2D>("fhRingHitsXY", "position of ring hits;X [cm];Y [cm];Entries", 37, -11.2, 11.2, 82, -24.115,
69 23.885);
70
71 fHM->Create2<TH2D>("fhRingCenterXY", "position of ring centers;X [cm];Y [cm];Entries", 37, -11.2, 11.2, 82, -24.115,
72 23.885);
73}
74
75void CbmRichMCbmDenoiseQa::Exec(Option_t* /*option*/)
76{
77 LOG(debug) << GetName() << " TS " << fTsNum;
78 fTsNum++;
79
80 if (fCbmEvents) {
81 Int_t nEvents = fCbmEvents->GetEntriesFast();
82 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
83 CbmEvent* event = static_cast<CbmEvent*>(fCbmEvents->At(iEvent));
84 if (!event) continue;
85 if (event->GetNofData(ECbmDataType::kRichRing) > 0) {
86 fHM->H1("fhRingsPerEvent")->Fill(static_cast<double>(event->GetNofData(ECbmDataType::kRichRing)));
87 if (fnSEDs < fMaxSEDs) {
88 DrawSED(event);
89 fnSEDs++;
90 }
91 }
92 Process(event);
93 fEventNum++;
94 }
95 }
96 else {
97 Process(nullptr);
98 }
99}
100
102{
103 // RICH Hits
104 Int_t nRichHits = event ? static_cast<Int_t>(event->GetNofData(ECbmDataType::kRichHit)) : fRichHits->GetEntriesFast();
105 for (int i = 0; i < nRichHits; i++) {
106 Int_t hitIndex = event ? static_cast<Int_t>(event->GetIndex(ECbmDataType::kRichHit, static_cast<uint32_t>(i))) : i;
107 CbmRichHit* hit = static_cast<CbmRichHit*>(fRichHits->At(hitIndex));
108 if (!hit) continue;
109 fHM->H2("fhRichHitsXY")->Fill(hit->GetX(), hit->GetY());
110 }
111
112 // RICH Rings
113 Int_t nRichRings =
114 event ? static_cast<Int_t>(event->GetNofData(ECbmDataType::kRichRing)) : fRichRings->GetEntriesFast();
115 for (int i = 0; i < nRichRings; i++) {
116 Int_t ringIndex =
117 event ? static_cast<Int_t>(event->GetIndex(ECbmDataType::kRichRing, static_cast<uint32_t>(i))) : i;
118 CbmRichRing* ring = static_cast<CbmRichRing*>(fRichRings->At(ringIndex));
119 if (!ring) continue;
120 fHM->H1("fhRingRadius")->Fill(ring->GetRadius());
121 fHM->H1(std::string("fhRingRadius") + std::string((ring->GetCenterY() > 0) ? "Up" : "Down"))
122 ->Fill(ring->GetRadius());
123 fHM->H1("fhRingHits")->Fill(ring->GetNofHits());
124 fHM->H1(std::string("fhRingHits") + std::string((ring->GetCenterY() > 0) ? "Up" : "Down"))
125 ->Fill(ring->GetNofHits());
126 fHM->H2("fhRingCenterXY")->Fill(ring->GetCenterX(), ring->GetCenterY());
127 // Ring hits
128 Int_t nRingHits = ring->GetNofHits();
129 for (int j = 0; j < nRingHits; j++) {
130 CbmRichHit* hit = static_cast<CbmRichHit*>(fRichHits->At(static_cast<Int_t>(ring->GetHit(j))));
131 if (!hit) continue;
132 fHM->H1("fhRingHitsTimeDiff")->Fill(hit->GetTime() - ring->GetTime());
133 fHM->H2("fhRingHitsXY")->Fill(hit->GetX(), hit->GetY());
134 }
135 }
136}
137
139{
140 std::string hName = "fhSED" + std::to_string(fEventNum);
141 auto c = fHMSed->CreateCanvas(hName, hName, 800, 1505);
142 c->SetMargin(0.15, 0.15, 0.1, 0.1);
143 c->cd();
144 TH2D* h = new TH2D(hName.c_str(), hName.c_str(), 37, -11.2, 11.2, 82, -24.115, 23.885);
145 h->SetXTitle("X [cm]");
146 h->SetYTitle("Y [cm]");
147 h->SetZTitle("NN classification 0: noise, 1: signal");
148 h->GetZaxis()->SetTitleOffset(1.3);
149 fHMSed->Add(hName, h);
150
151 Int_t nHits = static_cast<Int_t>(event->GetNofData(ECbmDataType::kRichHit));
152 for (int i = 0; i < nHits; i++) {
153 Int_t hitIndex = static_cast<Int_t>(event->GetIndex(ECbmDataType::kRichHit, static_cast<uint32_t>(i)));
154 CbmRichHit* hit = static_cast<CbmRichHit*>(fRichHits->At(hitIndex));
155 if (!hit) continue;
156 if (h->GetBinContent(h->FindBin(hit->GetX(), hit->GetY())) > 0) {
157 LOG(info) << GetName() << "::DrawSED channel already filled for this event. Skipping.";
158 continue;
159 }
160 if (!hit->GetIsNoiseNN()) {
161 h->Fill(hit->GetX(), hit->GetY(), 1.0);
162 }
163 else {
164 h->Fill(hit->GetX(), hit->GetY(), 0.001);
165 }
166 }
167 h->Draw("colz");
168
169 Int_t nRings = static_cast<Int_t>(event->GetNofData(ECbmDataType::kRichRing));
170 for (int i = 0; i < nRings; i++) {
171 Int_t ringIndex = static_cast<Int_t>(event->GetIndex(ECbmDataType::kRichRing, static_cast<uint32_t>(i)));
172 CbmRichRing* ring = static_cast<CbmRichRing*>(fRichRings->At(ringIndex));
173 if (!ring) continue;
174 TEllipse* circle = new TEllipse(ring->GetCenterX(), ring->GetCenterY(), ring->GetAaxis(), ring->GetBaxis(), 0.,
175 360., ring->GetPhi() * 180. / TMath::Pi());
176 circle->SetFillStyle(0);
177 circle->SetLineWidth(2);
178 circle->Draw("same");
179 }
180}
181
183{
184 TDirectory* oldir = gDirectory;
185 TFile* outFile = FairRootManager::Instance()->GetOutFile();
186 if (outFile) {
187 outFile->mkdir(GetName());
188 outFile->cd(GetName());
189 fHM->WriteToFile();
190 outFile->cd();
191 outFile->mkdir((std::string(GetName()) + std::string("SEDs")).c_str());
192 outFile->cd((std::string(GetName()) + std::string("SEDs")).c_str());
193 fHMSed->WriteCanvasToFile();
194 fHM->Clear();
195 fHMSed->Clear();
196 }
197 gDirectory->cd(oldir->GetPath());
198}
Histogram manager.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
double GetTime() const
Definition CbmHit.h:76
double GetY() const
Definition CbmPixelHit.h:74
double GetX() const
Definition CbmPixelHit.h:73
bool GetIsNoiseNN() const
Definition CbmRichHit.h:68
void Exec(Option_t *option)
std::unique_ptr< CbmHistManager > fHMSed
void Process(CbmEvent *event)
Process data and fill histograms.
std::unique_ptr< CbmHistManager > fHM
void DrawSED(CbmEvent *event)
uint32_t GetHit(int32_t i) const
Definition CbmRichRing.h:39
float GetRadius() const
Definition CbmRichRing.h:79
double GetTime() const
Definition CbmRichRing.h:99
double GetAaxis() const
Definition CbmRichRing.h:80
int32_t GetNofHits() const
Definition CbmRichRing.h:37
float GetCenterX() const
Definition CbmRichRing.h:77
double GetPhi() const
Definition CbmRichRing.h:84
float GetCenterY() const
Definition CbmRichRing.h:78
double GetBaxis() const
Definition CbmRichRing.h:81
Data class with information on a STS local track.