CbmRoot
Loading...
Searching...
No Matches
CbmMuchDigitizerQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2020-2021 Facility for Antiproton and Ion Research in Europe, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Dominik Smith [committer], Sergey Gorbunov */
4
13
14#include "CbmMuchDigitizerQa.h"
15
16#include "CbmDigiManager.h"
17#include "CbmLink.h"
18#include "CbmMCDataArray.h"
19#include "CbmMCDataManager.h"
20#include "CbmMCTrack.h"
21#include "CbmMatch.h"
22#include "CbmMuchAddress.h"
23#include "CbmMuchDigi.h"
24#include "CbmMuchGeoScheme.h"
25#include "CbmMuchModule.h"
26#include "CbmMuchModuleGem.h"
27#include "CbmMuchPad.h"
28#include "CbmMuchPoint.h"
29#include "CbmMuchPointInfo.h"
30#include "CbmMuchRecoDefs.h"
31#include "CbmMuchStation.h"
32#include "CbmQaCanvas.h"
33#include "CbmTimeSlice.h"
34
35#include <FairRootManager.h>
36#include <FairSink.h>
37#include <FairTask.h>
38#include <Logger.h>
39
40#include "TClonesArray.h"
41#include "TDatabasePDG.h"
42#include "TF1.h"
43#include "TGraph.h"
44#include "TH1.h"
45#include "TH2.h"
46#include "TParameter.h"
47#include "TString.h"
48#include "TStyle.h"
49#include <TAxis.h>
50#include <TCanvas.h>
51#include <TDirectory.h>
52#include <TMath.h>
53#include <TParticlePDG.h>
54#include <TVector2.h>
55
56#include <iostream>
57#include <limits>
58#include <vector>
59
60#include <math.h>
61
62using std::endl;
63using std::vector;
64
65#define BINNING_CHARGE 100, 0, 300.0
66#define BINNING_LENGTH 100, 0, 2.5
67#define BINNING_CHARGE_LOG 100, 3, 10
68#define BINNING_ENERGY_LOG 100, -0.5, 4.5
69#define BINNING_ENERGY_LOG_EL 100, -0.5, 4.5
70
72
73// -------------------------------------------------------------------------
74CbmMuchDigitizerQa::CbmMuchDigitizerQa(const char* name, Int_t verbose)
75 : FairTask(name, verbose)
76 , fOutFolder("MuchDigiQA", "MuchDigitizerQA")
77 , fNevents("nEvents", 0)
78 , fvUsPadsFiredR()
79 , fvUsPadsFiredXY()
80 , fvTrackCharge()
81 , fvPadsTotalR()
82 , fvPadsFiredR()
83 , fvPadOccupancyR()
84{
85}
86
87// -------------------------------------------------------------------------
89
90// -------------------------------------------------------------------------
92{
93
94 histFolder = nullptr;
95 fGeoScheme = nullptr;
96 fDigiManager = nullptr;
97 fPoints = nullptr;
98 fDigis = nullptr;
99 fDigiMatches = nullptr;
100 fMCTracks = nullptr;
101 fOutFolder.Clear();
102 fMcPointInfoMap.clear();
103
104 for (uint i = 0; i < fvUsPadsFiredR.size(); i++) {
105 SafeDelete(fvUsPadsFiredR[i]);
106 }
107 for (uint i = 0; i < fvUsPadsFiredXY.size(); i++) {
108 SafeDelete(fvUsPadsFiredXY[i]);
109 }
110 for (uint i = 0; i < fvTrackCharge.size(); i++) {
111 SafeDelete(fvTrackCharge[i]);
112 }
113 for (uint i = 0; i < fvPadsTotalR.size(); i++) {
114 SafeDelete(fvPadsTotalR[i]);
115 }
116 for (uint i = 0; i < fvPadsFiredR.size(); i++) {
117 SafeDelete(fvPadsFiredR[i]);
118 }
119 for (uint i = 0; i < fvPadOccupancyR.size(); i++) {
120 SafeDelete(fvPadOccupancyR[i]);
121 }
122 fvUsPadsFiredR.clear();
123 fvUsPadsFiredXY.clear();
124 fvTrackCharge.clear();
125 fvPadsTotalR.clear();
126 fvPadsFiredR.clear();
127 fvPadOccupancyR.clear();
128
129 SafeDelete(fhTrackLength);
130 SafeDelete(fhTrackLengthPi);
131 SafeDelete(fhTrackLengthPr);
132 SafeDelete(fhTrackLengthEl);
137 SafeDelete(fhTrackChargeVsTrackLength);
141 SafeDelete(fhNpadsVsS);
142 SafeDelete(fCanvCharge);
143 SafeDelete(fCanvStationCharge);
144 SafeDelete(fCanvChargeVsEnergy);
145 SafeDelete(fCanvChargeVsLength);
146 SafeDelete(fCanvTrackLength);
147 SafeDelete(fCanvNpadsVsArea);
148 SafeDelete(fCanvUsPadsFiredXY);
149 SafeDelete(fCanvPadOccupancyR);
150 SafeDelete(fCanvPadsTotalR);
151 SafeDelete(fFitEl);
152 SafeDelete(fFitPi);
153 SafeDelete(fFitPr);
154
155 fNevents.SetVal(0);
156 fNstations = 0;
157 fnPadSizesX = 0;
158 fnPadSizesY = 0;
159 fPadMinLx = 0.;
160 fPadMinLy = 0.;
161 fPadMaxLx = 0.;
162 fPadMaxLy = 0.;
163}
164
165// -------------------------------------------------------------------------
167{
168 DeInit();
169 fMcPointInfoMap.clear();
170 TDirectory* oldDirectory = gDirectory;
171
172 fManager = FairRootManager::Instance();
173 if (!fManager) {
174 LOG(fatal) << "Can not find FairRootManager";
175 return kFATAL;
176 }
177
179 if (!fGeoScheme) {
180 LOG(fatal) << "Can not find Much geo scheme";
181 return kFATAL;
182 }
183
185 LOG(info) << "CbmMuchDigitizerQa: N Stations = " << fNstations;
186
188 if (!fDigiManager) {
189 LOG(fatal) << "Can not find Much digi manager";
190 return kFATAL;
191 }
193
194 fMCTracks = nullptr;
195 fPoints = nullptr;
196
197 fMcManager = dynamic_cast<CbmMCDataManager*>(fManager->GetObject("MCDataManager"));
198
199 fTimeSlice = (CbmTimeSlice*) fManager->GetObject("TimeSlice.");
200 if (!fTimeSlice) { LOG(error) << GetName() << ": No time slice found"; }
201
202 if (fMcManager) {
203 // Get MCTrack array
204 fMCTracks = fMcManager->InitBranch("MCTrack");
205
206 // Get StsPoint array
207 fPoints = fMcManager->InitBranch("MuchPoint");
208 }
209
211 LOG(info) << " CbmMuchDigitizerQa: MC data read.";
212 }
213 else {
214 LOG(info) << " CbmMuchDigitizerQa: No MC data found.";
215 fMCTracks = nullptr;
216 fPoints = nullptr;
217 }
218 histFolder = fOutFolder.AddFolder("hist", "Histogramms");
219 fNevents.SetVal(0);
220 histFolder->Add(&fNevents);
221
222 //fVerbose = 3;
223 InitCanvases();
228 InitFits();
230
231 gDirectory = oldDirectory;
232 return kSUCCESS;
233}
234
235// -------------------------------------------------------------------------
237{
238
239 fPadMinLx = std::numeric_limits<Double_t>::max();
240 fPadMinLy = std::numeric_limits<Double_t>::max();
241 fPadMaxLx = std::numeric_limits<Double_t>::min();
242 fPadMaxLy = std::numeric_limits<Double_t>::min();
243
244 if (fVerbose > 2) {
245 printf("=========================================================\n");
246 printf(" Station Nr.\t| Sectors\t| Channels\t| Pad min size\t\t| Pad max"
247 "length\t \n");
248 printf("---------------------------------------------------------\n");
249 }
250
251 // Int_t nTotSectors = 0; // not used FU 23.03.23
252 Int_t nTotChannels = 0;
253 for (Int_t iStation = 0; iStation < fNstations; iStation++) {
254 Int_t nChannels = 0;
255 Int_t nSectors = 0;
256 Double_t padMinLx = std::numeric_limits<Double_t>::max();
257 Double_t padMinLy = std::numeric_limits<Double_t>::max();
258 Double_t padMaxLx = std::numeric_limits<Double_t>::min();
259 Double_t padMaxLy = std::numeric_limits<Double_t>::min();
260 vector<CbmMuchModule*> modules = fGeoScheme->GetModules(iStation);
261 for (UInt_t im = 0; im < modules.size(); im++) {
262 CbmMuchModule* mod = modules[im];
263 if (!mod) {
264 LOG(fatal) << "module pointer is 0";
265 return -1;
266 }
267 if (mod->GetDetectorType() != 4 && mod->GetDetectorType() != 3) {
268 LOG(fatal) << "unknown detector type " << mod->GetDetectorType();
269 return -1;
270 }
271 CbmMuchModuleGem* module = dynamic_cast<CbmMuchModuleGem*>(mod);
272 if (!module) {
273 LOG(fatal) << "module is not a Gem module";
274 return -1;
275 }
276 nChannels += module->GetNPads();
277 nSectors += module->GetNSectors();
278 vector<CbmMuchPad*> pads = module->GetPads();
279 for (UInt_t ip = 0; ip < pads.size(); ip++) {
280 CbmMuchPad* pad = pads[ip];
281 if (!pad) {
282 LOG(fatal) << "pad pointer is 0";
283 return -1;
284 }
285 if (pad->GetDx() < padMinLx) padMinLx = pad->GetDx();
286 if (pad->GetDy() < padMinLy) padMinLy = pad->GetDy();
287 if (pad->GetDx() > padMaxLx) padMaxLx = pad->GetDx();
288 if (pad->GetDy() > padMaxLy) padMaxLy = pad->GetDy();
289 }
290 }
291
292 if (fPadMinLx > padMinLx) fPadMinLx = padMinLx;
293 if (fPadMinLy > padMinLy) fPadMinLy = padMinLy;
294 if (fPadMaxLx < padMaxLx) fPadMaxLx = padMaxLx;
295 if (fPadMaxLy < padMaxLy) fPadMaxLy = padMaxLy;
296 // nTotSectors += nSectors; // not used FU 23.03.23
297 nTotChannels += nChannels;
298
299 if (fVerbose > 2) {
300 printf("%i\t\t| %i\t\t| %i\t| %5.4fx%5.4f\t\t| %5.4fx%5.4f\n", iStation + 1, nSectors, nChannels, padMinLx,
301 padMinLy, padMaxLx, padMaxLy);
302 printf("%i\t\t| %i\t\t\n", iStation + 1, nChannels);
303 }
304 }
305 if (fVerbose > 2) {
306 printf("-----------------------------------------------------------\n");
307 printf(" Much total channels:\t\t| %i\t\t\n", nTotChannels);
308 printf("===========================================================\n");
309 }
310 return 0;
311}
312
313// -------------------------------------------------------------------------
315{
316
317 /***** charge canvases ****/
318 if (fMCTracks && fPoints) {
319 fCanvCharge = new CbmQaCanvas("cTrackCharge", "Charge collected per track", 2 * 800, 2 * 400);
321
322 fCanvStationCharge = new CbmQaCanvas("cTrackChargeVsStation", "Track charge per station", 2 * 800, 2 * 400);
324
325 fCanvChargeVsEnergy = new CbmQaCanvas("cTrackChargeVsEnergy", "Track charge vs particle Energy", 2 * 800, 2 * 400);
327
328 fCanvChargeVsLength = new CbmQaCanvas("cTrackChargeVsLength", "Track charge vs track length", 2 * 800, 2 * 400);
330
335 }
336
337 /***** length canvas ****/
338 if (fMCTracks && fPoints) {
339 fCanvTrackLength = new CbmQaCanvas("cTrackLength", "track length", 2 * 800, 2 * 400);
342 }
343
344 /***** pad canvases ****/
345 fCanvUsPadsFiredXY = new CbmQaCanvas("cPadsFiredXY", "Number of pads fired vs XY", 2 * 400, 2 * 400);
347
348 fCanvPadOccupancyR = new CbmQaCanvas("cPadOccupancyR", "Pad occupancy [%] vs radius", 2 * 800, 2 * 400);
350
351 fCanvPadsTotalR = new CbmQaCanvas("cPadsTotalR", "Total pads vs radius", 2 * 800, 2 * 400);
353
357
358 /***** pad canvas (MC) ****/
359 if (fMCTracks && fPoints) {
360 fCanvNpadsVsArea = new CbmQaCanvas("cNpadsVsArea", "N pads Vs Area", 2 * 800, 2 * 400);
362 }
363}
364
365// -------------------------------------------------------------------------
367{
368
369 if (!fMCTracks || !fPoints) { return; }
370
372 for (Int_t i = 0; i < fNstations; i++) {
373 fvTrackCharge[i] = new TH1F(Form("hTrackCharge%i", i + 1),
374 Form("Track charge : Station %i; Charge [10^4 e]; Count", i + 1), BINNING_CHARGE);
375 histFolder->Add(fvTrackCharge[i]);
376 }
377
378 fhTrackCharge = new TH1F("hCharge", "Charge distribution from tracks", BINNING_CHARGE);
379 fhTrackCharge->GetXaxis()->SetTitle("Charge [10^{4} electrons]");
380
381 fhTrackChargeLog = new TH1F("hChargeLog", "Charge (log.) distribution from tracks", BINNING_CHARGE_LOG);
382 fhTrackChargeLog->GetXaxis()->SetTitle("Charge [Lg(Number of electrons)]");
383
384 fhTrackChargePr_1GeV_3mm = new TH1F("hChargePr_1GeV_3mm", "Charge for 1 GeV protons", BINNING_CHARGE);
385 fhTrackChargePr_1GeV_3mm->GetXaxis()->SetTitle("Charge [10^{4} electrons]");
386
388 new TH2F("hChargeEnergyLog", "Charge vs energy (log.) for all tracks", BINNING_ENERGY_LOG, BINNING_CHARGE);
389
391 new TH2F("hChargeEnergyLogPi", "Charge vs energy (log.) for pions", BINNING_ENERGY_LOG, BINNING_CHARGE);
392
394 new TH2F("hChargeEnergyLogPr", "Charge vs energy (log.) for protons", BINNING_ENERGY_LOG, BINNING_CHARGE);
395
397 new TH2F("hChargeEnergyLogEl", "Charge vs energy (log.) for electrons", BINNING_ENERGY_LOG_EL, BINNING_CHARGE);
398
400 new TH2F("hChargeTrackLength", "Charge vs length for all tracks", BINNING_LENGTH, BINNING_CHARGE);
401
403 new TH2F("hChargeTrackLengthPi", "Charge vs length for pions", BINNING_LENGTH, BINNING_CHARGE);
404
406 new TH2F("hChargeTrackLengthPr", "Charge vs length for proton", BINNING_LENGTH, BINNING_CHARGE);
407
409 new TH2F("hChargeTrackLengthEl", "Charge vs length for electrons", BINNING_LENGTH, BINNING_CHARGE);
410 std::vector<TH2F*> vChargeHistos;
411 vChargeHistos.push_back(fhTrackChargeVsTrackEnergyLog);
412 vChargeHistos.push_back(fhTrackChargeVsTrackEnergyLogPi);
413 vChargeHistos.push_back(fhTrackChargeVsTrackEnergyLogPr);
414 vChargeHistos.push_back(fhTrackChargeVsTrackEnergyLogEl);
415 vChargeHistos.push_back(fhTrackChargeVsTrackLength);
416 vChargeHistos.push_back(fhTrackChargeVsTrackLengthPi);
417 vChargeHistos.push_back(fhTrackChargeVsTrackLengthPr);
418 vChargeHistos.push_back(fhTrackChargeVsTrackLengthEl);
419
420 for (UInt_t i = 0; i < vChargeHistos.size(); i++) {
421 TH2F* histo = vChargeHistos[i];
422 histo->SetStats(0);
423 histo->GetYaxis()->SetDecimals(kTRUE);
424 histo->GetYaxis()->SetTitleOffset(1.4);
425 histo->GetYaxis()->CenterTitle();
426 histo->GetYaxis()->SetTitle("Charge [10^{4} electrons]");
427 if (i < 4) { histo->GetXaxis()->SetTitle("Lg(E_{kin} [MeV])"); }
428 else {
429 histo->GetXaxis()->SetTitle("Track length [cm]");
430 }
431 histFolder->Add(histo);
432 }
433}
434
435// -------------------------------------------------------------------------
437{
438
439 if (!fMCTracks || !fPoints) { return; }
440
441 fhTrackLength = new TH1F("hTrackLength", "Track length", BINNING_LENGTH);
442
443 fhTrackLengthPi = new TH1F("hTrackLengthPi", "Track length for pions", BINNING_LENGTH);
444
445 fhTrackLengthPr = new TH1F("hTrackLengthPr", "Track length for protons", BINNING_LENGTH);
446
447 fhTrackLengthEl = new TH1F("hTrackLengthEl", "Track length for electrons", BINNING_LENGTH);
448
449 std::vector<TH1F*> vLengthHistos;
450 vLengthHistos.push_back(fhTrackLength);
451 vLengthHistos.push_back(fhTrackLengthPi);
452 vLengthHistos.push_back(fhTrackLengthPr);
453 vLengthHistos.push_back(fhTrackLengthEl);
454
455 for (UInt_t i = 0; i < vLengthHistos.size(); i++) {
456 TH1F* histo = vLengthHistos[i];
457 histo->GetXaxis()->SetTitle("Track length [cm]");
458 histFolder->Add(histo);
459 }
460}
461
462// -------------------------------------------------------------------------
464{
465 // non-MC
466 fvPadsTotalR.resize(fNstations);
469 fvPadsFiredR.resize(fNstations);
471
472 for (Int_t i = 0; i < fNstations; i++) {
473 CbmMuchStation* station = fGeoScheme->GetStation(i);
474 Double_t rMax = station->GetRmax();
475 Double_t rMin = station->GetRmin();
476
477 fvPadsTotalR[i] =
478 new TH1F(Form("hPadsTotal%i", i + 1), Form("Number of pads vs radius: Station %i;Radius [cm]", i + 1), 100,
479 0.6 * rMin, 1.2 * rMax);
480
481 fvUsPadsFiredR[i] =
482 new TH1F(Form("hUsPadsFired%i", i + 1), Form("Number of fired pads vs radius: Station %i;Radius [cm]", i + 1),
483 100, 0.6 * rMin, 1.2 * rMax);
484
485 fvUsPadsFiredXY[i] = new TH2F(Form("hUsPadsFiredXY%i", i + 1), Form("Pads fired XY : Station %i; X; Y", i + 1), 100,
486 -1.2 * rMax, 1.2 * rMax, 100, -1.2 * rMax, 1.2 * rMax);
487
488 fvPadsFiredR[i] =
489 new TH1F(Form("hPadsFired%i", i + 1), Form("Number of fired pads vs radius: Station %i;Radius [cm]", i + 1), 100,
490 0.6 * rMin, 1.2 * rMax);
491
492 fvPadOccupancyR[i] = new TH1F(Form("hOccupancy%i", i + 1),
493 Form("Pad occupancy vs radius: Station %i;Radius [cm];Occupancy, %%", i + 1), 100,
494 0.6 * rMin, 1.2 * rMax);
495
496 histFolder->Add(fvPadsTotalR[i]);
498 histFolder->Add(fvPadsFiredR[i]);
500 }
501
502 // MC below
503 if (fMCTracks && fPoints) {
504 fhNpadsVsS =
505 new TH2F("hNpadsVsS", "Number of fired pads per particle vs average pad area", 50, 0, 5, 15, 0.5, 15.5);
506 fhNpadsVsS->SetYTitle("N fired pads");
507 fhNpadsVsS->SetXTitle("pad area [cm^2]");
509 }
510}
511
512// -------------------------------------------------------------------------
514{
515
516 vector<CbmMuchModule*> modules = fGeoScheme->GetModules();
517 for (vector<CbmMuchModule*>::iterator im = modules.begin(); im != modules.end(); im++) {
518 CbmMuchModule* mod = (*im);
519 if (mod->GetDetectorType() == 4 || mod->GetDetectorType() == 3) { // modified for rpc
520 CbmMuchModuleGem* module = (CbmMuchModuleGem*) mod;
521 vector<CbmMuchPad*> pads = module->GetPads();
522 for (UInt_t ip = 0; ip < pads.size(); ip++) {
523 CbmMuchPad* pad = pads[ip];
524 Int_t stationId = CbmMuchAddress::GetStationIndex(pad->GetAddress());
525 Double_t x0 = pad->GetX();
526 Double_t y0 = pad->GetY();
527 Double_t r0 = TMath::Sqrt(x0 * x0 + y0 * y0);
528 fvPadsTotalR[stationId]->Fill(r0);
529 }
530 }
531 }
532
533 Double_t errors[100];
534 for (Int_t i = 0; i < 100; i++) {
535 errors[i] = 0;
536 }
537 for (Int_t i = 0; i < fNstations; i++) {
538 //fvPadsTotalR[i]->Sumw2();
539 fvPadsTotalR[i]->SetError(errors);
540 }
541}
542
543// -------------------------------------------------------------------------
545{
546
547 fFitEl = new TF1("fit_el", LandauMPV, -0.5, 4.5, 1);
548 fFitEl->SetParameter(0, 0.511);
549 fFitEl->SetLineWidth(2);
550 fFitEl->SetLineColor(kBlack);
551
552 fFitPi = new TF1("fit_pi", LandauMPV, -0.5, 4.5, 1);
553 fFitPi->SetParameter(0, 139.57);
554 fFitPi->SetLineWidth(2);
555 fFitPi->SetLineColor(kBlack);
556
557 fFitPr = new TF1("fit_pr", LandauMPV, -0.5, 4.5, 1);
558 fFitPr->SetParameter(0, 938.27);
559 fFitPr->SetLineWidth(2);
560 fFitPr->SetLineColor(kBlack);
561}
562
563// -------------------------------------------------------------------------
565{
566 // Get run and runtime database
567
568 // The code currently does not work,
569 // CbmMuchGeoScheme::Instance() must be initialised outside.
570 // - Sergey
571
572 // FairRuntimeDb* db = FairRuntimeDb::instance();
573 // if ( ! db ) Fatal("SetParContainers", "No runtime database");
574 // Get MUCH geometry parameter container
575 // CbmGeoMuchPar *fGeoPar = (CbmGeoMuchPar*)
576 // db->getContainer("CbmGeoMuchPar"); TObjArray *stations =
577 // fGeoPar->GetStations(); cout<<"\n\n SG: stations:
578 // "<<stations->GetEntriesFast()<<endl;
579 // TString geoTag;
580 // CbmSetup::Instance()->GetGeoTag(ECbmModuleId::kMuch, geoTag);
581 // bool mcbmFlag = geoTag.Contains("mcbm", TString::kIgnoreCase);
582 // CbmMuchGeoScheme::Instance()->Init(stations, mcbmFlag);
583}
584
585// -------------------------------------------------------------------------x
587{
588 fNevents.SetVal(fNevents.GetVal() + 1);
589 LOG(debug) << "Event: " << fNevents.GetVal();
590
591 if (CheckConsistency() != 0) { return; }
592
593 TDirectory* oldDirectory = gDirectory;
594
595 OccupancyQa();
596
599
601
602 if (fVerbose > 1) {
605 }
606
607 gDirectory = oldDirectory;
608 return;
609}
610
611
612// -------------------------------------------------------------------------
614{
615 // check consistency of geometry & data
616 if (!fDigiManager) {
617 LOG(error) << "Can not find Much digi manager";
618 return -1;
619 }
620 if (!fGeoScheme) {
621 LOG(error) << "Can not find Much geo scheme";
622 return -1;
623 }
624
625 for (Int_t i = 0; i < fDigiManager->GetNofDigis(ECbmModuleId::kMuch); i++) {
627 if (!digi) {
628 LOG(error) << " Much digi " << i << " out of " << fDigiManager->GetNofDigis(ECbmModuleId::kMuch)
629 << " doesn't exist";
630 return -1;
631 }
632 UInt_t address = digi->GetAddress();
633
634 int ista = CbmMuchAddress::GetStationIndex(address);
635 if (ista < 0 || ista >= fNstations) {
636 LOG(error) << " Much station " << ista << " for adress " << address << " is out of the range";
637 return -1;
638 }
639
640 CbmMuchModule* module = fGeoScheme->GetModuleByDetId(address);
641 if (!module) {
642 LOG(error) << " Much module " << address << " doesn't exist";
643 return -1;
644 }
645 if (module->GetDetectorType() != 4 && module->GetDetectorType() != 3) {
646 LOG(error) << " Much module: unknown detector type " << module->GetDetectorType();
647 return -1;
648 }
649 CbmMuchModuleGem* moduleGem = dynamic_cast<CbmMuchModuleGem*>(module);
650 if (!moduleGem) {
651 LOG(error) << " Unexpected Much module type: module " << address << " is not a Gem module";
652 return -1;
653 }
654 CbmMuchPad* pad = moduleGem->GetPad(address);
655 if (!pad) {
656 LOG(error) << " Much pad for adress " << address << " doesn't exist";
657 return -1;
658 }
659
662 if (!match) {
663 LOG(error) << " Much MC match for digi " << i << " out of " << fDigiManager->GetNofDigis(ECbmModuleId::kMuch)
664 << "doesn't exist";
665 return -1;
666 }
667 }
668 }
669 return 0;
670}
671
672// -------------------------------------------------------------------------
673const CbmMuchPad* CbmMuchDigitizerQa::GetPad(UInt_t address) const
674{
675 // get Much pad from the digi address
676 CbmMuchModuleGem* moduleGem = dynamic_cast<CbmMuchModuleGem*>(fGeoScheme->GetModuleByDetId(address));
677 CbmMuchPad* pad = moduleGem->GetPad(address);
678 return pad;
679}
680
681
682// -------------------------------------------------------------------------
684{
685 // Filling occupancy plots
686 for (Int_t i = 0; i < fDigiManager->GetNofDigis(ECbmModuleId::kMuch); i++) {
688 UInt_t address = digi->GetAddress();
689 int ista = CbmMuchAddress::GetStationIndex(address);
690 const CbmMuchPad* pad = GetPad(address);
691 Double_t x0 = pad->GetX();
692 Double_t y0 = pad->GetY();
693 double r0 = TMath::Sqrt(x0 * x0 + y0 * y0);
694 fvUsPadsFiredR[ista]->Fill(r0);
695 fvUsPadsFiredXY[ista]->Fill(x0, y0);
696 }
697}
698
699// -------------------------------------------------------------------------
701{
702
703 if (!fMCTracks || !fPoints || !fTimeSlice) {
704 LOG(debug) << " CbmMuchDigitizerQa::DigitizerQa(): Found no MC data. Skipping.";
705 return 0;
706 }
707 TVector3 vIn; // in position of the track
708 TVector3 vOut; // out position of the track
709 TVector3 p; // track momentum
710
711 fMcPointInfoMap.clear();
712
713 std::vector<CbmLink> events = fTimeSlice->GetMatch().GetLinks();
714 std::sort(events.begin(), events.end());
715
716 for (uint iLink = 0; iLink < events.size(); iLink++) {
717 CbmLink link = events[iLink];
718 int nMuchPoints = fPoints->Size(link);
719 int nMcTracks = fMCTracks->Size(link);
720 for (Int_t ip = 0; ip < nMuchPoints; ip++) {
721 link.SetIndex(ip);
722 CbmMuchPoint* point = (CbmMuchPoint*) fPoints->Get(link);
723 if (!point) {
724 LOG(error) << " Much MC point " << ip << " out of " << nMuchPoints << " doesn't exist";
725 return -1;
726 }
727
728 Int_t stId = CbmMuchAddress::GetStationIndex(point->GetDetectorID());
729
730 if (stId < 0 || stId > fNstations) {
731 LOG(error) << " Much MC point station " << stId << " is out of the range";
732 return -1;
733 }
734
735 // Check if the point corresponds to a certain MC Track
736 Int_t trackId = point->GetTrackID();
737 if (trackId < 0 || trackId >= nMcTracks) {
738 LOG(error) << " Much MC point track Id " << trackId << " is out of the range";
739 return -1;
740 }
741
742 CbmMCTrack* mcTrack = (CbmMCTrack*) fMCTracks->Get(link.GetFile(), link.GetEntry(), trackId);
743 if (!mcTrack) {
744 LOG(error) << " MC track" << trackId << " is not found";
745 return -1;
746 }
747
748 Int_t pdgCode = mcTrack->GetPdgCode();
749 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle(pdgCode);
750
751 point->PositionIn(vIn);
752 point->PositionOut(vOut);
753 point->Momentum(p);
754 Double_t length = (vOut - vIn).Mag(); // Track length
755 Double_t kine = 0; // Kinetic energy
756 if (particle) {
757 Double_t mass = particle->Mass();
758 kine = sqrt(p.Mag2() + mass * mass) - mass;
759 }
760
761 // Get mother pdg code
762 Int_t motherPdg = 0;
763 Int_t motherId = mcTrack->GetMotherId();
764 if (motherId >= nMcTracks) {
765 LOG(error) << " MC track mother Id" << trackId << " is out of the range";
766 return -1;
767 }
768 if (motherId >= 0) {
769 CbmMCTrack* motherTrack = (CbmMCTrack*) fMCTracks->Get(link.GetFile(), link.GetEntry(), motherId);
770 if (!motherTrack) {
771 LOG(error) << " MC track" << trackId << " is not found";
772 return -1;
773 }
774 motherPdg = motherTrack->GetPdgCode();
775 }
776 fMcPointInfoMap.insert({link, CbmMuchPointInfo(pdgCode, motherPdg, kine, length, stId)});
777 } // points
778 } // events
779 return 0;
780}
781
782// -------------------------------------------------------------------------
784{
785
786 if (!fMCTracks || !fPoints) {
787 LOG(debug) << " CbmMuchDigitizerQa::FillChargePerPoint()): Found no MC "
788 "data. Skipping.";
789 return;
790 }
791 for (Int_t i = 0; i < fDigiManager->GetNofDigis(ECbmModuleId::kMuch); i++) {
794 if (!match) {
795 LOG(error) << "digi has no match";
796 return;
797 }
798 const CbmMuchPad* pad = GetPad(digi->GetAddress());
799 Double_t area = pad->GetDx() * pad->GetDy();
800 Int_t nofLinks = match->GetNofLinks();
801 for (Int_t pt = 0; pt < nofLinks; pt++) {
802 const CbmLink& link = match->GetLink(pt);
803 Int_t charge = link.GetWeight();
804 CbmMuchPointInfo& info = getPointInfo(link);
805 info.AddCharge(charge);
806 info.AddArea(area);
807 }
808 }
809}
810
811// -------------------------------------------------------------------------
813{
814
815 if (!fMCTracks || !fPoints) {
816 LOG(debug) << " CbmMuchDigitizerQa::FillDigitizerPerformancePlots(): Found "
817 "no MC data. Skipping.";
818 return;
819 }
820 Bool_t verbose = (fVerbose > 2);
821 int iPoint = 0;
822 for (auto it = fMcPointInfoMap.begin(); it != fMcPointInfoMap.end(); it++, ++iPoint) {
823 const CbmMuchPointInfo& info = it->second;
824 if (verbose) {
825 printf("%i", iPoint);
826 info.Print();
827 }
828 Double_t length = info.GetLength();
829 Double_t kine = 1000 * (info.GetKine());
830 Double_t charge = info.GetCharge();
831 Int_t pdg = info.GetPdgCode();
832 if (charge <= 0) continue;
833
834 if (pdg == 22 || // photons
835 pdg == 2112) // neutrons
836 {
837 LOG(error) << "Particle with pdg code " << pdg << " produces signal in Much";
838 }
839
840 // special entry at -0.2 for the particles that are not known for TDataBasePDG
841 Double_t log_kine = (kine > 0) ? TMath::Log10(kine) : -0.2;
842 Double_t log_charge = TMath::Log10(charge);
843 charge = charge / 1.e+4; // measure charge in 10^4 electrons
844
845 Int_t nPads = info.GetNPads();
846 Double_t area = info.GetArea() / nPads;
847 // printf("%f %i\n",TMath::Log2(area),nPads);
848 //fhNpadsVsS->Fill(TMath::Log2(area), nPads);
849 fhNpadsVsS->Fill(area, nPads);
850 fhTrackCharge->Fill(charge);
851 fvTrackCharge[info.GetStationId()]->Fill(charge);
852 fhTrackChargeLog->Fill(log_charge);
853 fhTrackChargeVsTrackEnergyLog->Fill(log_kine, charge);
854 fhTrackChargeVsTrackLength->Fill(length, charge);
855 fhTrackLength->Fill(length);
856
857 if (pdg == 2212) {
858 fhTrackChargeVsTrackEnergyLogPr->Fill(log_kine, charge);
859 fhTrackChargeVsTrackLengthPr->Fill(length, charge);
860 fhTrackLengthPr->Fill(length);
861 }
862 else if (pdg == 211 || pdg == -211) {
863 fhTrackChargeVsTrackEnergyLogPi->Fill(log_kine, charge);
864 fhTrackChargeVsTrackLengthPi->Fill(length, charge);
865 fhTrackLengthPi->Fill(length);
866 }
867 else if (pdg == 11 || pdg == -11) {
868 fhTrackChargeVsTrackEnergyLogEl->Fill(log_kine, charge);
869 fhTrackChargeVsTrackLengthEl->Fill(length, charge);
870 fhTrackLengthEl->Fill(length);
871 }
872 if (pdg == 2212 && length > 0.3 && length < 0.32 && kine > 1000 && kine < 1020)
873 fhTrackChargePr_1GeV_3mm->Fill(charge);
874 }
875}
876
877
878// -------------------------------------------------------------------------
880{
881
882 if (!fMCTracks || !fPoints) { return; }
883
884 for (Int_t i = 0; i < fNstations; i++) {
885 fCanvStationCharge->cd(i + 1);
886 fvTrackCharge[i]->DrawCopy("", "");
887 }
888 fCanvCharge->cd(1);
889 fhTrackCharge->DrawCopy("", "");
890 fCanvCharge->cd(2);
891 fhTrackChargeLog->DrawCopy("", "");
892 fCanvCharge->cd(3);
893 fhTrackChargePr_1GeV_3mm->DrawCopy("", "");
894
895 for (Int_t i = 0; i < 4; i++) {
896 fCanvChargeVsEnergy->cd(i + 1);
897 gPad->Range(0, 0, 200, 200);
898 gPad->SetBottomMargin(0.11);
899 gPad->SetRightMargin(0.14);
900 gPad->SetLeftMargin(0.12);
901 gPad->SetLogz();
902 }
903 fCanvChargeVsEnergy->cd(1);
904 fhTrackChargeVsTrackEnergyLog->DrawCopy("colz", "");
905 fCanvChargeVsEnergy->cd(2);
906 fhTrackChargeVsTrackEnergyLogPi->DrawCopy("colz", "");
907 fFitPi->DrawCopy("same");
908 fCanvChargeVsEnergy->cd(3);
909 fhTrackChargeVsTrackEnergyLogPr->DrawCopy("colz", "");
910 fFitPr->DrawCopy("same");
911 fCanvChargeVsEnergy->cd(4);
912 fhTrackChargeVsTrackEnergyLogEl->DrawCopy("colz", "");
913 fFitEl->DrawCopy("same");
914
915 for (Int_t i = 0; i < 4; i++) {
916 fCanvChargeVsLength->cd(i + 1);
917 gPad->Range(0, 0, 200, 200);
918 gPad->SetBottomMargin(0.11);
919 gPad->SetRightMargin(0.14);
920 gPad->SetLeftMargin(0.12);
921 gPad->SetLogz();
922 }
923 fCanvChargeVsLength->cd(1);
924 fhTrackChargeVsTrackLength->DrawCopy("colz", "");
925 fCanvChargeVsLength->cd(2);
926 fhTrackChargeVsTrackLengthPi->DrawCopy("colz", "");
927 fCanvChargeVsLength->cd(3);
928 fhTrackChargeVsTrackLengthPr->DrawCopy("colz", "");
929 fCanvChargeVsLength->cd(4);
930 fhTrackChargeVsTrackLengthEl->DrawCopy("colz", "");
931}
932
933// -------------------------------------------------------------------------
935{
936 //non-MC
937 for (Int_t i = 0; i < fNstations; i++) {
939 //fvPadsFiredR[i]->Sumw2();
940 fvPadsFiredR[i]->Scale(1. / fNevents.GetVal());
941 fvPadOccupancyR[i]->Divide(fvPadsFiredR[i], fvPadsTotalR[i]);
942 fvPadOccupancyR[i]->Scale(100.);
943
944 fCanvPadOccupancyR->cd(i + 1);
945 fvPadOccupancyR[i]->DrawCopy("", "");
946 fCanvPadsTotalR->cd(i + 1);
947 fvPadsTotalR[i]->DrawCopy("", "");
948 fCanvUsPadsFiredXY->cd(i + 1);
949 fvUsPadsFiredXY[i]->DrawCopy("colz", "");
950 }
951 //MC below
952 if (fMCTracks && fPoints) {
953 fCanvNpadsVsArea->cd();
954 fhNpadsVsS->DrawCopy("colz", "");
955 }
956}
957
958// -------------------------------------------------------------------------
960{
961
962 if (!fMCTracks || !fPoints) { return; }
963
964 for (Int_t i = 0; i < 4; i++) {
965 fCanvTrackLength->cd(i + 1);
966 gPad->SetLogy();
967 gStyle->SetOptStat(1110);
968 }
969 fCanvTrackLength->cd(1);
970 fhTrackLength->DrawCopy("", "");
971 fCanvTrackLength->cd(2);
972 fhTrackLengthPi->DrawCopy("", "");
973 fCanvTrackLength->cd(3);
974 fhTrackLengthPr->DrawCopy("", "");
975 fCanvTrackLength->cd(4);
976 fhTrackLengthEl->DrawCopy("", "");
977}
978
979// -------------------------------------------------------------------------
981{
982
983 if (!fMCTracks || !fTimeSlice) { return; }
984
985 const CbmMatch& match = fTimeSlice->GetMatch();
986 for (int iLink = 0; iLink < match.GetNofLinks(); iLink++) {
987 CbmLink link = match.GetLink(iLink);
988 int nMuchPoints = fPoints->Size(link);
989 for (Int_t ip = 0; ip < nMuchPoints; ip++) {
990 link.SetIndex(ip);
991 CbmMuchPoint* point = (CbmMuchPoint*) fPoints->Get(link);
992 if (!point) {
993 LOG(error) << " Much MC point " << ip << " out of " << nMuchPoints << " doesn't exist";
994 return;
995 }
996 Int_t stId = CbmMuchAddress::GetStationIndex(point->GetDetectorID());
997 if (stId != 0) continue;
998 Int_t layerId = CbmMuchAddress::GetLayerIndex(point->GetDetectorID());
999 if (layerId != 0) continue;
1000 printf("point %4i xin=%6.2f yin=%6.2f xout=%6.2f yout=%6.2f zin=%6.2f\n", ip, point->GetXIn(), point->GetYIn(),
1001 point->GetXOut(), point->GetYOut(), point->GetZIn());
1002 }
1003 }
1004}
1005
1006// -------------------------------------------------------------------------
1008{
1009 for (Int_t i = 0; i < fDigiManager->GetNofDigis(ECbmModuleId::kMuch); i++) {
1011 UInt_t address = digi->GetAddress();
1012 Int_t stId = CbmMuchAddress::GetStationIndex(address);
1013 if (stId != 0) continue;
1014 Int_t layerId = CbmMuchAddress::GetLayerIndex(address);
1015 if (layerId != 0) continue;
1017 if (!module) continue;
1018 CbmMuchPad* pad = module->GetPad(address);
1019 Double_t x0 = pad->GetX();
1020 Double_t y0 = pad->GetY();
1021 UInt_t charge = digi->GetAdc();
1022 printf("digi %4i x0=%5.1f y0=%5.1f charge=%3i\n", i, x0, y0, charge);
1023 }
1024}
1025
1026// -------------------------------------------------------------------------
1028{
1029 TDirectory* oldDirectory = gDirectory;
1033 gDirectory = oldDirectory;
1034 return fOutFolder;
1035}
1036
1037// -------------------------------------------------------------------------
1039{
1040
1041 if (!FairRootManager::Instance() || !FairRootManager::Instance()->GetSink()) {
1042 LOG(error) << "No sink found";
1043 return;
1044 }
1045 FairSink* sink = FairRootManager::Instance()->GetSink();
1046 sink->WriteObject(&GetQa(), nullptr);
1047}
1048
1049// -------------------------------------------------------------------------
1051{
1052 TCanvas* c = new TCanvas("nMeanVsS", "nMeanVsS", 2 * 800, 2 * 400);
1053 printf("===================================\n");
1054 printf("DigitizerQa:\n");
1055
1056 Double_t nMean[11];
1057 Double_t s[11];
1058 for (Int_t iS = 1; iS <= 10; iS++) {
1059 nMean[iS] = 0;
1060 s[iS] = -5.25 + 0.5 * iS;
1061 Double_t total = 0;
1062 for (Int_t iN = 1; iN <= 10; iN++) {
1063 nMean[iS] += iN * fhNpadsVsS->GetBinContent(iS, iN);
1064 total += fhNpadsVsS->GetBinContent(iS, iN);
1065 }
1066 if (total > 0) nMean[iS] /= total;
1067 printf("%f %f\n", s[iS], nMean[iS]);
1068 }
1069 c->cd();
1070
1071 TGraph* gNvsS = new TGraph(11, s, nMean);
1072 //gNvsS->SetLineColor(2);
1073 //gNvsS->SetLineWidth(4);
1074 gNvsS->SetMarkerColor(4);
1075 gNvsS->SetMarkerSize(1.5);
1076 gNvsS->SetMarkerStyle(21);
1077 gNvsS->SetTitle("nMeanVsS");
1078 gNvsS->GetYaxis()->SetTitle("nMean");
1079 gNvsS->GetYaxis()->SetTitle("nMean");
1080 //gNvsS->DrawClone("ALP");
1081 gNvsS->DrawClone("AP");
1082 fOutFolder.Add(c);
1083}
1084
1085// -------------------------------------------------------------------------
1086Double_t CbmMuchDigitizerQa::LandauMPV(Double_t* lg_x, Double_t* par)
1087{
1088 Double_t gaz_gain_mean = 1.7e+4;
1089 Double_t scale = 1.e+6;
1090 gaz_gain_mean /= scale;
1091 Double_t mass = par[0]; // mass in MeV
1092 Double_t x = TMath::Power(10, lg_x[0]);
1093 return gaz_gain_mean * MPV_n_e(x, mass);
1094}
1095
1096// -------------------------------------------------------------------------
1097Double_t CbmMuchDigitizerQa::MPV_n_e(Double_t Tkin, Double_t mass)
1098{
1099 Double_t logT;
1100 TF1 fPol6("fPol6", "pol6", -5, 10);
1101 if (mass < 0.1) {
1102 logT = log(Tkin * 0.511 / mass);
1103 if (logT > 9.21034) logT = 9.21034;
1104 if (logT < min_logT_e) logT = min_logT_e;
1105 return fPol6.EvalPar(&logT, mpv_e);
1106 }
1107 else if (mass >= 0.1 && mass < 0.2) {
1108 logT = log(Tkin * 105.658 / mass);
1109 if (logT > 9.21034) logT = 9.21034;
1110 if (logT < min_logT_mu) logT = min_logT_mu;
1111 return fPol6.EvalPar(&logT, mpv_mu);
1112 }
1113 else {
1114 logT = log(Tkin * 938.272 / mass);
1115 if (logT > 9.21034) logT = 9.21034;
1116 if (logT < min_logT_p) logT = min_logT_p;
1117 return fPol6.EvalPar(&logT, mpv_p);
1118 }
1119}
@ kMuch
Muon detection system.
#define BINNING_ENERGY_LOG
#define BINNING_LENGTH
ClassImp(CbmMuchDigitizerQa)
#define BINNING_CHARGE_LOG
#define BINNING_ENERGY_LOG_EL
#define BINNING_CHARGE
constexpr double min_logT_e
constexpr double mpv_e[]
constexpr double mpv_p[]
constexpr double min_logT_p
constexpr double mpv_mu[]
constexpr double min_logT_mu
Definition of the CbmQaCanvas class.
friend fvec sqrt(const fvec &a)
friend fvec log(const fvec &a)
static Int_t GetNofDigis(ECbmModuleId systemId)
static Bool_t IsMatchPresent(ECbmModuleId systemId)
Presence of a digi match branch.
InitStatus Init()
Initialisation.
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
const CbmMatch * GetMatch(ECbmModuleId systemId, UInt_t index) const
Get a match object.
TObject * Get(const CbmLink *lnk)
Int_t Size(Int_t fileNumber, Int_t eventNumber)
Task class creating and managing CbmMCDataArray objects.
CbmMCDataArray * InitBranch(const char *name)
int32_t GetMotherId() const
Definition CbmMCTrack.h:69
int32_t GetPdgCode() const
Definition CbmMCTrack.h:68
const CbmLink & GetLink(int32_t i) const
Definition CbmMatch.h:39
int32_t GetNofLinks() const
Definition CbmMatch.h:42
const std::vector< CbmLink > & GetLinks() const
Definition CbmMatch.h:40
static int32_t GetLayerIndex(int32_t address)
static int32_t GetStationIndex(int32_t address)
uint16_t GetAdc() const
Definition CbmMuchDigi.h:90
int32_t GetAddress() const
Definition CbmMuchDigi.h:93
virtual void SetParContainers()
CbmTimeSlice * fTimeSlice
TH1F * fhTrackChargeLog
MC point charge.
std::vector< TH1F * > fvUsPadsFiredR
number of processed events
CbmQaCanvas * fCanvStationCharge
CbmQaCanvas * fCanvChargeVsEnergy
static Double_t LandauMPV(Double_t *x, Double_t *par)
CbmMCDataManager * fMcManager
const CbmMuchPad * GetPad(UInt_t address) const
get pad from the digi address
std::vector< TH2F * > fvUsPadsFiredXY
FairRootManager * fManager
folder wich contains histogramms
CbmQaCanvas * fCanvPadOccupancyR
CbmQaCanvas * fCanvPadsTotalR
virtual void Exec(Option_t *option)
std::vector< TH1F * > fvPadsFiredR
CbmMuchGeoScheme * fGeoScheme
CbmMuchPointInfo & getPointInfo(const CbmLink &link)
map point link -> point info
CbmMuchDigitizerQa(const char *name="MuchHitFinderQa", Int_t verbose=1)
TParameter< int > fNevents
output folder with histos and canvases
CbmQaCanvas * fCanvChargeVsLength
CbmQaCanvas * fCanvNpadsVsArea
CbmQaCanvas * fCanvTrackLength
CbmMCDataArray * fMCTracks
TH1F * fhTrackLength
MC point charge for selected protons.
CbmDigiManager * fDigiManager
std::vector< TH1F * > fvPadOccupancyR
static Double_t MPV_n_e(Double_t Tkin, Double_t mass)
std::map< CbmLink, CbmMuchPointInfo > fMcPointInfoMap
CbmMCDataArray * fPoints
TH1F * fhTrackChargePr_1GeV_3mm
MC point charge log scale.
CbmQaCanvas * fCanvUsPadsFiredXY
virtual InitStatus Init()
TClonesArray * fDigiMatches
std::vector< TH1F * > fvPadsTotalR
std::vector< TH1F * > fvTrackCharge
CbmMuchStation * GetStation(Int_t iStation) const
std::vector< CbmMuchModule * > GetModules() const
static CbmMuchGeoScheme * Instance()
Int_t GetNStations() const
CbmMuchModule * GetModuleByDetId(Int_t detId) const
CbmMuchPad * GetPad(Int_t address)
Int_t GetDetectorType() const
Double_t GetX() const
Definition CbmMuchPad.h:32
Double_t GetDy() const
Definition CbmMuchPad.h:35
Int_t GetAddress() const
Definition CbmMuchPad.h:31
Double_t GetDx() const
Definition CbmMuchPad.h:34
Double_t GetY() const
Definition CbmMuchPad.h:33
Double_t GetArea() const
void Print(Option_t *="") const
Double_t GetLength() const
Int_t GetCharge() const
Int_t GetPdgCode() const
void AddCharge(Int_t charge)
Int_t GetNPads() const
Double_t GetKine() const
Int_t GetStationId() const
void AddArea(Double_t s)
double GetYOut() const
double GetYIn() const
void PositionIn(TVector3 &pos) const
double GetZIn() const
void PositionOut(TVector3 &pos) const
double GetXOut() const
double GetXIn() const
Double_t GetRmin() const
Double_t GetRmax() const
void Divide2D(int nPads)
Divide canvas into nPads in 2D in a nice way.
Bookkeeping of time-slice content.
const CbmMatch & GetMatch() const