CbmRoot
Loading...
Searching...
No Matches
CbmRichMCbmHitProducer.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2024 UGiessen/JINR-LIT, Giessen/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Adrian Amatus Weber, Semen Lebedev [committer], Martin Beyer */
4
6
7#include "CbmDigiManager.h"
8#include "CbmEvent.h"
9#include "CbmRichDetectorData.h" // for CbmRichPmtData, CbmRichPixelData
10#include "CbmRichDigi.h"
12#include "CbmRichGeoManager.h"
13#include "CbmRichHit.h"
14#include "CbmRichPoint.h"
15#if HAVE_ONNXRUNTIME
17#endif
18
19#include <FairRootManager.h>
20#include <Logger.h>
21
22#include <TClonesArray.h>
23#include <TStopwatch.h>
24
25#include <iomanip>
26#include <iostream>
27
28using std::fixed;
29using std::right;
30using std::setprecision;
31using std::setw;
32
33using namespace std;
34
35
37 : FairTask("CbmRichMCbmHitProducer")
38 , fRichHits(NULL)
39 , fNofTs(0)
40 , fHitError(0.6 / sqrt(12.))
41 ,
42 //fMappingFile("mRICH_Mapping_vert_20190318.geo")
43 fMappingFile("mRICH_Mapping_vert_20190318_elView.geo")
44 , fICD_offset_read()
45{
46}
47
49{
50 FairRootManager* manager = FairRootManager::Instance();
51 manager->Write();
52}
53
55
57{
58 FairRootManager* manager = FairRootManager::Instance();
59
60 fCbmEvents = dynamic_cast<TClonesArray*>(manager->GetObject("CbmEvent"));
61 if (fCbmEvents == nullptr) {
62 LOG(info) << GetName() << "::Init() CbmEvent NOT found \n";
63 }
64 else {
65 LOG(info) << GetName() << "::Init() CbmEvent found";
66 }
67
69 fDigiMan->Init();
70 if (!fDigiMan->IsPresent(ECbmModuleId::kRich)) Fatal("CbmRichMCbmHitProducer::Init", "No RichDigi array!");
71
72 fRichHits = new TClonesArray("CbmRichHit");
73 manager->Register("RichHit", "RICH", fRichHits, IsOutputBranchPersistent("RichHit"));
74
76
77 for (auto& a : fICD_offset_read)
78 a = 0.;
80
81#if HAVE_ONNXRUNTIME
82 if (fUseDenoiseNN) {
83 fDenoiseCnn = std::make_unique<CbmRichMCbmDenoiseCnn>();
84 fDenoiseCnn->Init();
85 fDenoiseCnn->SetClassificationThreshold(fDenoiseCnnThreshold);
86 LOG(info) << "CbmRichMCbmHitProducer: Initialized CbmRichMCbmDenoiseCnn";
87 }
88#endif
89
90 return kSUCCESS;
91}
92
94{
95
96 string line;
97 ifstream file(fMappingFile);
98 if (!file.is_open()) {
99 std::cout << "<CbmRichMCbmHitProducer::InitMapping>: Unable to open mapping file:" << fMappingFile.c_str()
100 << std::endl;
101 }
102
103 fRichMapping.clear();
104
105 while (getline(file, line)) {
106
107 istringstream iss(line);
108 vector<std::string> results(istream_iterator<string>{iss}, std::istream_iterator<string>());
109 if (results.size() != 8) continue;
110
112 data.fTrbId = stoi(results[0], nullptr, 16);
113 data.fChannel = stoi(results[1]);
114 data.fX = stod(results[6]);
115 data.fY = stod(results[7]);
116 data.fZ = 348.;
117
118 data.fX -= 6.3; //Shift by 1Pmt + PmtGap + 1cm
119
120 Int_t adr = ((data.fTrbId << 16) | (data.fChannel & 0x00FF));
121
122 // cout << data.fTrbId << " " << data.fChannel << " " << data.fX << " " << data.fY << " " << adr << endl;
123
124 fRichMapping[adr] = data;
125 }
126 file.close();
127
128 //cout << "Mapping size:" << fRichMapping.size() <<endl;
129}
130
131void CbmRichMCbmHitProducer::Exec(Option_t* /*option*/)
132{
133 fNofDigis = 0;
134 TStopwatch timer;
135 timer.Start();
136 fRichHits->Delete();
137
138 // if CbmEvent does not exist then process standard event.
139 // if CbmEvent exists then proceed all events in time slice.
140 Int_t nUnits = fCbmEvents ? fCbmEvents->GetEntriesFast() : 1;
141 if (fCbmEvents) fNofEvents += nUnits;
142 for (Int_t iUnit = 0; iUnit < nUnits; iUnit++) {
143 CbmEvent* event = fCbmEvents ? static_cast<CbmEvent*>(fCbmEvents->At(iUnit)) : nullptr;
144 ProcessData(event);
145 }
146 timer.Stop();
147 fTotalTime += timer.RealTime();
148
149 std::stringstream logOut;
150 logOut << setw(20) << left << GetName() << "[";
151 logOut << fixed << setw(8) << setprecision(1) << right << timer.RealTime() * 1000. << " ms] ";
152 logOut << "TS " << fNofTs;
153 if (fCbmEvents) logOut << ", events " << nUnits;
154 logOut << ", digis " << fNofDigis;
155 logOut << ", hits " << fRichHits->GetEntriesFast();
156 LOG(info) << logOut.str();
158 fTotalNofHits += fRichHits->GetEntriesFast();
159 fNofTs++;
160}
161
163{
164 TStopwatch timer;
165 timer.Start();
166 if (event != NULL) {
167 LOG(debug) << "CbmRichMCbmHitProducer CbmEvent mode. CbmEvent # " << event->GetNumber();
168 Int_t nofDigis = event->GetNofData(ECbmDataType::kRichDigi);
169 //LOG(info) << "nofDigis: " << nofDigis;
170
171 fNofHits = 0;
172 for (Int_t iDigi = 0; iDigi < nofDigis; iDigi++) {
173 Int_t digiIndex = event->GetIndex(ECbmDataType::kRichDigi, iDigi);
174 ProcessDigi(event, digiIndex);
175 }
176 LOG(debug) << "nofDigis: " << nofDigis << "\t\t "
177 << "nofHits : " << fNofHits;
178 fNofHits = 0;
179 }
180 else {
181 for (Int_t iDigi = 0; iDigi < fDigiMan->GetNofDigis(ECbmModuleId::kRich); iDigi++) {
182 ProcessDigi(event, iDigi);
183 }
184 }
185 timer.Stop();
186 fHitProducerTime += timer.RealTime();
187
188// Process all hits with DenoiseNN, setting value RichHit.fIsNoiseNN
189#if HAVE_ONNXRUNTIME
190 if (fUseDenoiseNN) {
191 timer.Start();
192 fDenoiseCnn->Process(event, fRichHits);
193 timer.Stop();
194 fDenoiseCnnTime += timer.RealTime();
195 }
196#endif
197}
198
200{
201 const CbmRichDigi* digi = fDigiMan->Get<CbmRichDigi>(digiIndex);
202 if (digi == nullptr) return;
203 if (digi->GetAddress() < 0) return;
204 Int_t DiRICH_Add = (digi->GetAddress() >> 16) & 0xFFFF;
205 if (DiRICH_Add == 0x7901 || DiRICH_Add == 0x7902) return; // TRBaddresses 0x7901 and 0x7902 are for FSD/NCAL
206 fNofDigis++;
207 if (isInToT(digi->GetToT())) {
208 TVector3 posPoint;
210
211 if (data != nullptr) {
212 posPoint.SetXYZ(data->fX, data->fY, data->fZ);
213 }
214 else {
215 LOG(info) << "CbmRichMCbmHitProducer: No Node for 0x" << std::hex << digi->GetAddress() << std::dec
216 << " found. Using ASCII File! ";
218
219 /* if (dataAscii == NULL) {
220 LOG(error) << "CbmRichMCbmHitProducer: No Position for "<<digi->GetAddress()<<"found! ";
221 return;
222 }*/
223 posPoint.SetXYZ(dataAscii.fX, dataAscii.fY, dataAscii.fZ);
224 }
225
226 //std::cout<<std::hex<<digi->GetAddress()<<std::dec<<" "<<data->fX<<" "<<data->fY<<" "<<data->fZ<<std::endl;
227 if (!RestrictToAcc(posPoint)) return;
228 if (!RestrictToFullAcc(posPoint)) return;
229 if (!RestrictToAerogelAccDec2019(posPoint)) return;
230 AddHit(event, posPoint, digi, digiIndex, data->fPmtId);
231 }
232}
233
234
235void CbmRichMCbmHitProducer::AddHit(CbmEvent* event, TVector3& posHit, const CbmRichDigi* digi, Int_t index,
236 Int_t PmtId)
237{
238 Int_t nofHits = fRichHits->GetEntriesFast();
239
240 int32_t DiRICH_Addr = digi->GetAddress();
241 unsigned int addr = (((DiRICH_Addr >> 24) & 0xF) * 18 * 32) + (((DiRICH_Addr >> 20) & 0xF) * 2 * 32)
242 + (((DiRICH_Addr >> 16) & 0xF) * 32) + ((DiRICH_Addr & 0xFFFF) - 0x1);
243
244 new ((*fRichHits)[nofHits]) CbmRichHit();
245 CbmRichHit* hit = (CbmRichHit*) fRichHits->At(nofHits);
246 hit->SetPosition(posHit);
247 hit->SetDx(fHitError);
248 hit->SetDy(fHitError);
249 hit->SetRefId(index);
250 if (fDoICD) {
251 hit->SetTime(digi->GetTime() - fICD_offset_read.at(addr));
252 }
253 else {
254 hit->SetTime(digi->GetTime());
255 }
256
257 hit->SetToT(digi->GetToT());
258 hit->SetAddress(digi->GetAddress());
259 hit->SetPmtId(PmtId);
260 // Not Implemented in Digi: hit->SetPmtId(digi->GetPmtId());
261
262 if (event != NULL) {
263 event->AddData(ECbmDataType::kRichHit, nofHits);
264 fNofHits++;
265 }
266}
267
268
270{
271 fRichHits->Clear();
272 std::cout << std::endl;
273 LOG(info) << "=====================================";
274 LOG(info) << GetName() << ": Run summary";
275 LOG(info) << "Time slices : " << fNofTs;
276 LOG(info) << "Digis / TS : " << fixed << setprecision(2) << fTotalNofDigis / Double_t(fNofTs);
277 LOG(info) << "Hits / TS : " << fixed << setprecision(2) << fTotalNofHits / Double_t(fNofTs);
278 LOG(info) << "Time / TS : " << fixed << setprecision(2) << 1000. * fTotalTime / Double_t(fNofTs) << " ms";
279 if (fCbmEvents) {
280 LOG(info) << "Events : " << fNofEvents;
281 LOG(info) << "Events / TS : " << fixed << setprecision(2) << fNofEvents / Double_t(fNofTs);
282 if (fNofEvents > 0) {
283 LOG(info) << "Digis / ev : " << fixed << setprecision(2) << fTotalNofDigis / Double_t(fNofEvents);
284 LOG(info) << "Hits / ev : " << fixed << setprecision(2) << fTotalNofHits / Double_t(fNofEvents);
285 LOG(info) << "Time / ev : " << fixed << setprecision(2) << 1000. * fTotalTime / Double_t(fNofEvents)
286 << " ms";
287 }
288 }
289
290 TString eventOrTsStr = fCbmEvents && fNofEvents > 0 ? "ev: " : "TS: ";
291 Double_t eventOrTsValue = Double_t(fCbmEvents && fNofEvents > 0 ? fNofEvents : fNofTs);
292 if (fTotalTime != 0.) {
293 LOG(info) << "===== Time by task (real time) ======";
294 LOG(info) << "HitProducer / " << eventOrTsStr << fixed << setprecision(2) << setw(9) << right
295 << 1000. * fHitProducerTime / eventOrTsValue << " ms [" << setw(5) << right
296 << 100 * fHitProducerTime / fTotalTime << " %]";
297#if HAVE_ONNXRUNTIME
298 if (fUseDenoiseNN)
299 LOG(info) << "Denoise NN / " << eventOrTsStr << fixed << setprecision(2) << setw(9) << right
300 << 1000. * fDenoiseCnnTime / eventOrTsValue << " ms [" << setw(5) << right
301 << 100 * fDenoiseCnnTime / fTotalTime << " %]";
302#endif
303 }
304 LOG(info) << "=====================================\n";
305}
306
307
309{
310
311 if (!fDoToT) return true;
312
313 if ((ToT > fToTLimitLow) && (ToT < fToTLimitHigh)) {
314 return true;
315 }
316 else {
317 return false;
318 }
319}
320
321
323{
324 Double_t x = pos.X();
325 Double_t y = pos.Y();
326
327 return this->RestrictToAcc(x, y);
328}
329
331{ //check if Track is in mRICH acceptance
332 if (fRestrictToAcc == false) return true;
333 bool inside = false;
334 if (x >= -6.25 && x <= -1.05) {
335 // left part of mRICH
336 if (y >= 8 && y <= 15.9) {
337 inside = true;
338 }
339 }
340 else if (x > -1.05 && x <= 4.25) {
341 //right part
342 if (y >= 8 && y <= 13.2) {
343 inside = true;
344 }
345 }
346
347 return inside;
348}
349
351{
352 Double_t x = pos.X();
353 Double_t y = pos.Y();
354
355 return this->RestrictToFullAcc(x, y);
356}
357
359{ //check if Track is in mRICH acceptance
360 if (fRestrictToFullAcc == false) return true;
361 bool inside = false;
362 if (x >= -16.85 && x <= 4.25) { //TODO:check values
363 // left part of mRICH
364 if (y >= -23.8 && y <= 23.8) {
365 inside = true;
366 }
367 }
368
369 return inside;
370}
371
373{
374 Double_t x = pos.X();
375 Double_t y = pos.Y();
376
377 return this->RestrictToAerogelAccDec2019(x, y);
378}
379
381{ //check if Track is in mRICH acceptance
382 if (fRestrictToAerogelAccDec2019 == false) return true;
383 //bool inside = false;
384
385 if (y > 13.5) return false;
386 if (y > 8.0 && x < 12.5) return false;
387
388 return true;
389}
390
391void CbmRichMCbmHitProducer::read_ICD(std::array<Double_t, 2304>& ICD_offsets, unsigned int iteration)
392{
393 std::string line;
394 std::string filename = ("" == fIcdFilenameBase ? "icd_offset_it" : fIcdFilenameBase);
395 filename += Form("_%u.data", iteration);
396 std::ifstream icd_file(filename);
397 unsigned int lineCnt = 0;
398 if (icd_file.is_open()) {
399 while (getline(icd_file, line)) {
400 if (line[0] == '#') continue; // just a comment
401 std::istringstream iss(line);
402 unsigned int addr = 0;
403 Double_t value;
404 if (!(iss >> addr >> value)) {
405 LOG(info) << "A Problem accured in line " << lineCnt;
406 break;
407 } // error
408 lineCnt++;
409 ICD_offsets.at(addr) += value;
410 }
411 icd_file.close();
412 LOG(info) << "Loaded inter channel delay file " << filename << " for RICH.";
413 }
414 else {
415 LOG(info) << "Unable to open inter channel delay file " << filename;
416 }
417}
418
ClassImp(CbmConverterManager)
@ kRich
Ring-Imaging Cherenkov Detector.
static Int_t GetNofDigis(ECbmModuleId systemId)
static Bool_t IsPresent(ECbmModuleId systemId)
Presence of a digi branch.
InitStatus Init()
Initialisation.
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
void SetAddress(int32_t address)
Definition CbmHit.h:83
void SetRefId(int32_t refId)
Definition CbmHit.h:82
void SetTime(double time)
Definition CbmHit.h:85
void SetDx(double dx)
Definition CbmPixelHit.h:94
void SetDy(double dy)
Definition CbmPixelHit.h:95
void SetPosition(const TVector3 &pos)
Sets position of the hit.
static CbmRichDigiMapManager & GetInstance()
Return Instance of CbmRichGeoManager.
CbmRichPixelData * GetPixelDataByAddress(Int_t address)
Return CbmRichDataPixel by digi address.
int32_t GetAddress() const
Definition CbmRichDigi.h:43
double GetToT() const
Definition CbmRichDigi.h:78
double GetTime() const
Definition CbmRichDigi.h:72
virtual void SetPmtId(int32_t det)
Definition CbmRichHit.h:57
void SetToT(double tot)
Definition CbmRichHit.h:60
virtual ~CbmRichMCbmHitProducer()
Destructor.
std::map< Int_t, CbmRichMCbmMappingData > fRichMapping
bool isInToT(const double ToT)
void read_ICD(std::array< Double_t, 2304 > &offsets, unsigned int iteration)
void ProcessDigi(CbmEvent *event, Int_t digiIndex)
bool RestrictToAerogelAccDec2019(TVector3 &pos)
bool RestrictToFullAcc(TVector3 &pos)
void ProcessData(CbmEvent *event)
void AddHit(CbmEvent *event, TVector3 &posHit, const CbmRichDigi *digi, Int_t index, Int_t PmtId)
Add hit to the output array (and) CbmEvent if it is not NULL.
CbmRichMCbmHitProducer()
Default constructor.
virtual void SetParContainers()
Inherited from FairTask.
virtual void Finish()
Inherited from FairTask.
virtual InitStatus Init()
Inherited from FairTask.
std::array< Double_t, 2304 > fICD_offset_read
virtual void Exec(Option_t *option)
Inherited from FairTask.
Hash for CbmL1LinkKey.