BmnRoot
Loading...
Searching...
No Matches
BmnQaOffline.cxx
Go to the documentation of this file.
1#include "BmnQaOffline.h"
2#include "BmnGlobalTrack.h"
3#include <UniRun.h>
4
5Int_t BmnQaOffline::fCurrentEvent = 0;
6
8: fBmnHeader(nullptr),
9 fDstHeader(nullptr),
10 fChainDst(nullptr),
11 fSiliconHits(nullptr),
12 fSiliconTracks(nullptr),
13 fGemHits(nullptr),
14 fGemTracks(nullptr),
15 fCscHits(nullptr),
16 fTof400Hits(nullptr),
17 fTof700Hits(nullptr),
18 fDchTracks(nullptr),
19 fMwpcTracks(nullptr),
20 fDchHits(nullptr),
21 fMwpcHits(nullptr),
22 fVertex(nullptr),
23 fGlobalTracks(nullptr),
24 period(0),
25 fRunId(0),
26 fSteering(new BmnOfflineQaSteering()),
27 fDetGem(nullptr),
28 fDetCsc(nullptr),
29 fDetSilicon(nullptr)
30{
31 isDstRead = ReadDstTree(file);
32}
33
34Bool_t BmnQaOffline::ReadDstTree(TString fileDst) {
35 fChainDst = new TChain("bmndata");
36 fChainDst->Add(fileDst.Data());
37
38 // Get dst header ...
39 fChainDst->SetBranchAddress("DstEventHeader.", &fDstHeader);
40 fChainDst->GetEntry(0);
41
42 fChainDst->SetBranchAddress("BmnDchHit", &fDchHits);
43 fChainDst->SetBranchAddress("BmnMwpcHit", &fMwpcHits);
44 fChainDst->SetBranchAddress("BmnTof400Hit", &fTof400Hits);
45 fChainDst->SetBranchAddress("BmnTof700Hit", &fTof700Hits);
46 fChainDst->SetBranchAddress("BmnCSCHit", &fCscHits);
47 fChainDst->SetBranchAddress("BmnSiliconHit", &fSiliconHits);
48 fChainDst->SetBranchAddress("BmnSiliconTrack", &fSiliconTracks);
49 fChainDst->SetBranchAddress("BmnGemStripHit", &fGemHits);
50 fChainDst->SetBranchAddress("BmnGemTrack", &fGemTracks);
51 fChainDst->SetBranchAddress("BmnDchTrack", &fDchTracks);
52 fChainDst->SetBranchAddress("BmnMwpcTrack", &fMwpcTracks);
53 fChainDst->SetBranchAddress("BmnGlobalTrack", &fGlobalTracks);
54 fChainDst->SetBranchAddress("BmnVertex", &fVertex);
55
56 if (fileDst.IsNull() || fChainDst->GetEntries() == 0)
57 return kFALSE;
58 else
59 return kTRUE;
60}
61
62InitStatus BmnQaOffline::Init() {
63 cout << " BmnQaOffline::Init() " << endl;
64
65 ioman = FairRootManager::Instance();
66
67 fBmnHeader = static_cast<BmnEventHeader*>(ioman->GetObject("BmnEventHeader."));
68 if (fBmnHeader){
69 printf("BmnEventHeader. found\n");
70 fRunId = fBmnHeader->GetRunId();
71 } else
72 {
73 fDstHeader = static_cast<DstEventHeader*>(ioman->GetObject("DstEventHeader."));
74 printf("DstEventHeader. found\n");
75 fRunId = fDstHeader->GetRunId();
76
77 }
78 // Get general info on period and exp. setup, detectors and trigger been
79 period = fSteering->GetRunAndSetupByRunId(fRunId).first;
80 setup = fSteering->GetRunAndSetupByRunId(fRunId).second;
81 nDets = fSteering->GetDetectors(period, setup).size(); // Number of detectors should be set for a certain run / setup extension
82 Int_t nTrigs = fSteering->GetTriggers(period, setup).size(); // Number of triggers should be set for a certain run / setup extension
83
84 prefix = TString::Format("RUN%d_SETUP_%s_", period, setup.Data());
85
86 // Initialize some detector geometries ...
87 fSteering->SetGeometriesByRunId(fRunId, fDetGem, fDetSilicon, fDetCsc);
88
89 // Read input arrays with det. and trig. info ...
90 DETECTORS = new TClonesArray*[nDets];
91 TRIGGERS = new TClonesArray*[nTrigs];
92
93 for (Int_t iDet = 0; iDet < nDets; iDet++)
94 DETECTORS[iDet] = (TClonesArray*) ioman->GetObject(fSteering->GetDetectors(period, setup)[iDet].Data());
95
96 for (Int_t iTrigger = 0; iTrigger < nTrigs; iTrigger++)
97 TRIGGERS[iTrigger] = (TClonesArray*) ioman->GetObject(fSteering->GetTriggers(period, setup)[iTrigger].Data());
98
99 nCoordinate = fSteering->GetNumberOfDets(period, setup, "coordinate");
100 nTime = fSteering->GetNumberOfDets(period, setup, "time");
101 nCalorimeter = fSteering->GetNumberOfDets(period, setup, "calorimeter");
102
103 coordinate = new BmnCoordinateDetQa*[nCoordinate];
104 for (Int_t iDet = 0; iDet < nCoordinate; iDet++)
105 coordinate[iDet] = new BmnCoordinateDetQa(fSteering->GetDetectors(period, setup)[iDet], fRunId);
106
107 time = new BmnTimeDetQa*[nTime];
108 for (Int_t iDet = nCoordinate; iDet < nCoordinate + nTime; iDet++)
109 time[iDet - nCoordinate] = new BmnTimeDetQa(fSteering->GetDetectors(period, setup)[iDet], fRunId);
110
111 calorimeter = new BmnCalorimeterDetQa*[nCalorimeter];
112 for (Int_t iDet = nCoordinate + nTime; iDet < nCoordinate + nTime + nCalorimeter; iDet++)
113 calorimeter[iDet - nCoordinate - nTime] = new BmnCalorimeterDetQa(fSteering->GetDetectors(period, setup)[iDet], fRunId);
114
115 for (Int_t iDet = 0; iDet < nTrigs; iDet++)
116 fTrigCorr[TRIGGERS[iDet]] = fSteering->GetTriggers(period, setup)[iDet];
117
118 if (fTrigCorr.size() != 0)
119 triggers = new BmnTrigDetQa(fTrigCorr, fRunId);
120
121 // Dst
122 dst = new BmnDstQa(fRunId);
123
124 Char_t* geoFileName = (Char_t*) "current_geo_file.root";
125 Int_t res_code = UniRun::ReadGeometryFile(period, fRunId, geoFileName);
126 if (res_code != 0) {
127 cout << "Geometry file can't be read from the database" << endl;
128 exit(-1);
129 }
130 TGeoManager::Import(geoFileName);
131 UniRun* db = UniRun::GetRun(period, fRunId);
132 isField = (*db->GetFieldVoltage() > 10.) ? kTRUE : kFALSE;
133
134 return kSUCCESS;
135}
136
138 vector <TString> detectors;
139
140 // TRIGGERS are considered as a whole detector when saving to output file
141 for (Int_t iDet = 0; iDet < nDets; iDet++)
142 detectors.push_back(fSteering->GetDetectors(period, setup)[iDet]);
143 detectors.push_back("TRIGGERS");
144
145 // Fill a vector with histo managers
146 vector <BmnQaHistoManager*> managers;
147 for (Int_t iDet = 0; iDet < nCoordinate; iDet++)
148 managers.push_back(coordinate[iDet]->GetManager());
149
150 for (Int_t iDet = 0; iDet < nTime; iDet++)
151 managers.push_back(time[iDet]->GetManager());
152
153 for (Int_t iDet = 0; iDet < nCalorimeter; iDet++)
154 managers.push_back(calorimeter[iDet]->GetManager());
155
156 managers.push_back(triggers->GetManager());
157
158 TDirectory** directories = new TDirectory*[nDets + 1];
159 for (Int_t iDet = 0; iDet < nDets + 1; iDet++) {
160 directories[iDet] = ioman->GetOutFile()->mkdir(detectors[iDet].Data());
161 directories[iDet]->cd();
162 managers[iDet]->WriteToFile();
163 }
164
165 ioman->GetOutFile()->mkdir("DST")->cd();
166 dst->GetManager()->WriteToFile();
167
168 // Delete detector classes and its histo classes
169 for (Int_t iDet = 0; iDet < nCoordinate; iDet++)
170 delete coordinate[iDet];
171
172 for (Int_t iDet = 0; iDet < nTime; iDet++)
173 delete time[iDet];
174
175 for (Int_t iDet = 0; iDet < nCalorimeter; iDet++)
176 delete calorimeter[iDet];
177
178 delete triggers;
179 delete dst;
180}
181
182void BmnQaOffline::Exec(Option_t* opt) {
183 if (isDstRead) {
184 if (fCurrentEvent > fChainDst->GetEntries() - 1)
185 return;
186 else
187 fChainDst->GetEntry(fCurrentEvent);
188 }
189 fCurrentEvent++;
190 if (fCurrentEvent % 1000 == 0)
191 cout << "Event# = " << fCurrentEvent << endl;
192
193 if (isDstRead) {
194 vector <TString> dets = {"MWPC", "SILICON", "GEM", "CSC", "TOF400", "TOF700", "DCH"};
195 vector <TClonesArray*> hits = {fMwpcHits, fSiliconHits, fGemHits, fCscHits, fTof400Hits, fTof700Hits, fDchHits};
196
197 for (size_t iDet = 0; iDet < dets.size(); iDet++)
198 for (size_t iHit = 0; iHit < dets.size(); iHit++) {
199 if (iDet != iHit)
200 continue;
201
202 doHitsDistributions(hits[iHit], dst, dets[iDet]);
203 }
204
205 GetBasicTrackDistributions();
206 doMatchingAnal();
207 doResidualsPullsAnal();
208 doAverageStripValuePerTrack();
209
210 for (auto hit : hits) {
211 doOccupancy(hit);
212 doAverageStripValuePerHit(hit);
213 }
214 }
215
216 // Coord. dets
217 for (Int_t iDet = 0; iDet < nCoordinate; iDet++) {
218 TString detName = fSteering->GetDetectors(period, setup)[iDet];
219 TClonesArray* detData = DETECTORS[iDet];
220 BmnCoordinateDetQa* detQa = coordinate[iDet];
221
222 if (detName.Contains("GEM")) {
223 GetDistributionOfFiredStrips <BmnGemStripDigit> (detData, detQa, prefix + detName); // histos 1
224 GetDistributionOfFiredStripsVsSignal <BmnGemStripDigit> (detData, detQa, prefix + detName); // histos 2
225 } else if (detName.Contains("SILICON")) {
226 GetDistributionOfFiredStrips <BmnSiliconDigit> (detData, detQa, prefix + detName);
227 } else if (detName.Contains("CSC")) {
228 GetDistributionOfFiredStrips <BmnCSCDigit> (detData, detQa, prefix + detName);
229 }
230 }
231
232 // Time dets
233 for (Int_t iDet = nCoordinate; iDet < nCoordinate + nTime; iDet++) {
234 TString detName = fSteering->GetDetectors(period, setup)[iDet];
235 TClonesArray* detData = DETECTORS[iDet];
236 BmnTimeDetQa* detQa = time[iDet - nCoordinate];
237
238 if (detName.Contains("TOF400")) {
239 GetCommonInfo <BmnTof1Digit> (detData, detQa, prefix + detName);
240 GetTofInfo <BmnTof1Digit> (detData, detQa, prefix + detName);
241 } else if (detName.Contains("TOF700")) {
242 GetCommonInfo <BmnTof2Digit> (detData, detQa, prefix + detName);
243 GetTofInfo <BmnTof2Digit> (detData, detQa, prefix + detName);
244 } else if (detName.Contains("DCH")) {
245 GetCommonInfo <BmnDchDigit> (detData, detQa, prefix + detName);
246 GetMwpcDchInfo <BmnDchDigit> (detData, detQa, prefix + detName);
247 } else if (detName.Contains("MWPC")) {
248 GetCommonInfo <BmnMwpcDigit> (detData, detQa, prefix + detName);
249 GetMwpcDchInfo <BmnMwpcDigit> (detData, detQa, prefix + detName);
250 }
251 }
252
253 // Calorim. dets
254 for (Int_t iDet = nCoordinate + nTime; iDet < nCoordinate + nTime + nCalorimeter; iDet++) {
255 TString detName = fSteering->GetDetectors(period, setup)[iDet];
256 TClonesArray* detData = DETECTORS[iDet];
257 BmnCalorimeterDetQa* detQa = calorimeter[iDet - nCoordinate - nTime];
258
259 if (detName.Contains("ECAL")) {
260 GetCommonInfo <BmnECALDigit> (detData, detQa, prefix + detName);
261 } else if (detName.Contains("ZDC")) {
262 GetCommonInfo <BmnZDCDigit> (detData, detQa, prefix + detName);
263 }
264 }
265
266 // Trig. dets
267 for (auto it : fTrigCorr)
268 GetCommonInfo <BmnTrigDigit> (it.first, triggers, it.second);
269
270 if (fCurrentEvent == fChainDst->GetEntries() - 1) {
271 doEfficiencyAnal();
272 doOccupancyAnal();
273 }
274}
275
276// Functions to be used when getting previously filled histos
277// Coordinate detectors
278
279// Histos 1
280
281template <class T> void BmnQaOffline::GetDistributionOfFiredStrips(TClonesArray* digiArray, BmnCoordinateDetQa* detHistoClass, TString detName) {
282 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
283 T* dig = (T*) digiArray->UncheckedAt(iDig);
284 Int_t module = dig->GetModule();
285 Int_t station = dig->GetStation();
286 Int_t layer = dig->GetStripLayer();
287 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of fired strips, Stat %d Mod %d Lay %d", detName.Data(), station, module, layer))->Fill(dig->GetStripNumber());
288 }
289}
290
291// Histos 2
292
293template <class T> void BmnQaOffline::GetDistributionOfFiredStripsVsSignal(TClonesArray* digiArray, BmnCoordinateDetQa* detHistoClass, TString detName) {
294 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
295 T* dig = (T*) digiArray->UncheckedAt(iDig);
296 Int_t module = dig->GetModule();
297 Int_t station = dig->GetStation();
298 Int_t layer = dig->GetStripLayer();
299 detHistoClass->GetManager()->H2(TString::Format("%s_2d, Distribution of fired strips vs. signal, Stat %d Mod %d Lay %d", detName.Data(), station, module, layer))->Fill(dig->GetStripNumber(), dig->GetStripSignal());
300 }
301}
302
303// Time detectors
304
305template <class T> void BmnQaOffline::GetCommonInfo(TClonesArray* digiArray, BmnTimeDetQa* detHistoClass, TString detName) {
306 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
307 T* dig = (T*) digiArray->UncheckedAt(iDig);
308 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of planes", detName.Data()))->Fill(dig->GetPlane());
309 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of times", detName.Data()))->Fill(dig->GetTime());
310 }
311}
312
313template <class T> void BmnQaOffline::GetMwpcDchInfo(TClonesArray* digiArray, BmnTimeDetQa* detHistoClass, TString detName) {
314 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
315 T* dig = (T*) digiArray->UncheckedAt(iDig);
316 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of wires", detName.Data()))->Fill(dig->GetWireNumber());
317 }
318}
319
320template <class T> void BmnQaOffline::GetTofInfo(TClonesArray* digiArray, BmnTimeDetQa* detHistoClass, TString detName) {
321 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
322 T* dig = (T*) digiArray->UncheckedAt(iDig);
323 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of strips", detName.Data()))->Fill(dig->GetStrip());
324 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of amplitudes", detName.Data()))->Fill(dig->GetAmplitude());
325 }
326}
327
328// Calorim. detectors
329
330template <class T> void BmnQaOffline::GetCommonInfo(TClonesArray* digiArray, BmnCalorimeterDetQa* detHistoClass, TString detName) {
331 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
332 T* dig = (T*) digiArray->UncheckedAt(iDig);
333 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of iX", detName.Data()))->Fill(dig->GetIX());
334 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of iY", detName.Data()))->Fill(dig->GetIY());
335 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of X", detName.Data()))->Fill(dig->GetX());
336 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of Y", detName.Data()))->Fill(dig->GetY());
337 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of channels", detName.Data()))->Fill(dig->GetChannel());
338 detHistoClass->GetManager()->H1(TString::Format("%s_1d, Distribution of amplitudes", detName.Data()))->Fill(dig->GetAmp());
339 }
340}
341
342// Trigger detectors
343
344template <class T> void BmnQaOffline::GetCommonInfo(TClonesArray* digiArray, BmnTrigDetQa* detHistoClass, TString detName) {
345 for (Int_t iDig = 0; iDig < digiArray->GetEntriesFast(); iDig++) {
346 T* dig = (T*) digiArray->UncheckedAt(iDig);
347 detHistoClass->GetManager()->H1(TString::Format("%sTRIGGERS_1d, %s, Distribution of inn. channels", prefix.Data(), detName.Data()))->Fill(dig->GetMod());
348 detHistoClass->GetManager()->H1(TString::Format("%sTRIGGERS_1d, %s, Distribution of times", prefix.Data(), detName.Data()))->Fill(dig->GetTime());
349 detHistoClass->GetManager()->H1(TString::Format("%sTRIGGERS_1d, %s, Distribution of amplitudes", prefix.Data(), detName.Data()))->Fill(dig->GetAmp());
350 }
351}
352
353// Dst
354
355void BmnQaOffline::GetBasicTrackDistributions() {
356 if (!fGlobalTracks)
357 return;
358
359 const Int_t nDims = 3;
360 TString dim[nDims] = {"X", "Y", "Z"};
361
362 dst->GetManager()->H1(Form("%sDST_1d, Distribution of total multiplicity", prefix.Data()))->Fill(fGlobalTracks->GetEntriesFast());
363
364 for (Int_t iTrack = 0; iTrack < fGlobalTracks->GetEntriesFast(); iTrack++) {
365 BmnGlobalTrack* track = (BmnGlobalTrack*) fGlobalTracks->UncheckedAt(iTrack);
366 dst->GetManager()->H1(Form("%sDST_1d, Distribution of momenta", prefix.Data()))->Fill(1. / track->GetParamFirst()->GetQp());
367 dst->GetManager()->H1(Form("%sDST_1d, Distribution of Nhits", prefix.Data()))->Fill(track->GetNHits());
368
369 FairTrackParam* parFirst = track->GetParamFirst();
370 FairTrackParam* parLast = track->GetParamLast();
371
372 Double_t xyzFirst[nDims] = {parFirst->GetX(), parFirst->GetY(), parFirst->GetZ()};
373 Double_t xyzLast[nDims] = {parLast->GetX(), parLast->GetY(), parLast->GetZ()};
374
375 Double_t txtyFirst[] = {parFirst->GetTx(), parFirst->GetTy()};
376 Double_t txtyLast[] = {parLast->GetTx(), parLast->GetTy()};
377
378 for (Int_t iDim = 0; iDim < nDims; iDim++) {
379 dst->GetManager()->H1(Form("%sDST_1d, Distribution of start%s", prefix.Data(), dim[iDim].Data()))->Fill(xyzFirst[iDim]);
380 if (!dim[iDim].Contains("Z"))
381 dst->GetManager()->H1(Form("%sDST_1d, Distribution of start T%s", prefix.Data(), dim[iDim].Data()))->Fill(txtyFirst[iDim]);
382 dst->GetManager()->H1(Form("%sDST_1d, Distribution of last%s", prefix.Data(), dim[iDim].Data()))->Fill(xyzLast[iDim]);
383 if (!dim[iDim].Contains("Z"))
384 dst->GetManager()->H1(Form("%sDST_1d, Distribution of last T%s", prefix.Data(), dim[iDim].Data()))->Fill(txtyLast[iDim]);
385 }
386 }
387
388 if (fSiliconTracks)
389 for (Int_t iTrack = 0; iTrack < fSiliconTracks->GetEntriesFast(); iTrack++) {
390 BmnSiliconTrack* track = (BmnSiliconTrack*) fSiliconTracks->UncheckedAt(iTrack);
391 dst->GetManager()->H1(Form("%sDST_1d, Distribution of Nhits, %s track", prefix.Data(), "SILICON"))->Fill(track->GetNHits());
392 }
393
394 if (fGemTracks)
395 for (Int_t iTrack = 0; iTrack < fGemTracks->GetEntriesFast(); iTrack++) {
396 BmnGemTrack* track = (BmnGemTrack*) fGemTracks->UncheckedAt(iTrack);
397 dst->GetManager()->H1(Form("%sDST_1d, Distribution of Nhits, %s track", prefix.Data(), "GEM"))->Fill(track->GetNHits());
398 }
399
400 if (fVertex)
401 dst->GetManager()->H2(Form("%sDST_2d, Vp_{z} vs. Ntracks", prefix.Data()))->Fill(fGlobalTracks->GetEntriesFast(), ((CbmVertex*) fVertex->UncheckedAt(0))->GetZ());
402}
403
404void BmnQaOffline::doHitsDistributions(TClonesArray* hitsArray, BmnDstQa* detHistoClass, TString detName) {
405 const Int_t nDims = 3;
406 TString dim[nDims] = {"X", "Y", "Z"};
407
408 if (!hitsArray)
409 return;
410
411 // Basic distributions: coordinates and errors (all detectors)
412 for (Int_t iHit = 0; iHit < hitsArray->GetEntriesFast(); iHit++) {
413 BmnHit* hit = (BmnHit*) hitsArray->UncheckedAt(iHit);
414
415 for (Int_t iDim = 0; iDim < nDims; iDim++) {
416 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of %s", prefix.Data(), detName.Data(), dim[iDim].Data()))->Fill(
417 (iDim == 0) ? hit->GetX() : (iDim == 1) ? hit->GetY() : hit->GetZ());
418 if (iDim == 2)
419 continue;
420 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of d%s", prefix.Data(), detName.Data(), dim[iDim].Data()))->Fill(
421 (iDim == 0) ? hit->GetDx() : hit->GetDy());
422 }
423 }
424
425 // Specific distributions for SILICON, GEM and CSC
426 if (detName.Contains("SILICON"))
427 for (Int_t iHit = 0; iHit < hitsArray->GetEntriesFast(); iHit++) {
428 BmnSiliconHit* hit = (BmnSiliconHit*) hitsArray->UncheckedAt(iHit);
429
430 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X cluster size", prefix.Data(), detName.Data()))->Fill(hit->GetClusterSizeInLowerLayer());
431 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X' cluster size", prefix.Data(), detName.Data()))->Fill(hit->GetClusterSizeInUpperLayer());
432
433 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X cluster signals", prefix.Data(), detName.Data()))->Fill(hit->GetStripTotalSignalInLowerLayer());
434 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X' cluster signals", prefix.Data(), detName.Data()))->Fill(hit->GetStripTotalSignalInUpperLayer());
435 } else if (detName.Contains("GEM"))
436 for (Int_t iHit = 0; iHit < hitsArray->GetEntriesFast(); iHit++) {
437 BmnGemStripHit* hit = (BmnGemStripHit*) hitsArray->UncheckedAt(iHit);
438
439 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X cluster size", prefix.Data(), detName.Data()))->Fill(hit->GetClusterSizeInLowerLayer());
440 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X' cluster size", prefix.Data(), detName.Data()))->Fill(hit->GetClusterSizeInUpperLayer());
441
442 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X cluster signals", prefix.Data(), detName.Data()))->Fill(hit->GetStripTotalSignalInLowerLayer());
443 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X' cluster signals", prefix.Data(), detName.Data()))->Fill(hit->GetStripTotalSignalInUpperLayer());
444 } else if (detName.Contains("CSC"))
445 for (Int_t iHit = 0; iHit < hitsArray->GetEntriesFast(); iHit++) {
446 BmnCSCHit* hit = (BmnCSCHit*) hitsArray->UncheckedAt(iHit);
447
448 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X cluster size", prefix.Data(), detName.Data()))->Fill(hit->GetClusterSizeInLowerLayer());
449 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X' cluster size", prefix.Data(), detName.Data()))->Fill(hit->GetClusterSizeInUpperLayer());
450
451 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X cluster signals", prefix.Data(), detName.Data()))->Fill(hit->GetStripTotalSignalInLowerLayer());
452 detHistoClass->GetManager()->H1(Form("%sDST_1d, %s, Distribution of X' cluster signals", prefix.Data(), detName.Data()))->Fill(hit->GetStripTotalSignalInUpperLayer());
453 }
454}
455
456BmnHit* BmnQaOffline::MatchDetector(FairTrackParam* parTrack, TClonesArray* hits, Bool_t doUpdate, Int_t detNum) {
457 FairTrackParam minPar;
458 Double_t minDX = DBL_MAX;
459 Double_t minDY = DBL_MAX;
460 BmnHit* minHit = nullptr;
461 Double_t xCut = 5.0;
462 Double_t yCut = 5.0;
463 Int_t pdg = 2212;
464
465 Double_t minZ = DBL_MAX;
466 BmnKalmanFilter* kalman = new BmnKalmanFilter();
467
468 map<Double_t, FairTrackParam> zParamMap;
469 for (Int_t iHit = 0; iHit < hits->GetEntriesFast(); ++iHit) {
470 BmnHit* hit = (BmnHit*) hits->At(iHit);
471 if (!hit) continue;
472 if (hit->GetStation() != detNum) continue;
473 minZ = Min(minZ, hit->GetZ());
474 zParamMap[hit->GetZ()] = FairTrackParam();
475 }
476 FairTrackParam par(*parTrack);
477 kalman->TGeoTrackPropagate(&par, minZ, pdg, NULL, NULL, isField);
478
479 for (map<Double_t, FairTrackParam>::iterator it = zParamMap.begin(); it != zParamMap.end(); ++it) {
480 (*it).second = par;
481 kalman->TGeoTrackPropagate(&(*it).second, (*it).first, pdg, NULL, NULL, isField);
482 }
483
484 for (Int_t iHit = 0; iHit < hits->GetEntriesFast(); ++iHit) {
485 BmnHit* hit = (BmnHit*) hits->At(iHit);
486 if (!hit) continue;
487 if (hit->GetStation() != detNum) continue;
488 if (zParamMap.find(hit->GetZ()) == zParamMap.end()) {
489 std::cout << "-E- CbmLitNearestHitTofMerger::DoMerge: Z position " << hit->GetZ() << " not found in map. Something is wrong.\n";
490 }
491 FairTrackParam tpar(zParamMap[hit->GetZ()]);
492 Double_t dX = tpar.GetX() - hit->GetX();
493 Double_t dY = tpar.GetY() - hit->GetY();
494 if (Abs(dX) < xCut && Abs(dY) < yCut && Abs(dX) < minDX && Abs(dY) < minDY) {
495 minHit = hit;
496 minDX = dX;
497 minDY = dY;
498 minPar = tpar;
499 }
500 }
501
502 *(parTrack) = par;
503
504 if (minHit) {
505 *(parTrack) = minPar;
506 if (doUpdate) {
507 Double_t chi = 0;
508 kalman->Update(parTrack, minHit, chi);
509 }
510 }
511 delete kalman;
512
513 return minHit;
514}
515
516void BmnQaOffline::doMatchingAnal() {
517 for (Int_t iTrack = 0; iTrack < fGemTracks->GetEntriesFast(); iTrack++) {
518 BmnGemTrack* track = (BmnGemTrack*) fGemTracks->UncheckedAt(iTrack);
519
520 const Int_t ndets = 4;
521 const Int_t nSeqs = 4;
522 TString dets[ndets] = {"CSC", "DCH1", "TOF700", "DCH2"};
523 TClonesArray * detSequence[ndets][nSeqs] = {
524 {fDchHits, fTof700Hits, fDchHits, fCscHits},
525 {fCscHits, fTof700Hits, fDchHits, fDchHits},
526 {fCscHits, fDchHits, fDchHits, fTof700Hits},
527 {fCscHits, fDchHits, fTof700Hits, fDchHits}
528 };
529
530 Int_t idx[ndets][nSeqs] = {
531 {0, -1, 1, 0},
532 {0, -1, 1, 0},
533 {0, 0, 1, -1},
534 {0, 0, -1, 1}
535 };
536 Bool_t doUpdate[ndets][nSeqs] = {
537 {kTRUE, kFALSE, kTRUE, kFALSE},
538 {kTRUE, kFALSE, kTRUE, kFALSE},
539 {kTRUE, kTRUE, kTRUE, kFALSE},
540 {kTRUE, kTRUE, kFALSE, kFALSE}
541 };
542
543 for (Int_t iDet = 0; iDet < ndets; iDet++) {
544 FairTrackParam parLast(*(track->GetParamLast()));
545 if (MatchDetector(&parLast, detSequence[iDet][0], doUpdate[iDet][0], idx[iDet][0]))
546 if (MatchDetector(&parLast, detSequence[iDet][1], doUpdate[iDet][1], idx[iDet][1]))
547 if (MatchDetector(&parLast, detSequence[iDet][2], doUpdate[iDet][2], idx[iDet][2])) {
548 BmnHit* hit = MatchDetector(&parLast, detSequence[iDet][3], doUpdate[iDet][3], idx[iDet][3]);
549 dst->GetManager()->H1(Form("%sDST_1d, hAllNo%s", prefix.Data(), dets[iDet].Data()))->Fill(1.0 / parLast.GetQp());
550 dst->GetManager()->H2(Form("%sDST_2d, hAllNo%sXY", prefix.Data(), dets[iDet].Data()))->Fill(parLast.GetX(), parLast.GetY());
551 if (hit) {
552 dst->GetManager()->H1(Form("%sDST_1d, hAllYes%s", prefix.Data(), dets[iDet].Data()))->Fill(1.0 / parLast.GetQp());
553 dst->GetManager()->H2(Form("%sDST_2d, hAllYes%sXY", prefix.Data(), dets[iDet].Data()))->Fill(parLast.GetX(), parLast.GetY());
554 }
555 }
556 }
557 }
558}
559
560void BmnQaOffline::doEfficiencyAnal() {
561 const Int_t ndets = 4;
562 TString dets[ndets] = {"CSC", "DCH1", "TOF700", "DCH2"};
563
564 for (Int_t iDet = 0; iDet < ndets; iDet++) {
565 TH1F* denom1 = (TH1F*) dst->GetManager()->H1(Form("%sDST_1d, hAllNo%s", prefix.Data(), dets[iDet].Data()));
566 denom1->Sumw2();
567
568 TH1F* nom1 = (TH1F*) dst->GetManager()->H1(Form("%sDST_1d, hAllYes%s", prefix.Data(), dets[iDet].Data()));
569 nom1->Sumw2();
570
571 TH1F* eff1 = (TH1F*) dst->GetManager()->H1(Form("%sDST_1d, hEff%s", prefix.Data(), dets[iDet].Data()));
572 eff1->Divide(nom1, denom1);
573 eff1->Scale(100);
574 eff1->SetMaximum(110);
575 eff1->SetMinimum(0);
576
577 TH2F* denom2 = (TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hAllNo%sXY", prefix.Data(), dets[iDet].Data()));
578 denom2->Sumw2();
579
580 TH2F* nom2 = (TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hAllYes%sXY", prefix.Data(), dets[iDet].Data()));
581 nom2->Sumw2();
582
583 TH2F* eff2 = (TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hEff%sXY", prefix.Data(), dets[iDet].Data()));
584 eff2->Divide(nom2, denom2);
585 }
586}
587
588void BmnQaOffline::doResidualsPullsAnal() {
589 const Int_t ndets = 8;
590 TString dets[ndets] = {"MWPC", "SILICON", "GEM", "CSC", "TOF400", "DCH1", "TOF700", "DCH2"};
591
592 TClonesArray * detSequence[ndets] = {fMwpcHits, fSiliconHits, fGemHits, fCscHits, fTof400Hits, fDchHits, fTof700Hits, fDchHits};
593
594 // Residuals ...
595 for (Int_t iDet = 0; iDet < ndets; iDet++) {
596 if (!detSequence[iDet])
597 continue;
598 for (Int_t iHit = 0; iHit < detSequence[iDet]->GetEntriesFast(); iHit++) {
599 BmnHit* hit = (BmnHit*) detSequence[iDet]->UncheckedAt(iHit);
600
601 if (dets[iDet].Contains("DCH")) {
602 dets[iDet] = (hit->GetStation() == 0) ? "DCH1" : "DCH2";
603 dst->GetManager()->H1(Form("%sDST_1d, hResX%s", prefix.Data(), dets[iDet].Data()))->Fill(hit->GetResX());
604 dst->GetManager()->H1(Form("%sDST_1d, hResY%s", prefix.Data(), dets[iDet].Data()))->Fill(hit->GetResY());
605 } else {
606 dst->GetManager()->H1(Form("%sDST_1d, hResX%s", prefix.Data(), dets[iDet].Data()))->Fill(hit->GetResX());
607 dst->GetManager()->H1(Form("%sDST_1d, hResY%s", prefix.Data(), dets[iDet].Data()))->Fill(hit->GetResY());
608 }
609 }
610 }
611
612 // Pulls ...
613 for (Int_t iGlob = 0; iGlob < fGlobalTracks->GetEntriesFast(); iGlob++) {
614 BmnGlobalTrack* glTrack = (BmnGlobalTrack*) fGlobalTracks->UncheckedAt(iGlob);
615
616 BmnKalmanFilter* kalman = new BmnKalmanFilter();
617 Int_t pdg = 2212;
618
619 FairTrackParam* first = glTrack->GetParamFirst(); // param. estimation
620
621 // DCH pulls ...
622 for (Int_t iHit = 0; iHit < fDchHits->GetEntriesFast(); iHit++) {
623 BmnHit* hit = (BmnHit*) fDchHits->UncheckedAt(iHit);
624
625 // BmnTrack* track = (BmnTrack*) fDchTracks->UncheckedAt(hit->GetIndex());
626
627 Double_t X = hit->GetX();
628 Double_t dX = hit->GetDx();
629
630 Double_t Y = hit->GetY();
631 Double_t dY = hit->GetDy();
632
633 FairTrackParam par(*first);
634 kalman->TGeoTrackPropagate(&par, hit->GetZ(), pdg, nullptr, nullptr, isField);
635
636 TString det = (hit->GetStation() == 0) ? "DCH1" : "DCH2";
637 dst->GetManager()->H1(Form("%sDST_1d, hPullX%s", prefix.Data(), det.Data()))->Fill(((par.GetX() - X) / dX));
638 dst->GetManager()->H1(Form("%sDST_1d, hPullY%s", prefix.Data(), det.Data()))->Fill(((par.GetY() - Y) / dY));
639 }
640
641 vector <TString> trackDets = {"MWPC", "SILICON", "GEM"};
642 vector <TClonesArray*> tracks = {fMwpcTracks, fSiliconTracks, fGemTracks};
643 vector <TClonesArray*> trackHits = {fMwpcHits, fSiliconHits, fGemHits};
644
645 for (auto detector : trackDets) {
646 if (detector.Contains("MWPC")) // FIXME
647 continue;
648
649 Int_t idx = detector.Contains("SILICON") ? glTrack->GetSilTrackIndex() : glTrack->GetGemTrackIndex();
650 Int_t idxArr = detector.Contains("MWPC") ? 0 : detector.Contains("SILICON") ? 1 : 2;
651
652 if (idx != -1) {
653 BmnTrack* track = (BmnTrack*) tracks[idxArr]->UncheckedAt(idx);
654
655 for (Int_t iHit = 0; iHit < track->GetNHits(); iHit++) {
656 BmnHit* hit = (BmnHit*) trackHits[idxArr]->UncheckedAt(track->GetHitIndex(iHit));
657
658 Double_t X = hit->GetX();
659 Double_t dX = hit->GetDx();
660
661 Double_t Y = hit->GetY();
662 Double_t dY = hit->GetDy();
663
664 FairTrackParam par(*first);
665 kalman->TGeoTrackPropagate(&par, hit->GetZ(), pdg, nullptr, nullptr, isField);
666
667 dst->GetManager()->H1(Form("%sDST_1d, hPullX%s", prefix.Data(), detector.Data()))->Fill(((par.GetX() - X) / dX));
668 dst->GetManager()->H1(Form("%sDST_1d, hPullY%s", prefix.Data(), detector.Data()))->Fill(((par.GetY() - Y) / dY));
669 }
670 }
671 }
672
673 vector <TString> hitDets = {"CSC", "TOF400", "TOF700"};
674 vector <TClonesArray*> hits = {fCscHits, fTof400Hits, fTof700Hits};
675
676 for (auto detector : hitDets) {
677 Int_t idx = detector.Contains("CSC") ? glTrack->GetCscHitIndex(0) :
678 detector.Contains("TOF400") ? glTrack->GetTof1HitIndex() :
679 glTrack->GetTof2HitIndex();
680 Int_t idxArr = detector.Contains("CSC") ? 0 : detector.Contains("TOF400") ? 1 : 2;
681
682 if (idx != -1) {
683 BmnHit* hit = (BmnHit*) hits[idxArr]->UncheckedAt(idx);
684
685 Double_t X = hit->GetX();
686 Double_t dX = hit->GetDx();
687
688 Double_t Y = hit->GetY();
689 Double_t dY = hit->GetDy();
690
691 FairTrackParam par(*first);
692 kalman->TGeoTrackPropagate(&par, hit->GetZ(), pdg, nullptr, nullptr, isField);
693
694 dst->GetManager()->H1(Form("%sDST_1d, hPullX%s", prefix.Data(), detector.Data()))->Fill(((par.GetX() - X) / dX));
695 dst->GetManager()->H1(Form("%sDST_1d, hPullY%s", prefix.Data(), detector.Data()))->Fill(((par.GetY() - Y) / dY));
696 }
697 }
698 delete kalman;
699 }
700}
701
702void BmnQaOffline::doOccupancy(TClonesArray* hits) {
703 if (!hits)
704 return;
705
706 const Double_t zBorder = 500.; // FIXME
707 for (Int_t iHit = 0; iHit < hits->GetEntriesFast(); iHit++) {
708 BmnHit* hit = (BmnHit*) hits->UncheckedAt(iHit);
709
710 Double_t x = hit->GetX();
711 Double_t y = hit->GetY();
712 Double_t z = hit->GetZ();
713
714 TString className = (TString) hit->GetName();
715 TString detName = (className.Contains("Gem")) ? "GEM" :
716 (className.Contains("Silicon")) ? "SILICON" :
717 (className.Contains("CSC")) ? "CSC" : "";
718
719 if (detName.IsNull()) {
720 if (className.Contains("Tof") && (z < zBorder))
721 detName = "TOF400";
722 else if (className.Contains("Tof") && (z > zBorder))
723 detName = "TOF700";
724 }
725
726 if (detName.Contains("GEM")) {
727 for (Int_t iStat = 0; iStat < fDetGem->GetNStations(); iStat++)
728 if (hit->GetStation() == iStat)
729 dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s, stat %d", prefix.Data(), detName.Data(), iStat))->Fill(x, y);
730 } else if (detName.Contains("SILICON")) {
731 for (Int_t iStat = 0; iStat < fDetSilicon->GetNStations(); iStat++)
732 if (hit->GetStation() == iStat)
733 dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s, stat %d", prefix.Data(), detName.Data(), iStat))->Fill(x, y);
734 } else if (detName.Contains("CSC")) {
735 for (Int_t iStat = 0; iStat < fDetCsc->GetNStations(); iStat++)
736 if (hit->GetStation() == iStat)
737 dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s, stat %d", prefix.Data(), detName.Data(), iStat))->Fill(x, y);
738 } else if (detName.Contains("TOF400") || detName.Contains("TOF700"))
739 dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s", prefix.Data(), detName.Data()))->Fill(x, y);
740 }
741}
742
743void BmnQaOffline::doOccupancyAnal() {
744 const Int_t ndets = 5;
745 TString dets[ndets] = {"SILICON", "GEM", "CSC", "TOF400", "TOF700"};
746
747 for (Int_t iDet = 0; iDet < ndets; iDet++) {
748 if (dets[iDet].Contains("GEM"))
749 for (Int_t iStat = 0; iStat < fDetGem->GetNStations(); iStat++)
750 changeHistoContent((TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s, stat %d", prefix.Data(), dets[iDet].Data(), iStat)));
751
752 else if (dets[iDet].Contains("SILICON"))
753 for (Int_t iStat = 0; iStat < fDetSilicon->GetNStations(); iStat++)
754 changeHistoContent((TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s, stat %d", prefix.Data(), dets[iDet].Data(), iStat)));
755
756 else if (dets[iDet].Contains("CSC"))
757 for (Int_t iStat = 0; iStat < fDetCsc->GetNStations(); iStat++)
758 changeHistoContent((TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s, stat %d", prefix.Data(), dets[iDet].Data(), iStat)));
759
760 else
761 changeHistoContent((TH2F*) dst->GetManager()->H2(Form("%sDST_2d, hOccupancy%s", prefix.Data(), dets[iDet].Data())));
762 }
763}
764
765void BmnQaOffline::changeHistoContent(TH2F* h) {
766 for (Int_t iBinX = 0; iBinX < h->GetNbinsX() + 1; iBinX++)
767 for (Int_t iBinY = 1; iBinY < h->GetNbinsY() + 1; iBinY++) {
768 Double_t wX = h->GetXaxis()->GetBinWidth(iBinX);
769 Double_t wY = h->GetYaxis()->GetBinWidth(iBinY);
770 Double_t content = h->GetBinContent(iBinX, iBinY);
771
772 content /= (wX * wY * fCurrentEvent);
773 h->SetBinContent(iBinX, iBinY, content);
774 }
775}
776
777void BmnQaOffline::doAverageStripValuePerHit(TClonesArray* hits) {
778 if (!hits)
779 return;
780
781 const Int_t nLayers = 2; // X and X'
782 TString layers[nLayers] = {"X", "X'"};
783
784 for (Int_t iHit = 0; iHit < hits->GetEntriesFast(); iHit++) {
785 BmnHit* hit = (BmnHit*) hits->UncheckedAt(iHit);
786
787 TString className = (TString) hit->GetName();
788 TString detName = (className.Contains("Gem")) ? "GEM" :
789 (className.Contains("Silicon")) ? "SILICON" :
790 (className.Contains("CSC")) ? "CSC" : "";
791
792 if (detName.IsNull())
793 continue;
794
795 if (detName.Contains("GEM")) {
796 BmnGemStripHit* gemHit = (BmnGemStripHit*) hit;
797 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per hit in %s (%s)", prefix.Data(), detName.Data(), layers[0].Data()))->Fill(gemHit->GetClusterSizeInLowerLayer());
798 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per hit in %s (%s)", prefix.Data(), detName.Data(), layers[1].Data()))->Fill(gemHit->GetClusterSizeInUpperLayer());
799 } else if (detName.Contains("SILICON")) {
800 BmnSiliconHit* silHit = (BmnSiliconHit*) hit;
801 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per hit in %s (%s)", prefix.Data(), detName.Data(), layers[0].Data()))->Fill(silHit->GetClusterSizeInLowerLayer());
802 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per hit in %s (%s)", prefix.Data(), detName.Data(), layers[1].Data()))->Fill(silHit->GetClusterSizeInUpperLayer());
803 } else if (detName.Contains("CSC")) {
804 BmnCSCHit* cscHit = (BmnCSCHit*) hit;
805 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per hit in %s (%s)", prefix.Data(), detName.Data(), layers[0].Data()))->Fill(cscHit->GetClusterSizeInLowerLayer());
806 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per hit in %s (%s)", prefix.Data(), detName.Data(), layers[1].Data()))->Fill(cscHit->GetClusterSizeInUpperLayer());
807 }
808 }
809}
810
811void BmnQaOffline::doAverageStripValuePerTrack() {
812 if (!fGlobalTracks)
813 return;
814
815 const Int_t nLayers = 2; // X and X'
816 TString layers[nLayers] = {"X", "X'"};
817
818 vector <TString> trackDets = {"SILICON", "GEM"};
819 vector <TClonesArray*> tracks = {fSiliconTracks, fGemTracks};
820 vector <TClonesArray*> trackHits = {fSiliconHits, fGemHits};
821
822 for (Int_t iTrack = 0; iTrack < fGlobalTracks->GetEntriesFast(); iTrack++) {
823 BmnGlobalTrack* glTrack = (BmnGlobalTrack*) fGlobalTracks->UncheckedAt(iTrack);
824
825 // Get CSC ...
826 Int_t idxCsc = glTrack->GetCscHitIndex(0);
827 if (idxCsc != -1) {
828 BmnCSCHit* hit = (BmnCSCHit*) fCscHits->UncheckedAt(idxCsc);
829 for (Int_t iLayer = 0; iLayer < nLayers; iLayer++) {
830 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per track in %s (%s)", prefix.Data(), "CSC", layers[0].Data()))->Fill(hit->GetClusterSizeInLowerLayer());
831 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per track in %s (%s)", prefix.Data(), "CSC", layers[1].Data()))->Fill(hit->GetClusterSizeInUpperLayer());
832 }
833 }
834
835 for (auto detector : trackDets) {
836 Int_t idx = (detector.Contains("SILICON")) ? glTrack->GetSilTrackIndex() : glTrack->GetGemTrackIndex();
837 Int_t idxArr = detector.Contains("SILICON") ? 0 : 1;
838
839 if (idx == -1)
840 continue;
841
842 BmnTrack* track = (BmnTrack*) tracks[idxArr]->UncheckedAt(idx);
843
844 Double_t sumLower = 0.;
845 Double_t sumUpper = 0.;
846
847 if (detector.Contains("SILICON")) {
848 for (Int_t iHit = 0; iHit < track->GetNHits(); iHit++) {
849 BmnSiliconHit* hit = (BmnSiliconHit*) trackHits[idxArr]->UncheckedAt(track->GetHitIndex(iHit));
850 sumLower += hit->GetClusterSizeInLowerLayer();
851 sumUpper += hit->GetClusterSizeInUpperLayer();
852 }
853 for (Int_t iLayer = 0; iLayer < nLayers; iLayer++) {
854 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per track in %s (%s)", prefix.Data(), detector.Data(), layers[0].Data()))->Fill(sumLower / track->GetNHits());
855 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per track in %s (%s)", prefix.Data(), detector.Data(), layers[1].Data()))->Fill(sumUpper / track->GetNHits());
856 }
857 } else {
858 for (Int_t iHit = 0; iHit < track->GetNHits(); iHit++) {
859 BmnGemStripHit* hit = (BmnGemStripHit*) trackHits[idxArr]->UncheckedAt(track->GetHitIndex(iHit));
860 sumLower += hit->GetClusterSizeInLowerLayer();
861 sumUpper += hit->GetClusterSizeInUpperLayer();
862 }
863 for (Int_t iLayer = 0; iLayer < nLayers; iLayer++) {
864 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per track in %s (%s)", prefix.Data(), detector.Data(), layers[0].Data()))->Fill(sumLower / track->GetNHits());
865 dst->GetManager()->H1(Form("%sDST_1d, Ave. strip value per track in %s (%s)", prefix.Data(), detector.Data(), layers[1].Data()))->Fill(sumUpper / track->GetNHits());
866 }
867 }
868 }
869 }
870}
Double_t GetStripTotalSignalInLowerLayer()
Definition BmnCSCHit.h:48
Double_t GetStripTotalSignalInUpperLayer()
Definition BmnCSCHit.h:52
Int_t GetClusterSizeInLowerLayer()
Definition BmnCSCHit.h:64
Int_t GetClusterSizeInUpperLayer()
Definition BmnCSCHit.h:68
BmnQaHistoManager * GetManager()
BmnQaHistoManager * GetManager()
BmnQaHistoManager * GetManager()
Definition BmnDstQa.h:29
Int_t GetClusterSizeInUpperLayer()
Double_t GetStripTotalSignalInUpperLayer()
Int_t GetClusterSizeInLowerLayer()
Double_t GetStripTotalSignalInLowerLayer()
Int_t GetTof2HitIndex() const
Int_t GetCscHitIndex(Int_t idx) const
Int_t GetGemTrackIndex() const
Int_t GetSilTrackIndex() const
Int_t GetTof1HitIndex() const
TH2 * H2(const TString &name) const
Return pointer to TH2 histogram.
void WriteToFile()
Write all histograms to current opened file.
TH1 * H1(const TString &name) const
Return pointer to TH1 histogram.
Double_t GetResX()
Definition BmnHit.h:98
Double_t GetResY()
Definition BmnHit.h:102
Short_t GetStation() const
Definition BmnHit.h:45
BmnStatus TGeoTrackPropagate(FairTrackParam *par, Double_t zOut, Int_t pdg, vector< Double_t > *F, Double_t *length, Bool_t isField)
BmnStatus Update(FairTrackParam *par, const BmnHit *hit, Double_t &chiSq)
void SetGeometriesByRunId(Int_t id, BmnGemStripStationSet *&gem, BmnSiliconStationSet *&silicon, BmnCSCStationSet *&csc)
Int_t GetNumberOfDets(Int_t period, TString setup, TString type)
vector< TString > GetDetectors(Int_t period, TString setup="BM@N")
pair< Int_t, TString > GetRunAndSetupByRunId(Int_t id)
vector< TString > GetTriggers(Int_t period, TString setup="BM@N")
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
virtual void Finish()
Double_t GetStripTotalSignalInLowerLayer()
Double_t GetStripTotalSignalInUpperLayer()
Int_t GetClusterSizeInLowerLayer()
Int_t GetClusterSizeInUpperLayer()
BmnQaHistoManager * GetManager()
FairTrackParam * GetParamFirst()
Definition BmnTrack.h:72
Int_t GetNHits() const
Definition BmnTrack.h:44
Int_t GetHitIndex(Int_t iHit) const
Definition BmnTrack.h:48
FairTrackParam * GetParamLast()
Definition BmnTrack.h:76
BmnQaHistoManager * GetManager()
static UniRun * GetRun(int period_number, int run_number)
get run from the database
Definition UniRun.cxx:177
static int ReadGeometryFile(int period_number, int run_number, char *geo_file_path, bool usePrevGeometryIfMissing=false)
Definition UniRun.cxx:1105
double * GetFieldVoltage()
get field voltage of the current run
Definition UniRun.h:125
const UInt_t dim
Definition setup.py:1