CbmRoot
Loading...
Searching...
No Matches
CbmMuchFindHitsGem.cxx
Go to the documentation of this file.
1/* Copyright (C) 2008-2021 St. Petersburg Polytechnic University, St. Petersburg
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Vikas Singhal, Ajit Kumar, Ekata Nandy, Evgeny Kryshen [committer] , Ajay Kumar, Florian Uhlig */
4
5/*
6 * CbmMuchFindHitsGem.cxx
7 *
8 * Modified on 08/08/2019 : Hit reconstruction in Event (in time slice) and Time slice mode
9 * Default is time slice (kCbmTimeSlice) and it will run in event mode (kCbmEvent) if find event branch in the tree
10 * @authors Vikas Singhal and Ajit Kumar
11 *
12 //@author : Ekata Nandy (ekata@vecc.gov.in) since 21-06-19- Drift time correction for GEM and RPC have been done separately.
13 */
14#include "CbmMuchFindHitsGem.h"
15
16#include "CbmDigiManager.h"
17#include "CbmMuchCluster.h"
18#include "CbmMuchGeoScheme.h"
19#include "CbmMuchModule.h"
20#include "CbmMuchModuleGem.h"
21#include "CbmMuchPad.h"
22#include "CbmMuchPixelHit.h"
23#include "FairRootManager.h"
24#include "TClonesArray.h"
25#include "TFile.h"
26#include "TMath.h"
27#include "TStopwatch.h"
28//#include "CbmTimeSlice.h"
29#include "CbmMuchAddress.h"
30#include "CbmMuchBeamTimeDigi.h"
31#include "CbmMuchDigi.h"
32
33#include <algorithm>
34#include <iomanip>
35#include <iostream>
36
37using std::cout;
38using std::endl;
39using std::fixed;
40using std::left;
41using std::multimap;
42using std::right;
43using std::setw;
44using std::vector;
45
46// -------------------------------------------------------------------------
47CbmMuchFindHitsGem::CbmMuchFindHitsGem(const char* digiFileName, Int_t flag)
48 : FairTask("MuchFindHitsGem", 1)
49 , fDigiFile(digiFileName)
50 , fFlag(flag)
51 , fAlgorithm(3)
52 , fClusterSeparationTime(100.)
53 , fThresholdRatio(0.1)
54 , fEvent(0)
55 , fNofTimeslices(0)
56 ,
57 //fDigis(NULL),
58 fEvents(NULL)
59 , fClusterCharges()
60 , fLocalMax()
61 , fClusterPads()
62 , fNeighbours()
63 , fClusters(new TClonesArray("CbmMuchCluster", 1000))
64 , fHits(new TClonesArray("CbmMuchPixelHit", 1000))
65 , fGeoScheme(CbmMuchGeoScheme::Instance())
66 , fDigiIndices()
67 , fFiredPads()
68 ,
69 // fDaq(),
70 // fTimeSlice(NULL),
71 // fDigiData(),
72 fuClusters(0)
73{
74}
75
76// ----- Private method Init -------------------------------------------
78{
79 FairRootManager* ioman = FairRootManager::Instance();
80 //if (fDaq) fTimeSlice = (CbmTimeSlice*) ioman->GetObject("TimeSlice.");
81 //else fDigis = (TClonesArray*) ioman->GetObject("MuchDigi");
82
83 // --- Digi Manager for reading digis which were stored in vector
87
88 // fDigis will not be used now. Just for checking. Need to remove
89 /*fDigis = (TClonesArray*) ioman->GetObject("MuchDigi");
90 if (! fDigis)
91 fDigis = (TClonesArray*) ioman->GetObject("MuchBeamTimeDigi");
92 if (! fDigis)
93 fDigis = (TClonesArray*) ioman->GetObject("CbmMuchBeamTimeDigi");
94 if (! fDigis)
95 fDigis = (TClonesArray*) ioman->GetObject("CbmMuchDigi");
96 if (! fDigis)
97 LOG(info) << "MuchFindHitsGem: No MuchDigi or MuchBeamTimeDigi or CbmMuchDigi or CbmMuchBeamTimeDigi exist";
98 */
99
100 // Implementation of Event by event execution after TimeSlice Mode and Event Building should have one Event branch.
101
102 fEvents = dynamic_cast<TClonesArray*>(FairRootManager::Instance()->GetObject("Event"));
103 if (!fEvents) fEvents = dynamic_cast<TClonesArray*>(FairRootManager::Instance()->GetObject("CbmEvent"));
104
105 if (!fEvents) {
106 LOG(info) << GetName() << ": No event branch present.";
107 // return kfatal;
108 }
109 else {
110 fEventMode = kTRUE;
111 //fMode=kCbmEvent;
112 LOG(info) << GetName() << "TimeSlice: Event-by-event mode after Event Building selected. ";
113 }
114
115 ioman->Register("MuchCluster", "Cluster in MUCH", fClusters, IsOutputBranchPersistent("MuchCluster"));
116 ioman->Register("MuchPixelHit", "Hit in MUCH", fHits, IsOutputBranchPersistent("MuchPixelHit"));
117
118 // Initialize GeoScheme
120 TFile* oldFile = gFile;
121 TDirectory* oldDir = gDirectory;
122
123 TFile* file = new TFile(fDigiFile);
124 LOG_IF(fatal, !file) << "Could not open file " << fDigiFile;
125 TObjArray* stations = file->Get<TObjArray>("stations");
126 LOG_IF(fatal, !stations) << "TObjArray stations not found in file " << fDigiFile;
127
128 //For X, Y position correction according to different geometry file and its tag.
129
130 if (fDigiFile.Contains("mcbm")) {
131 if (fDigiFile.Contains("v19a")) {
132 fGemTX = 18.5;
133 fGemTY = 80.5;
134 }
135 else if (fDigiFile.Contains("v20a")) {
136 fGemTX = 18.5;
137 fGemTY = 80.0;
138 }
139 else if (fDigiFile.Contains("v22j")) { // for high intensity runs during June 2022 mMuCh in acceptance
140 fGemTX = 8.0;
141 fGemTY = 81.5;
142 fRpcTX = 66.25; //RPC introduced during 2022 data taking
143 fRpcTY = -70.5;
144 }
145 else if (fDigiFile.Contains("v22k")) { // during benchmark runs. mMuCh is out of acceptance
146 fGemTX = -44.5;
147 fGemTY = 81.5;
148 fRpcTX = 48.0;
149 fRpcTY = -70.0;
150 }
151 }
152
153 file->Close();
154 file->Delete();
156 gFile = oldFile;
157 gDirectory = oldDir;
158
159 fGeoScheme->Init(stations, fFlag);
160 return kSUCCESS;
161}
162// -------------------------------------------------------------------------
163
164// ----- Task execution ------------------------------------------------
166{
167 TStopwatch timer;
168 timer.Start();
169 // fDigiData.clear();
170 // Removing SetDaq functionality as Cluster and Hit Finder algorithm is same for both the Time Based and Event Based mode.
171 //if (fDaq) ;
172 //fDigiData = fTimeSlice->GetMuchData();
173 // else {
174 //LOG(debug)<<"Start Reading digi from a module ";
175
176 /*for (Int_t iDigi = 0; iDigi < fDigiManager->GetNofDigis(ECbmModuleId::kMuch); iDigi++) {
177 //Reading digi from CbmDigiManager which stors digis in vector
178 const auto * digi;
179 if(!bBeamTimeDigi) digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(iDigi);
180 else digi = (CbmMuchBeamTimeDigi*) fDigiManager->Get<CbmMuchBeamTimeDigi>(iDigi);
181 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigis->At(iDigi);
182 CbmMuchModule* module = fGeoScheme->GetModuleByDetId(digi->GetAddress()); //AZ
183 //std::cout << module->GetDetectorType() << std::endl; //AZ
184 if (module->GetDetectorType() == 2) continue; //AZ - skip 2-D straws
185 fDigiData.push_back(*digi);
186 }*/
187 //}
188
189 // Clear output array
190 if (fHits) fHits->Clear();
191 if (fClusters) fClusters->Delete(); //Clear(); // Delete because of memory escape
192
193 fuClusters = 0;
194 // --- Time-slice mode: process entire array
195
196 if (!fEventMode) {
197 //if ( fMode == kCbmTimeslice ){
198 ProcessData(nullptr);
199 LOG(info) << setw(20) << left << GetName() << ": processing time is " << timer.RealTime()
200 << " Time Slice Number is " << fNofTimeslices << " digis "
202 //<< "s digis " << fDigis->GetEntriesFast()
203 << " clusters " << fClusters->GetEntriesFast() << " total hits " << fHits->GetEntriesFast();
204 }
205 // --- Event mode: loop over events
206 else {
207 assert(fEvents);
208 Int_t nEvents = fEvents->GetEntriesFast();
209 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
210 CbmEvent* event = dynamic_cast<CbmEvent*>(fEvents->At(iEvent));
211 assert(event);
212 Int_t nDigis =
214 //: fDigis->GetEntriesFast() );
215 //if (event) LOG(debug)<<" Timeslice "<< fNofTimeslices <<" event : " << event->GetNumber() <<" nDigi : " << nDigis;
216 ProcessData(event);
217 LOG(debug) << setw(20) << left
218 << GetName()
219 //<< ": Processing Time for an event is " << timer.RealTime()
220 << ": Time slice " << fNofTimeslices << " with " << nEvents << (nEvents == 1 ? " event" : " events")
221 << " and processing event nubmer " << iEvent << " digis "
222 << nDigis
223 //<< "s digis " << fDigis->GetEntriesFast()
224 << " and created cluster " << event->GetNofData(ECbmDataType::kMuchCluster) << " and created hit "
225 << event->GetNofData(ECbmDataType::kMuchPixelHit);
226 } //# events
227 LOG(info) << setw(20) << left << GetName() << ": Processing Time is " << timer.RealTime() << ": Time slice "
228 << fNofTimeslices << " with " << nEvents << (nEvents == 1 ? " event" : " events") << "s digis "
230 //<< "s digis " << fDigis->GetEntriesFast()
231 << " and event wise total "
232 << " clusters " << fClusters->GetEntriesFast() << " total hits " << fHits->GetEntriesFast();
233
234 } //? event mode
236 /* double *pos = IgnoredAddresses.data();
237 for (int i=0;i<IgnoredAddresses.size();i++)
238 LOG(info) << "Ignored Address " << *pos++ ; */
239 LOG(info) << "Total ignored address " << IgnoredAddresses.size();
240 IgnoredAddresses.clear();
241}
242// -------------------------------------------------------------------------
243
244// ----- Public method Exec --------------------------------------------
246{
247 TStopwatch EventTimer;
248 EventTimer.Start();
249
250 fEvent++;
251 //LOG(debug3)<<" Start creating cluster ";
252 // Find clusters
253 FindClusters(event);
254 Int_t NuOfClusterInEvent = (event ? event->GetNofData(ECbmDataType::kMuchCluster) : fClusters->GetEntriesFast());
255
256 for (Int_t clusterIndex = 0; clusterIndex < NuOfClusterInEvent; ++clusterIndex) {
257 UInt_t iCluster = (event ? event->GetIndex(ECbmDataType::kMuchCluster, clusterIndex) : clusterIndex);
258 CbmMuchCluster* cluster = (CbmMuchCluster*) fClusters->At(iCluster);
259 switch (fAlgorithm) {
260 // Hit
261 case 0:
262 // One hit per pad
263 case 1: {
264 // One hit per cluster
265 CreateHits(cluster, iCluster, event);
266 break;
267 }
268 case 2: {
269 // Simple cluster deconvolution
270 ExecClusteringSimple(cluster, iCluster, event);
271 break;
272 }
273 case 3: {
274 ExecClusteringPeaks(cluster, iCluster, event);
275 break;
276 }
277 default: {
278 LOG(fatal) << GetName() << "::Exec: The algorithm index does not exist.";
279 break;
280 }
281 }
282 }
283 fDigiIndices.clear();
284 fFiredPads.clear();
285}
286// -------------------------------------------------------------------------
287
288
289// ----- Private method FindClusters ------------------------------------
291{
292
293 Int_t nDigis = (event ? event->GetNofData(ECbmDataType::kMuchDigi) : fDigiManager->GetNofDigis(ECbmModuleId::kMuch));
294 //: fDigis->GetEntriesFast() );
295 if (event)
296 LOG(debug2) << " Timeslice " << fNofTimeslices << " event : " << event->GetNumber() << " nDigi : " << nDigis;
297 if (nDigis < 0) return;
298 if (fAlgorithm == 0) {
299 for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
300 UInt_t digiIndex = (event ? event->GetIndex(ECbmDataType::kMuchDigi, iDigi) : iDigi);
301 fDigiIndices.clear();
302 fDigiIndices.push_back(digiIndex);
303 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(digiIndex);
304 const CbmMuchDigi* digi;
305 if (!bBeamTimeDigi)
306 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(digiIndex));
307 else
308 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(digiIndex));
309 //const CbmMuchDigi* digi = static_cast<const CbmMuchDigi*>(fDigis->At(digiIndex));
310 CbmMuchCluster* cluster = new ((*fClusters)[fuClusters++]) CbmMuchCluster();
311 Int_t address = CbmMuchAddress::GetAddress(
314 cluster->SetAddress(address);
315 cluster->AddDigis(fDigiIndices);
316 // --- In event-by-event mode after event building: register clusters to event using ECbmDataType::kMuchCluster
317 //Uncomment below code
318 if (event) {
319 event->AddData(ECbmDataType::kMuchCluster, fuClusters - 1);
320 } //? Event object
321 }
322 return;
323 }
324
325 vector<CbmMuchModuleGem*> modules = fGeoScheme->GetGemModules();
326
327 // Clear array of digis in the modules
328 for (UInt_t m = 0; m < modules.size(); m++)
329 modules[m]->ClearDigis();
330
331 // Fill array of digis in the modules. Digis are automatically sorted in time
332
333 for (Int_t iDigi = 0; iDigi < nDigis; iDigi++) {
334 UInt_t digiIndex = (event ? event->GetIndex(ECbmDataType::kMuchDigi, iDigi) : iDigi);
335 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(digiIndex);
336 //const auto * digi;
337 const CbmMuchDigi* digi;
338 if (!bBeamTimeDigi)
339 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(digiIndex));
340 else
341 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(digiIndex));
342 //const CbmMuchDigi* digi =static_cast<const CbmMuchDigi*>(fDigis->At(digiIndex));
343 Double_t time = digi->GetTime();
344 // Double_t chanid = digi->GetChannelId();
345 UInt_t address = digi->GetAddress();
346 // UInt_t adc = digi->GetAdc();
347 // 22/02/23 for mCBM Data checking for GEM and RPC sector and channel limits in a module.
348 //
349 if (fFlag) {
351 || (CbmMuchAddress::GetStationIndex(digi->GetAddress()) == 1)) {
354 IgnoredAddresses.push_back(digi->GetAddress());
355 continue;
356 }
359 IgnoredAddresses.push_back(digi->GetAddress());
360 continue;
361 }
362 }
364 || (CbmMuchAddress::GetStationIndex(digi->GetAddress()) == 3)) {
367 IgnoredAddresses.push_back(digi->GetAddress());
368 continue;
369 }
372 IgnoredAddresses.push_back(digi->GetAddress());
373 continue;
374 }
375 }
376 if (fGeoScheme->GetModuleByDetId(address) == NULL) {
377 IgnoredAddresses.push_back(digi->GetAddress());
378 continue;
379 }
380 }
381 fGeoScheme->GetModuleByDetId(address)->AddDigi(time, digiIndex);
382 }
383
384 // Find clusters module-by-module
385 for (UInt_t m = 0; m < modules.size(); m++) {
386 CbmMuchModuleGem* module = modules[m];
387 assert(module);
388 multimap<Double_t, Int_t> digis = modules[m]->GetDigis();
389 multimap<Double_t, Int_t>::iterator it, itmin, itmax;
390
391 // Split module digis into time slices according to fClusterSeparationTime
392 vector<multimap<Double_t, Int_t>::iterator> slices;
393 Double_t tlast = -10000;
394 // slices.push_back(digis.begin());
395 for (it = digis.begin(); it != digis.end(); ++it) {
396 Double_t t = it->first;
397 if (t > tlast + fClusterSeparationTime) slices.push_back(it);
398 tlast = t;
399 }
400 slices.push_back(it);
401 for (UInt_t s = 1; s < slices.size(); s++) {
402 fFiredPads.clear();
403 for (it = slices[s - 1]; it != slices[s]; it++) {
404 Int_t iDigi = it->second;
405 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(iDigi);
406 const CbmMuchDigi* digi;
407 //const auto * digi;
408 if (!bBeamTimeDigi)
409 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(iDigi));
410 else
411 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(iDigi));
412 //const CbmMuchDigi* digi = static_cast<const CbmMuchDigi*>(fDigis->At(iDigi));
413
414 if (module->GetPad(digi->GetAddress()) == nullptr) {
415 IgnoredAddresses.push_back(digi->GetAddress());
416 continue;
417 }
418
419 CbmMuchPad* pad = module->GetPad(digi->GetAddress());
420
421 pad->SetDigiIndex(iDigi);
422 fFiredPads.push_back(pad);
423 }
424 // Loop over fired pads in a time slice of 100 ns
425 for (UInt_t p = 0; p < fFiredPads.size(); p++) {
426 fDigiIndices.clear();
428 if (fDigiIndices.size() == 0) continue;
429 //const CbmMuchDigi* digi = static_cast<const CbmMuchDigi*>(fDigis->At(fDigiIndices.front()));
430 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(fDigiIndices.front());
431 //const auto * digi;
432 const CbmMuchDigi* digi;
433 if (!bBeamTimeDigi)
434 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(fDigiIndices.front()));
435 else
436 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(fDigiIndices.front()));
437
438 CbmMuchCluster* cluster = new ((*fClusters)[fuClusters++]) CbmMuchCluster();
439 Int_t address = CbmMuchAddress::GetAddress(
442
443 cluster->SetAddress(address);
444 cluster->AddDigis(fDigiIndices);
445 // --- In event-by-event mode after event building: register clusters to event using ECbmDataType::kMuchCluster
446 if (event) {
447 event->AddData(ECbmDataType::kMuchCluster, fuClusters - 1);
448 } //? Event object
449 }
450 }
451 }
452}
453// -------------------------------------------------------------------------
454
455
456// ----- Private method CreateCluster -----------------------------------
458{
459 Int_t digiIndex = pad->GetDigiIndex();
460 if (digiIndex < 0) return;
461 fDigiIndices.push_back(digiIndex);
462 pad->SetDigiIndex(-1);
463 vector<CbmMuchPad*> neighbours = pad->GetNeighbours();
464 for (UInt_t i = 0; i < neighbours.size(); i++)
465 CreateCluster(neighbours[i]);
466}
467// -------------------------------------------------------------------------
468
469
470// ----- Private method ExecClusteringSimple ----------------------------
472{
473 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(cluster->GetDigi(0));
474 //const auto * digi;
475 const CbmMuchDigi* digi;
476 if (!bBeamTimeDigi)
477 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(cluster->GetDigi(0)));
478 else
479 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(cluster->GetDigi(0)));
480 //CbmMuchDigi* digi = static_cast<CbmMuchDigi*>(fDigis->At(cluster->GetDigi(0)));
482 CbmMuchModuleGem* module = (CbmMuchModuleGem*) m;
483 assert(module);
484 // Int_t iStation = CbmMuchAddress::GetStationIndex(digi->GetAddress());
485
486 Int_t maxCharge = 0;
487 for (Int_t iDigi = 0; iDigi < cluster->GetNofDigis(); iDigi++) {
488 Int_t digiIndex = cluster->GetDigi(iDigi);
489 //digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(digiIndex);
490 if (!bBeamTimeDigi)
491 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(digiIndex));
492 else
493 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(digiIndex));
494 //digi = static_cast<CbmMuchDigi*> (fDigis->At(digiIndex));
495 Int_t charge = digi->GetAdc();
496 if (charge > maxCharge) maxCharge = charge;
497 }
498
499 UInt_t threshold = UInt_t(fThresholdRatio * maxCharge);
500
501 // Fire pads which are higher than threshold in a cluster
502 fFiredPads.clear();
503 for (Int_t iDigi = 0; iDigi < cluster->GetNofDigis(); iDigi++) {
504 Int_t digiIndex = cluster->GetDigi(iDigi);
505 //digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(digiIndex);
506 if (!bBeamTimeDigi)
507 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(digiIndex));
508 else
509 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(digiIndex));
510 //digi = static_cast<CbmMuchDigi*> (fDigis->At(digiIndex));
511 if (digi->GetAdc() <= threshold) continue;
512 if (module->GetPad(digi->GetAddress()) == nullptr) {
513 IgnoredAddresses.push_back(digi->GetAddress());
514 continue;
515 }
516
517 CbmMuchPad* pad = module->GetPad(digi->GetAddress());
518 pad->SetDigiIndex(digiIndex);
519 fFiredPads.push_back(pad);
520 }
521 for (UInt_t p = 0; p < fFiredPads.size(); p++) {
522 fDigiIndices.clear();
524 if (fDigiIndices.size() == 0) continue;
527 CreateHits(&cl, iCluster, event);
528 }
529}
530// -------------------------------------------------------------------------
531
532
533// -------------------------------------------------------------------------
535{
536 Int_t nDigis = cluster->GetNofDigis();
537 if (nDigis <= 2) {
538 CreateHits(cluster, iCluster, event);
539 return;
540 }
541 fClusterCharges.clear();
542 fClusterPads.clear();
543 fLocalMax.clear();
544 // for (Int_t i=0;i<fNeighbours.size();i++) fNeighbours[i].clear();
545 fNeighbours.clear();
546
547 // Fill cluster map
548 for (Int_t i = 0; i < nDigis; i++) {
549 Int_t iDigi = cluster->GetDigi(i);
550 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(iDigi);
551 //const auto * digi;
552 const CbmMuchDigi* digi;
553 if (!bBeamTimeDigi)
554 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(iDigi));
555 else
556 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(iDigi));
557 //const CbmMuchDigi* digi = static_cast<const CbmMuchDigi*>(fDigis->At(iDigi));
558 UInt_t address = digi->GetAddress();
560 assert(module);
561 if (module->GetPad(digi->GetAddress()) == nullptr) {
562 IgnoredAddresses.push_back(digi->GetAddress());
563 continue;
564 }
565 CbmMuchPad* pad = module->GetPad(address);
566 Int_t adc = digi->GetAdc();
567 fClusterPads.push_back(pad);
568 fClusterCharges.push_back(adc);
569 fLocalMax.push_back(1);
570 }
571
572 // Fill neighbours
573 for (Int_t i = 0; i < nDigis; i++) {
574 CbmMuchPad* pad = fClusterPads[i];
575 vector<CbmMuchPad*> neighbours = pad->GetNeighbours();
576 vector<Int_t> selected_neighbours;
577 for (UInt_t ip = 0; ip < neighbours.size(); ip++) {
578 CbmMuchPad* p = neighbours[ip];
579 vector<CbmMuchPad*>::iterator it = find(fClusterPads.begin(), fClusterPads.end(), p);
580 if (it == fClusterPads.end()) continue;
581 selected_neighbours.push_back(it - fClusterPads.begin());
582 }
583 fNeighbours.push_back(selected_neighbours);
584 }
585
586 // Flag local maxima
587 for (Int_t i = 0; i < nDigis; i++) {
588 Int_t c = fClusterCharges[i];
589 for (UInt_t n = 0; n < fNeighbours[i].size(); n++) {
590 Int_t in = fNeighbours[i][n];
591 Int_t cn = fClusterCharges[in];
592 if (cn < c) fLocalMax[in] = 0;
593 }
594 }
595
596 // Fire pads corresponding to local maxima
597 fFiredPads.clear();
598 for (Int_t i = 0; i < nDigis; i++) {
599 if (fLocalMax[i] == 0) continue;
600 CbmMuchPad* pad = fClusterPads[i];
601 pad->SetDigiIndex(cluster->GetDigi(i));
602 fFiredPads.push_back(pad);
603 }
604
605 // Create clusters
606 for (UInt_t p = 0; p < fFiredPads.size(); p++) {
607 fDigiIndices.clear();
609 if (fDigiIndices.size() == 0) continue;
612 CreateHits(&cl, iCluster, event);
613 }
614}
615// -------------------------------------------------------------------------
616
617
618// ----- Private method CreateHits --------------------------------------
619void CbmMuchFindHitsGem::CreateHits(CbmMuchCluster* cluster, Int_t iCluster, CbmEvent* event)
620{
621 Int_t nDigis = cluster->GetNofDigis();
622 Double_t sumq = 0, sumx = 0, sumy = 0, sumdx2 = 0, sumdy2 = 0, sumdxy2 = 0,
623 sumdt2 = 0; // , sumt =0 // not used FU 22.03.23
624 Double_t q = 0, x = 0, y = 0, t = 0, z = 0, dx = 0, dy = 0, dxy = 0, dt = 0;
625 Double_t nX = 0, nY = 0, nZ = 0;
626 Int_t address = 0;
627 Int_t planeId = 0;
628 CbmMuchModuleGem* module = NULL;
629
630 Double_t tmin = -1;
631
632 for (Int_t i = 0; i < nDigis; i++) {
633 Int_t iDigi = cluster->GetDigi(i);
634 //const CbmMuchDigi* digi = (CbmMuchDigi*) fDigiManager->Get<CbmMuchDigi>(iDigi);
635 //const auto * digi;
636 const CbmMuchDigi* digi;
637 if (!bBeamTimeDigi)
638 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchDigi>(iDigi));
639 else
640 digi = static_cast<const CbmMuchDigi*>(fDigiManager->Get<CbmMuchBeamTimeDigi>(iDigi));
641 //const CbmMuchDigi* digi = static_cast<const CbmMuchDigi*>(fDigis->At(iDigi));
642 if (i == 0) {
644 planeId = fGeoScheme->GetLayerSideNr(address);
645 // Below address should be 32 bit CbmMuchAddress therefore changing it. Need to check for simulation. 22/02/23
646 module = (CbmMuchModuleGem*) fGeoScheme->GetModuleByDetId(digi->GetAddress());
647 assert(module);
648 z = module->GetPosition()[2];
649 }
650 // One cluster must be from single module. No need of below line.
651 // CbmMuchModule* moduleDet = fGeoScheme->GetModuleByDetId(digi->GetAddress()); /// added
652
653 if (module->GetPad(digi->GetAddress()) == nullptr) {
654 IgnoredAddresses.push_back(digi->GetAddress());
655 continue;
656 }
657
658 CbmMuchPad* pad = module->GetPad(digi->GetAddress());
659 x = pad->GetX();
660 y = pad->GetY();
661
663 Double_t gemDriftTimeCorrc = 15.0; // Drift time mean for GEM is 15 ns. Drift vel =100um/ns and drift width 3mm.
664 Double_t rpcDriftTimeCorrc = 8.33; // Drift time mean for RPC is 8.33 ns. Drift vel =120um/ns and drift width 2mm.
665 Double_t gemHitTimeError = 4.0; // Hit time error for GEM = 4.0 as residual dist width is 4, to make pull width 1.
666 Double_t rpcHitTimeError =
667 2.3; // Hit time error for RPC = 2.3 ns, as residual dist width is 2.3. This is made to make pull dist width ~1
668 Double_t timeShift = 0.5; // this is added because residual time dist is shifted by -0.5 from 0.
669
670 if (module->GetDetectorType() == 3)
671 {
672 if (fFlag == 0)
673 t = digi->GetTime() - gemDriftTimeCorrc + timeShift;
674 else
675 t = digi->GetTime(); // Not correcting Drift Time for mCBM data
676 dt = gemHitTimeError;
677 }
678 if (module->GetDetectorType() == 4)
679 {
680 t = digi->GetTime() - rpcDriftTimeCorrc + timeShift;
681 dt = rpcHitTimeError;
682 }
683
684 if (tmin < 0) tmin = t;
685 if (tmin < t) tmin = t;
686 q = digi->GetAdc();
687 dx = pad->GetDx();
688 dy = pad->GetDy();
689 dxy = pad->GetDxy();
690 //dt = 4.; // digi->GetDTime(); //TODO introduce time uncertainty determination
691 sumq += q;
692 sumx += q * x;
693 sumy += q * y;
694 //sumt += q * t; // not used FU 22.03.23
695 sumdx2 += q * q * dx * dx;
696 sumdy2 += q * q * dy * dy;
697 sumdxy2 += q * q * dxy * dxy;
698 sumdt2 += q * q * dt * dt;
699 //std::cout<<" i "<<i<<" q "<<q<<" sumq "<<sumq<<std::endl;
700 }
701
702 x = sumx / sumq;
703 y = sumy / sumq;
704 // t = sumt/sumq;
705 t = tmin;
706 dx = sqrt(sumdx2 / 12) / sumq;
707 dy = sqrt(sumdy2 / 12) / sumq;
708 dxy = sqrt(sumdxy2 / 12) / sumq;
709 dt = sqrt(sumdt2) / sumq;
710 Int_t iHit = fHits->GetEntriesFast();
711
712 //------------------------------Added by A. Agarwal 29.09.2022 for mCbm ---------------------------
713 // Double_t tX = 8.5, tY = 81.0;
714 Double_t tX = 8.5, tY = 81.0; // commented for testing
715 //Double_t tX =18.5 , tY = 80.0 ;
716 if (module->GetDetectorType() == 3) // GEM
717 {
718 tX = fGemTX;
719 tY = fGemTY;
720 }
721 else if (module->GetDetectorType() == 4) // RPC
722 {
723 tX = fRpcTX;
724 tY = fRpcTY;
725 //tX = -66.5;
726 //tY = 72.0;
727 }
728 else {
729 LOG(error) << "Unknown detector type";
730 }
731
732 nX = x + tX; // for miniCBM setup set these according to the geometry
733 // tag and in line number 125 during initialization.
734 nY = y + tY; //
735 nZ = z;
736
737 if (fFlag == 1) {
738 new ((*fHits)[iHit]) CbmMuchPixelHit(address, nX, nY, nZ, dx, dy, 0, dxy, iCluster, planeId, t, dt); //mCbm
739 }
740 else {
741 new ((*fHits)[iHit]) CbmMuchPixelHit(address, x, y, z, dx, dy, 0, dxy, iCluster, planeId, t, dt); //Cbm
742 }
743 //Adding CbmMuchPixelHit entries in the CbmEvent
744 if (event) {
745 event->AddData(ECbmDataType::kMuchPixelHit, iHit);
746 } //? Event object
747}
748// ---------------------------------------------------------------------------------
749
ClassImp(CbmConverterManager)
@ kMuch
Muon detection system.
@ kMuchModule
Module.
Data container for MUCH clusters.
Class for pixel hits in MUCH detector.
friend fvec sqrt(const fvec &a)
void SetAddress(int32_t address)
Definition CbmCluster.h:94
int32_t GetDigi(int32_t index) const
Get digi at position index.
Definition CbmCluster.h:76
void AddDigis(const std::vector< int32_t > &indices)
Add array of digi to cluster.
Definition CbmCluster.h:57
int32_t GetNofDigis() const
Number of digis in cluster.
Definition CbmCluster.h:69
static Int_t GetNofDigis(ECbmModuleId systemId)
void UseMuchBeamTimeDigi(Bool_t)
Use CbmMuchBeamTimeDigi instead of CbmMuchDigi for MUCH.
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
static int32_t GetModuleIndex(int32_t address)
static int32_t GetSectorIndex(int32_t address)
static int32_t GetChannelIndex(int32_t address)
static int32_t GetLayerIndex(int32_t address)
static int32_t GetLayerSideIndex(int32_t address)
static int32_t GetElementAddress(int32_t address, int32_t level)
static uint32_t GetAddress(int32_t station=0, int32_t layer=0, int32_t side=0, int32_t module=0, int32_t sector=0, int32_t channel=0)
static int32_t GetStationIndex(int32_t address)
Data container for MUCH clusters.
uint16_t GetAdc() const
Definition CbmMuchDigi.h:90
int32_t GetAddress() const
Definition CbmMuchDigi.h:93
double GetTime() const
Definition CbmMuchDigi.h:94
CbmDigiManager * fDigiManager
std::vector< Int_t > fClusterCharges
virtual void Exec(Option_t *opt)
void CreateHits(CbmMuchCluster *cluster, Int_t iCluster, CbmEvent *event)
void ExecClusteringPeaks(CbmMuchCluster *cluster, Int_t iCluster, CbmEvent *event)
std::vector< Bool_t > fLocalMax
void ProcessData(CbmEvent *)
std::vector< Double_t > IgnoredAddresses
std::vector< CbmMuchPad * > fClusterPads
std::vector< CbmMuchPad * > fFiredPads
void CreateCluster(CbmMuchPad *pad)
TClonesArray * fClusters
CbmMuchGeoScheme * fGeoScheme
CbmMuchFindHitsGem(const char *digiFileName, Int_t flag)
void FindClusters(CbmEvent *)
std::vector< Int_t > fDigiIndices
virtual InitStatus Init()
std::vector< std::vector< Int_t > > fNeighbours
TClonesArray * fEvents
Interface to digi branch.
void ExecClusteringSimple(CbmMuchCluster *cluster, Int_t iCluster, CbmEvent *event)
Int_t GetLayerSideNr(Int_t detId) const
std::vector< CbmMuchModuleGem * > GetGemModules() const
void Init(TObjArray *stations, Int_t flag)
CbmMuchModule * GetModuleByDetId(Int_t detId) const
CbmMuchPad * GetPad(Int_t address)
Int_t GetDetectorType() const
void AddDigi(Double_t time, Int_t id)
void SetDigiIndex(Int_t iDigi)
Definition CbmMuchPad.h:45
Double_t GetX() const
Definition CbmMuchPad.h:32
std::vector< CbmMuchPad * > GetNeighbours() const
Definition CbmMuchPad.h:42
Double_t GetDy() const
Definition CbmMuchPad.h:35
Double_t GetDxy() const
Definition CbmMuchPad.h:36
Double_t GetDx() const
Definition CbmMuchPad.h:34
Int_t GetDigiIndex() const
Definition CbmMuchPad.h:38
Double_t GetY() const
Definition CbmMuchPad.h:33