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