Belle II Software  release-05-01-25
ARICHDQMModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Kindo Haruki, Luka Santelj *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // Own include
12 #include <arich/modules/arichDQM/ARICHDQMModule.h>
13 
14 // ARICH
15 #include <arich/dbobjects/ARICHGeometryConfig.h>
16 #include <arich/dbobjects/ARICHChannelMapping.h>
17 #include <arich/dbobjects/ARICHMergerMapping.h>
18 #include <arich/dbobjects/ARICHGeoDetectorPlane.h>
19 #include <arich/dbobjects/ARICHGeoAerogelPlane.h>
20 #include <arich/dataobjects/ARICHHit.h>
21 #include <arich/dataobjects/ARICHDigit.h>
22 #include <arich/dataobjects/ARICHAeroHit.h>
23 #include <arich/dataobjects/ARICHTrack.h>
24 #include <arich/dataobjects/ARICHLikelihood.h>
25 #include <arich/dataobjects/ARICHPhoton.h>
26 
27 #include <mdst/dataobjects/Track.h>
28 #include <mdst/dataobjects/MCParticle.h>
29 
30 // framework - DataStore
31 #include <framework/datastore/StoreArray.h>
32 
33 // Dataobject classes
34 #include <framework/database/DBObjPtr.h>
35 
36 #include <TF1.h>
37 #include <TVector3.h>
38 
39 #include <fstream>
40 #include <math.h>
41 
42 using namespace std;
43 
44 namespace Belle2 {
50  //-----------------------------------------------------------------
51  // Register module
52  //-----------------------------------------------------------------
53 
54  REG_MODULE(ARICHDQM);
55 
56  ARICHDQMModule::ARICHDQMModule() : HistoModule()
57  {
58  // set module description (e.g. insert text)
59  setDescription("Make summary of data quality.");
61  addParam("debug", m_debug, "debug mode", false);
62  addParam("UpperMomentumLimit", m_momUpLim, "Upper momentum limit of tracks included in monitoring", 10.0);
63  addParam("LowerMomentumLimit", m_momDnLim, "Lower momentum limit of tracks included in monitoring", 2.5);
64  addParam("ArichEvents", m_arichEvents, "Include only hits from events where an extrapolated track to arich exists", false);
65  addParam("MaxHits", m_maxHits, "Include only events with less than MaxHits hits in ARICH (remove loud events)", 70000);
66  addParam("MinHits", m_minHits, "Include only events with more than MinHits hits in ARICH", 0);
67  }
68 
70  {
71  }
72 
74  {
75 
76  TDirectory* oldDir = gDirectory;
77  TDirectory* dirARICHDQM = oldDir->mkdir("ARICH");
78  dirARICHDQM->cd();
79 
80  //Histograms for analysis and statistics
81 
82  h_chStat = new TH1D("chStat", "Status of channels;Channel serial;Status", 420 * 144, -0.5, 420 * 144 - 0.5);
83  h_aeroStat = new TH1D("aeroStat", "Status of aerogels;Aerogel tile serial;Status", 160, -0.5, 160 - 0.5);
84  h_chHit = new TH1D("chHit", "Number of hits in each channel;Channel serial;Hits", 420 * 144, -0.5, 420 * 144 - 0.5);
85  h_chipHit = new TH1D("chipHit", "Number of hits in each chip;Chip serial;Hits", 420 * 4, -0.5, 420 * 4 - 0.5);
86  h_hapdHit = new TH1D("hapdHit", "Number of hits in each HAPD;HAPD serial;Hits", 420, 0.5, 421 - 0.5);
87  h_hapdHitPerEvent = new TH2D("hapdHitPerEvent", "Number of hits in each HAPD per Event;HAPD serial;Hits/event", 420, 0.5, 420 + 0.5,
88  144, -0.5, 143.5);
89  h_trackPerEvent = new TH1D("trackPerEvent", "Number of tracks in ARICH per event; # of tracks;Events", 6, -0.5, 5.5);
90 
91  h_mergerHit = new TH1D("mergerHit", "Number of hits in each merger board;MB ID;Hits", 72, 0.5, 72 + 0.5);
92  h_bitsPerMergerNorm = new TH2D("bitsPerMergerNorm", "Normalised number of hits in each bit in each Merger;Bit;MB ID;Hits", 5,
93  -1 - 0.5,
94  4 - 0.5, 72, 0.5, 72 + 0.5); // copy of h_bits, normalised to number of connected Hapds and to sum(bit1, bit2)
95  h_bitsPerHapdMerger = new TH2D("bitsPerHapdMerger",
96  "Number of hits in each bit in each Hapd sorted by mergers;Bit;HAPD unsorted;Hits", 5, -1 - 0.5,
97  4 - 0.5, 432, 1, 432);
98 
99  h_aerogelHit = new TH1D("aerogelHit", "Number track associated hits in each aerogel tile;Aerogel slot ID;Hits", 125, -0.5,
100  125 - 0.5);
101  h_bits = new TH1D("bits", "Number of hits in each bit;Bit;Hits", 5, -1 - 0.5,
102  4 - 0.5); //Bin at -1 is added to set minimum 0 without SetMinimum(0) due to bugs in jsroot on DQM server.
103  h_bitsPerChannel = new TH2D("bitsPerChannel", "Number of hits in each bit in each chanenel;Bit;channel #;Hits", 4, - 0.5,
104  4 - 0.5, 420 * 144, -0.5, 420 * 144 - 0.5); // copy of h_bits
105  h_hitsPerTrack2D = new TH2D("hitsPerTrack2D", "2D distribution of track associated hits;X[cm];Y[cm];Hits", 230, -115, 115, 230,
106  -115, 115);
107  h_tracks2D = new TH2D("tracks2D", "Distribution track positions;X[cm];Y[cm];Tracks", 230, -115, 115, 230, -115, 115);
108 
109  h_hitsPerEvent = new TH1D("hitsPerEvent", "Number of hit per event;Number of hits;Events", 150, -0.5, 150 - 0.5);
110  h_theta = new TH1D("theta", "Cherenkov angle distribution;Angle [rad];Events", 60, 0, M_PI / 6);
111  h_hitsPerTrack = new TH1D("hitsPerTrack", "Number of hit per track;Number of hits;Tracks", 41, -0.5, 40.5);
112 
113  for (int i = 0; i < 6; i++) {
114  h_secTheta[i] = new TH1D(Form("thetaSec%d", i + 1), Form("Cherenkov angle distribution in sector %d;Angle [rad];Events", i + 1),
115  60, 0, M_PI / 6);
116  h_secHitsPerTrack[i] = new TH1D(Form("hitsPerTrackSec%d", i + 1),
117  Form("Number of hit per track in sector %d;Number of hits;Tracks", i + 1), 40, 0, 40);
118  h_secHapdHit[i] = new TH1D(Form("hapdHit%d", i + 1), Form("Number of hits in each HAPD in sector %d;HAPD serial;Hits", i + 1), 70,
119  0.5, 71 - 0.5);
120  }
121 
122  h_chDigit = new TH1D("chDigit", "Number of raw digits in each channel;Channel serial;Hits", 420 * 144, -0.5, 420 * 144 - 0.5);
123  h_chipDigit = new TH1D("chipDigit", "Number of raw digits in each chip;Chip serial;Hits", 420 * 4, -0.5, 420 * 4 - 0.5);
124  h_hapdDigit = new TH1D("hapdDigit", "Number of raw digits in each HAPD;HAPD serial;Hits", 420, 0.5, 421 - 0.5);
125 
126  h_aerogelHits3D = new TH3D("aerogelHits3D", "Number of track associated hits for each aerogel tile; #phi section; r section", 125,
127  -0.5, 124.5, 20, 0, 20, 20, 0, 20);
128  h_mirrorThetaPhi = new TH3D("mirrorThetaPhi",
129  "Cherenkov theta vs Cherenkov phi for mirror reflected photons; mirroID; #phi_{c} [rad]; #theta_{c} [rad]", 18, 0.5, 18.5, 100,
130  -M_PI, M_PI, 100, 0, 0.5);
131  h_thetaPhi = new TH2D("thetaPhi", "Cherenkov theta vs phi;#phi [rad];#theta_{c} [rad]", 100, -M_PI, M_PI, 100, 0., 0.5);
132 
133  h_flashPerAPD = new TH1D("flashPerAPD", "Number of flashes per APD; APD serial; number of flash", 420 * 4, -0.5, 420 * 4 - 0.5);
134 
135  h_ARICHOccAfterInjLer = new TH1F("ARICHOccInjLER", " ARICHOccInjLER /Time;Time in #mus;Nhits/Time (#mus bins)", 4000, 0, 20000);
136  h_ARICHOccAfterInjHer = new TH1F("ARICHOccInjHER", " ARICHOccInjHER /Time;Time in #mus;Nhits/Time (#mus bins)", 4000, 0, 20000);
137  h_ARICHEOccAfterInjLer = new TH1F("ARICHEOccInjLER", "ARICHEOccInjLER/Time;Time in #mus;Triggers/Time (#mus bins)", 4000, 0, 20000);
138  h_ARICHEOccAfterInjHer = new TH1F("ARICHEOccInjHER", "ARICHEOccInjHER/Time;Time in #mus;Triggers/Time (#mus bins)", 4000, 0, 20000);
139 
140  //Select "LIVE" monitoring histograms
141  h_chStat->SetOption("LIVE");
142  h_aeroStat->SetOption("LIVE");
143 
144  h_chHit->SetOption("LIVE");
145  h_chipHit->SetOption("LIVE");
146  h_hapdHit->SetOption("LIVE");
147 
148  h_chDigit->SetOption("LIVE");
149  h_chipDigit->SetOption("LIVE");
150  h_hapdDigit->SetOption("LIVE");
151  h_mergerHit->SetOption("LIVE");
152  h_bitsPerMergerNorm->SetOption("LIVE");
153  h_bitsPerHapdMerger->SetOption("LIVE");
154 
155  h_aerogelHit->SetOption("LIVE");
156  h_bits->SetOption("LIVE");
157  h_hitsPerTrack2D->SetOption("LIVE");
158  h_tracks2D->SetOption("LIVE");
159 
160  h_hitsPerEvent->SetOption("LIVE");
161  h_theta->SetOption("LIVE");
162  h_hitsPerTrack->SetOption("LIVE");
163  h_trackPerEvent->SetOption("LIVE");
164  h_flashPerAPD->SetOption("LivE");
165 
166  for (int i = 0; i < 6; i++) {
167  h_secTheta[i]->SetOption("LIVE");
168  h_secHitsPerTrack[i]->SetOption("LIVE");
169  }
170  h_ARICHOccAfterInjLer->SetOption("LIVE");
171  h_ARICHEOccAfterInjLer->SetOption("LIVE");
172  h_ARICHOccAfterInjHer->SetOption("LIVE");
173  h_ARICHEOccAfterInjHer->SetOption("LIVE");
174 
175  //Set the minimum to 0
176  h_chDigit->SetMinimum(0);
177  h_chipDigit->SetMinimum(0);
178  h_hapdDigit->SetMinimum(0);
179  h_chHit->SetMinimum(0);
180  h_chipHit->SetMinimum(0);
181  h_hapdHit->SetMinimum(0);
182  h_mergerHit->SetMinimum(0);
183  h_bitsPerMergerNorm->SetMinimum(0);
184  h_bitsPerHapdMerger->SetMinimum(0);
185  h_aerogelHit->SetMinimum(0);
186  h_bits->SetMinimum(0);
187  h_bitsPerChannel->SetMinimum(0);
188  h_hitsPerTrack2D->SetMinimum(0);
189  h_tracks2D->SetMinimum(0);
190  h_flashPerAPD->SetMinimum(0);
191  h_hitsPerEvent->SetMinimum(0);
192  h_theta->SetMinimum(0);
193  h_hitsPerTrack->SetMinimum(0);
194  h_ARICHOccAfterInjLer->SetMinimum(0);
195  h_ARICHEOccAfterInjLer->SetMinimum(0);
196  h_ARICHOccAfterInjHer->SetMinimum(0);
197  h_ARICHEOccAfterInjHer->SetMinimum(0);
198 
199  for (int i = 0; i < 6; i++) {
200  h_secTheta[i]->SetMinimum(0);
201  h_secHitsPerTrack[i]->SetMinimum(0);
202  }
203 
204  oldDir->cd();
205  }
206 
208  {
209  REG_HISTOGRAM
210  StoreArray<MCParticle> MCParticles;
211  MCParticles.isOptional();
212  StoreArray<Track> tracks;
213  tracks.isOptional();
214  StoreArray<ARICHHit> arichHits;
215  arichHits.isOptional();
216  StoreArray<ARICHDigit> arichDigits;
217  arichDigits.isRequired();
218  StoreArray<ARICHTrack> arichTracks;
219  arichTracks.isOptional();
220  StoreArray<ARICHAeroHit> arichAeroHits;
221  arichAeroHits.isOptional();
222  StoreArray<ARICHLikelihood> likelihoods;
223  likelihoods.isOptional();
224  m_rawFTSW.isOptional();
225  }
226 
228  {
229 
230  h_chStat->Reset();
231  h_aeroStat->Reset();
232  h_chDigit->Reset();
233  h_chipDigit->Reset();
234  h_hapdDigit->Reset();
235 
236  h_chHit->Reset();
237  h_chipHit->Reset();
238  h_hapdHit->Reset();
239  h_mergerHit->Reset();
240  h_bitsPerMergerNorm->Reset();
241  h_bitsPerHapdMerger->Reset();
242 
243  h_aerogelHit->Reset();
244  h_bits->Reset();
245  h_bitsPerChannel->Reset();
246  h_hitsPerTrack2D->Reset();
247  h_tracks2D->Reset();
248  h_aerogelHits3D->Reset();
249 
250  h_hitsPerEvent->Reset();
251  h_theta->Reset();
252  h_hitsPerTrack->Reset();
253  h_trackPerEvent->Reset();
254  h_flashPerAPD->Reset();
255  h_mirrorThetaPhi->Reset();
256  h_thetaPhi->Reset();
257 
258  for (int i = 0; i < 6; i++) {
259  h_secTheta[i]->Reset();
260  h_secHitsPerTrack[i]->Reset();
261  }
262 
263  h_ARICHOccAfterInjLer->Reset();
264  h_ARICHEOccAfterInjLer->Reset();
265  h_ARICHOccAfterInjHer->Reset();
266  h_ARICHEOccAfterInjHer->Reset();
267  }
268 
270  {
271  StoreArray<MCParticle> MCParticles;
272  StoreArray<Track> tracks;
273  StoreArray<ARICHDigit> arichDigits;
274  StoreArray<ARICHHit> arichHits;
275  StoreArray<ARICHTrack> arichTracks;
276  StoreArray<ARICHAeroHit> arichAeroHits;
277  StoreArray<ARICHLikelihood> arichLikelihoods;
278  DBObjPtr<ARICHGeometryConfig> arichGeoConfig;
279  const ARICHGeoDetectorPlane& arichGeoDec = arichGeoConfig->getDetectorPlane();
280  const ARICHGeoAerogelPlane& arichGeoAero = arichGeoConfig->getAerogelPlane();
281  DBObjPtr<ARICHChannelMapping> arichChannelMap;
282  DBObjPtr<ARICHMergerMapping> arichMergerMap;
283 
284  setReturnValue(1);
285 
286  if (arichHits.getEntries() < m_minHits || arichHits.getEntries() > m_maxHits) { setReturnValue(0); return;}
287 
288  if (!arichLikelihoods.getEntries() && m_arichEvents) { setReturnValue(0); return;}
289 
290  std::vector<int> apds(420 * 4, 0);
291  for (const auto& digit : arichDigits) {
292  uint8_t bits = digit.getBitmap();
293  int moduleID = digit.getModuleID();
294  int channelID = digit.getChannelID();
295  int mergerID = arichMergerMap->getMergerID(moduleID);
296  int febSlot = arichMergerMap->getFEBSlot(moduleID);
297  unsigned binID = (mergerID - 1) * N_FEB2MERGER + (febSlot);
298 
299  for (int i = 0; i < 8; i++) {
300  if ((bits & (1 << i)) && !(bits & ~(1 << i))) {
301  h_bits->Fill(i);
302  h_bitsPerChannel->Fill(i, (moduleID - 1) * 144 + channelID);
303  h_bitsPerMergerNorm->Fill(i, mergerID);
304  h_bitsPerHapdMerger->Fill(i, binID);
305  } else if (!bits) {
306  h_bits->Fill(8);
307  h_bitsPerChannel->Fill(8, (moduleID - 1) * 144 + channelID);
308  h_bitsPerMergerNorm->Fill(8, mergerID);
309  h_bitsPerHapdMerger->Fill(8, binID);
310  }
311  }
312 
313  // fill occupancy histograms for raw data
314  h_chDigit ->Fill((moduleID - 1) * 144 + channelID);
315  int chip = (moduleID - 1) * 4 + channelID / 36;
316  h_chipDigit->Fill(chip);
317  apds[chip] += 1;
318  h_hapdDigit->Fill(moduleID);
319  }
320 
321  int iapd = 0;
322  for (auto apd : apds) { if (apd > 20) h_flashPerAPD->Fill(iapd); iapd++;}
323 
324  std::vector<int> hpd(420, 0);
325  int nHit = 0;
326 
327  for (int i = 0; i < arichHits.getEntries(); i++) {
328 
329  int moduleID = arichHits[i]->getModule();
330  int channelID = arichHits[i]->getChannel();
331  h_chHit->Fill((moduleID - 1) * 144 + channelID);
332  hpd[moduleID - 1]++;
333  h_chipHit->Fill((moduleID - 1) * 4 + channelID / 36);
334  h_hapdHit->Fill(moduleID);
335  if (moduleID > 420) B2INFO("Invalid hapd number " << LogVar("hapd ID", moduleID));
336 
337  for (int j = 1; j <= 7; j++) {
338  int ringStart = (j - 1) * (84 + (j - 2) * 6) / 2 + 1; // The smallest module ID in each ring
339  int ringEnd = j * (84 + (j - 1) * 6) / 2; // The biggest module ID in each ring
340  if (ringStart <= moduleID && moduleID <= ringEnd) {
341  h_secHapdHit[(moduleID - ringStart) / (6 + j)]->Fill((moduleID - ringStart) % (6 + j) + 1 + (ringStart - 1) / 6);
342  }
343  }
344 
345  int mergerID = arichMergerMap->getMergerID(moduleID);
346  h_mergerHit->Fill(mergerID);
347  nHit++;
348  }
349 
350  h_hitsPerEvent->Fill(nHit);
351  int mmid = 1;
352  for (auto hh : hpd) { h_hapdHitPerEvent->Fill(mmid, hh); mmid++;}
353 
354  int ntrk = 0;
355  for (const auto& arichTrack : arichTracks) {
356 
357 
358  //Momentum limits are applied
359  //if (arichTrack.getPhotons().size() == 0) continue;
360  if (arichTrack.getMomentum() > 0.5) ntrk++; // count tracks with momentum larger than 0.5 GeV
361  if (arichTrack.getMomentum() < m_momDnLim || arichTrack.getMomentum() > m_momUpLim) continue;
362 
363  TVector3 recPos = arichTrack.getPosition();
364  int trSector = 0;
365  double dPhi = 0;
366  if (recPos.Phi() >= 0) {
367  dPhi = recPos.Phi();
368  } else {
369  dPhi = 2 * M_PI + recPos.Phi();
370  }
371  while (dPhi > M_PI / 3 && trSector < 5) {
372  dPhi -= M_PI / 3;
373  trSector++;
374  }
375  double trR = recPos.XYvector().Mod();
376  double trPhi = 0;
377  if (recPos.Phi() >= 0) {
378  trPhi = recPos.Phi();
379  } else {
380  trPhi = 2 * M_PI + recPos.Phi();
381  }
382 
383  int iRing = 0;
384  int iAzimuth = 0;
385  while (arichGeoAero.getRingRadius(iRing + 1) < trR && iRing < 4) {
386  iRing++;
387  }
388  if (iRing == 0) continue;
389  while (arichGeoAero.getRingDPhi(iRing) * (iAzimuth + 1) < trPhi) {
390  iAzimuth++;
391  }
392 
393  h_tracks2D->Fill(recPos.X(), recPos.Y());
394 
395  std::vector<ARICHPhoton> photons = arichTrack.getPhotons();
396  for (auto& photon : photons) {
397  if (photon.getMirror() == 0) {
398  if (trR < 93.) {
399  h_thetaPhi->Fill(photon.getPhiCer(), photon.getThetaCer());
400  h_theta->Fill(photon.getThetaCer());
401  }
402  int hitSector = 0;
403  double hitPhi = arichGeoDec.getSlotPhi(arichHits[photon.getHitID()]->getModule());
404  if (hitPhi < 0) hitPhi += 2 * M_PI;
405  while (hitPhi > M_PI / 3 && hitSector < 5) {
406  hitPhi -= M_PI / 3;
407  hitSector++;
408  }
409  h_secTheta[hitSector]->Fill(photon.getThetaCer());
410  } else {
411  if (trR > 95.) h_mirrorThetaPhi->Fill(photon.getMirror(), photon.getPhiCer(), photon.getThetaCer());
412  }
413  }
414 
415  //Get ARICHLikelihood related to the ARICHTrack
416  const ExtHit* extHit = arichTrack.getRelated<ExtHit>();
417  const Track* mdstTrack = NULL;
418  if (extHit) mdstTrack = extHit->getRelated<Track>();
419  const ARICHLikelihood* lkh = NULL;
420  if (mdstTrack) lkh = mdstTrack->getRelated<ARICHLikelihood>();
421  else lkh = arichTrack.getRelated<ARICHLikelihood>();
422 
423  if (lkh) {
424  if (!lkh->getFlag()) continue; //Fill only when the number of expected photons is more than 0.
425  double nphoton = lkh->getDetPhot();
426  h_hitsPerTrack->Fill(nphoton);
427  h_secHitsPerTrack[trSector]->Fill(nphoton);
428  h_hitsPerTrack2D->Fill(recPos.X(), recPos.Y(), nphoton);
429  int aeroID = arichGeoAero.getAerogelTileID(recPos.X(), recPos.Y());
430  h_aerogelHits3D->Fill(aeroID, (trPhi - arichGeoAero.getRingDPhi(iRing)*iAzimuth) / (arichGeoAero.getRingDPhi(iRing) / 20) ,
431  (trR - arichGeoAero.getRingRadius(iRing)) / ((arichGeoAero.getRingRadius(iRing + 1) - arichGeoAero.getRingRadius(iRing)) / 20),
432  nphoton);
433  h_aerogelHit->Fill(aeroID, nphoton);
434  }
435  }
436 
437  h_trackPerEvent->Fill(ntrk);
438 
439  for (auto& it : m_rawFTSW) {
440  B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
441  (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
442  it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
443  auto difference = it.GetTimeSinceLastInjection(0);
444  if (difference != 0x7FFFFFFF) {
445  unsigned int nentries = arichDigits.getEntries();
446  float diff2 = difference / 127.; // 127MHz clock ticks to us, inexact rounding
447  if (it.GetIsHER(0)) {
448  h_ARICHOccAfterInjHer->Fill(diff2, nentries);
449  h_ARICHEOccAfterInjHer->Fill(diff2);
450  } else {
451  h_ARICHOccAfterInjLer->Fill(diff2, nentries);
452  h_ARICHEOccAfterInjLer->Fill(diff2);
453  }
454  }
455  }
456 
457  }
458 
460  {
461 
462  DBObjPtr<ARICHMergerMapping> arichMergerMap;
463  if (h_theta->GetEntries() < 200) return;
464  TF1* f1 = new TF1("arichFitFunc", "gaus(0)+pol1(3)", 0.25, 0.4);
465  f1->SetParameters(0.8 * h_theta->GetMaximum(), 0.323, 0.016, 0, 0);
466  f1->SetParName(0, "C");
467  f1->SetParName(1, "mean");
468  f1->SetParName(2, "sigma");
469  f1->SetParName(3, "p0");
470  f1->SetParName(4, "p1");
471  h_theta->Fit(f1, "R");
472 
473  //Normalise bins in histogram bitsPerMergerNorm
474  for (int mergerID = 1; mergerID < 73; ++mergerID) {
475  double NHapd = 0;
476  for (int febSlot = 1; febSlot < 7; ++febSlot) {
477  if (arichMergerMap->getModuleID(mergerID, febSlot) > 0) NHapd++;
478  }
479 
480  double bin_value[5];
481  for (int i = 1; i <= 5; ++i) {
482  // loop over bits and save values
483  bin_value[i - 1] = h_bitsPerMergerNorm->GetBinContent(i, mergerID);
484  }
485 
486  for (int i = 1; i <= 5; ++i) {
487  // loop over bits again and set bin content
488  h_bitsPerMergerNorm->SetBinContent(i, mergerID,
489  bin_value[i - 1] / (bin_value[3] + bin_value[2]) / NHapd); // normalise with sum of bit1 and bit2, and number of connected HAPDs
490  }
491  }
492 
493 
494  }
495 
497  {
498  }
500 }
Belle2::ARICHDQMModule::h_flashPerAPD
TH1 * h_flashPerAPD
Number of flashes in each APD.
Definition: ARICHDQMModule.h:125
Belle2::ARICHLikelihood
This is a class to store ARICH likelihoods in the datastore.
Definition: ARICHLikelihood.h:38
Belle2::ARICHDQMModule::h_hitsPerEvent
TH1 * h_hitsPerEvent
Ihe number of all hits in each event.
Definition: ARICHDQMModule.h:121
Belle2::ARICHDQMModule::terminate
virtual void terminate() override
Termination action.
Definition: ARICHDQMModule.cc:496
Belle2::ARICHDQMModule::h_hitsPerTrack2D
TH2 * h_hitsPerTrack2D
Sum of 2D hit/track map on each position of track.
Definition: ARICHDQMModule.h:116
Belle2::ARICHDQMModule::h_ARICHOccAfterInjLer
TH1 * h_ARICHOccAfterInjLer
Histogram Ndigits after LER injection.
Definition: ARICHDQMModule.h:130
Belle2::ARICHDQMModule::h_hapdHit
TH1 * h_hapdHit
The number of hits in each HAPD.
Definition: ARICHDQMModule.h:106
Belle2::ARICHDQMModule::~ARICHDQMModule
virtual ~ARICHDQMModule()
Destructor.
Definition: ARICHDQMModule.cc:69
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::ARICHDQMModule::h_hapdHitPerEvent
TH2 * h_hapdHitPerEvent
number of hits in each HAPD per event
Definition: ARICHDQMModule.h:112
Belle2::ARICHDQMModule::h_thetaPhi
TH2 * h_thetaPhi
cherenkov theta vs phi for non-mirror-reflected photons
Definition: ARICHDQMModule.h:120
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::Module::c_ParallelProcessingCertified
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:82
Belle2::ARICHDQMModule::h_theta
TH1 * h_theta
Reconstructed Cherenkov angles.
Definition: ARICHDQMModule.h:122
Belle2::ARICHLikelihood::getFlag
int getFlag() const
Get reconstruction flag.
Definition: ARICHLikelihood.h:76
Belle2::ARICHDQMModule::h_mergerHit
TH1 * h_mergerHit
The number of hits in each Merger Boards.
Definition: ARICHDQMModule.h:107
Belle2::ARICHDQMModule::m_momUpLim
double m_momUpLim
Upper momentum limit of tracks used in GeV (if set 0, no limit is applied)
Definition: ARICHDQMModule.h:139
Belle2::ARICHGeoAerogelPlane::getRingRadius
double getRingRadius(unsigned iRing) const
Get radius of i-th aluminum ring between aerogel tiles (inner radius of ring)
Definition: ARICHGeoAerogelPlane.h:273
Belle2::ARICHDQMModule::h_ARICHOccAfterInjHer
TH1 * h_ARICHOccAfterInjHer
Histogram Ndigits after HER injection.
Definition: ARICHDQMModule.h:131
Belle2::ARICHDQMModule::h_chStat
TH1 * h_chStat
Status of each channels.
Definition: ARICHDQMModule.h:97
Belle2::RelationsInterface::getRelated
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Definition: RelationsObject.h:280
Belle2::ARICHDQMModule::h_aeroStat
TH1 * h_aeroStat
Status of each aerogel tiles.
Definition: ARICHDQMModule.h:98
Belle2::ARICHDQMModule::m_rawFTSW
StoreArray< RawFTSW > m_rawFTSW
Input array for DAQ Status.
Definition: ARICHDQMModule.h:135
Belle2::ARICHDQMModule::h_secHitsPerTrack
TH1 * h_secHitsPerTrack[6]
Detailed average hits/track for each sector.
Definition: ARICHDQMModule.h:128
Belle2::ARICHDQMModule::event
virtual void event() override
Event processor.
Definition: ARICHDQMModule.cc:269
Belle2::ARICHDQMModule::h_chipHit
TH1 * h_chipHit
The number of hits in each ASIC chip.
Definition: ARICHDQMModule.h:105
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::ARICHDQMModule::h_ARICHEOccAfterInjHer
TH1 * h_ARICHEOccAfterInjHer
Histogram for Nr Entries (=Triggrs) for normalization after HER injection.
Definition: ARICHDQMModule.h:133
Belle2::ARICHDQMModule::h_aerogelHit
TH1 * h_aerogelHit
The number of reconstructed photons in each aerogel tiles.
Definition: ARICHDQMModule.h:113
Belle2::ARICHLikelihood::getDetPhot
float getDetPhot() const
Return number of detected photons for a given particle.
Definition: ARICHLikelihood.h:92
Belle2::ARICHDQMModule::h_chHit
TH1 * h_chHit
The number of hits in each channel.
Definition: ARICHDQMModule.h:104
Belle2::Module::setPropertyFlags
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:210
Belle2::ExtHit
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:40
Belle2::ARICHDQMModule::endRun
virtual void endRun() override
End-of-run action.
Definition: ARICHDQMModule.cc:459
Belle2::ARICHDQMModule::h_trackPerEvent
TH1 * h_trackPerEvent
Number of tracks in ARICH per event (with p>0.5 GeV)
Definition: ARICHDQMModule.h:124
Belle2::ARICHGeoAerogelPlane::getAerogelTileID
unsigned getAerogelTileID(double x, double y) const
Get ID of aerogel tile containing point (x,y) (actually this is tile slot ID, as it is the same for a...
Definition: ARICHGeoAerogelPlane.cc:122
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ARICHDQMModule::m_momDnLim
double m_momDnLim
Lower momentum limit of tracks used in GeV (if set 0, no limit is applied)
Definition: ARICHDQMModule.h:140
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::ARICHDQMModule::h_mirrorThetaPhi
TH3 * h_mirrorThetaPhi
cherenkov theta vs phi for mirror reflected photons (for each mirror plate)
Definition: ARICHDQMModule.h:119
Belle2::ARICHDQMModule::h_hitsPerTrack
TH1 * h_hitsPerTrack
Average hits/track calculated from h_hits2D and h_track2D.
Definition: ARICHDQMModule.h:123
Belle2::ARICHDQMModule::m_debug
bool m_debug
debug
Definition: ARICHDQMModule.h:90
Belle2::Module::setReturnValue
void setReturnValue(int value)
Sets the return value for this module as integer.
Definition: Module.cc:222
Belle2::ARICHDQMModule::initialize
virtual void initialize() override
Initialize the Module.
Definition: ARICHDQMModule.cc:207
Belle2::ARICHDQMModule::m_arichEvents
bool m_arichEvents
process only events that have extrapolated hit in arich
Definition: ARICHDQMModule.h:92
Belle2::Module::addParam
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:562
Belle2::ARICHGeoDetectorPlane
Geometry parameters of ARICH photon detector plane.
Definition: ARICHGeoDetectorPlane.h:37
Belle2::ARICHDQMModule::h_chDigit
TH1 * h_chDigit
The number of raw digits in each channel.
Definition: ARICHDQMModule.h:101
Belle2::ARICHDQMModule::m_maxHits
int m_maxHits
exclude events with very large number of hits in arich
Definition: ARICHDQMModule.h:93
Belle2::ARICHDQMModule::h_ARICHEOccAfterInjLer
TH1 * h_ARICHEOccAfterInjLer
Histogram for Nr Entries (=Triggrs) for normalization after LER injection.
Definition: ARICHDQMModule.h:132
Belle2::Track
Class that bundles various TrackFitResults.
Definition: Track.h:35
Belle2::ARICHDQMModule::h_hapdDigit
TH1 * h_hapdDigit
The number of raw digits in each HAPD.
Definition: ARICHDQMModule.h:103
Belle2::StoreArray< MCParticle >
Belle2::ARICHGeoAerogelPlane
Geometry parameters of HAPD.
Definition: ARICHGeoAerogelPlane.h:37
Belle2::ARICHDQMModule::h_secTheta
TH1 * h_secTheta[6]
Detailed view of Cherenkov angle for each sector.
Definition: ARICHDQMModule.h:127
Belle2::ARICHDQMModule::h_chipDigit
TH1 * h_chipDigit
The number of raw digits in each ASIC chip.
Definition: ARICHDQMModule.h:102
Belle2::ARICHDQMModule::h_bits
TH1 * h_bits
Timing bits.
Definition: ARICHDQMModule.h:114
Belle2::ARICHGeoAerogelPlane::getRingDPhi
double getRingDPhi(unsigned iRing) const
Get phi (angle) distance between "phi" aluminum wall between aerogel tiles in i-th tile ring.
Definition: ARICHGeoAerogelPlane.h:280
Belle2::StoreArray::getEntries
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:226
Belle2::ARICHGeoDetectorPlane::getSlotPhi
double getSlotPhi(unsigned modID) const
Get phi (angle) position of module with given module ID number.
Definition: ARICHGeoDetectorPlane.h:229
Belle2::ARICHDQMModule::m_minHits
int m_minHits
exclude events with number of hits lower than this
Definition: ARICHDQMModule.h:94
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Belle2::ARICHDQMModule::h_tracks2D
TH2 * h_tracks2D
2D track distribution of whole ARICH
Definition: ARICHDQMModule.h:117
Belle2::ARICHDQMModule::beginRun
virtual void beginRun() override
Called when entering a new run.
Definition: ARICHDQMModule.cc:227
Belle2::ARICHDQMModule::h_aerogelHits3D
TH3 * h_aerogelHits3D
3D histogram of
Definition: ARICHDQMModule.h:118
Belle2::ARICHDQMModule::defineHisto
virtual void defineHisto() override
Function to define histograms.
Definition: ARICHDQMModule.cc:73
Belle2::ARICHDQMModule::h_secHapdHit
TH1 * h_secHapdHit[6]
The number of hits in each HAPDs of each sector.
Definition: ARICHDQMModule.h:111