Belle II Software  release-06-02-00
arichBtestModule.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 #include <arich/modules/arichBtest/arichBtestModule.h>
10 
11 //include <boost/format.hpp>
12 //include <boost/foreach.hpp>
13 
14 
15 // Framework - DataStore
16 #include <framework/datastore/StoreObjPtr.h>
17 #include <framework/dataobjects/EventMetaData.h>
18 #include <framework/datastore/StoreArray.h>
19 #include <framework/logging/Logger.h>
20 #include <framework/gearbox/Unit.h>
21 
22 
23 #include <boost/algorithm/string.hpp>
24 #include <boost/foreach.hpp>
25 
26 
27 // Hit classes
28 #include <arich/dataobjects/ARICHDigit.h>
29 #include <arich/dataobjects/ARICHAeroHit.h>
30 #include "arich/geometry/ARICHGeometryPar.h"
31 #include "arich/geometry/ARICHBtestGeometryPar.h"
32 
33 #include "arich/modules/arichBtest/arichBtestData.h"
34 
35 #include <TH1F.h>
36 #include <TH2F.h>
37 #include <TFile.h>
38 #include <TNtuple.h>
39 #include <TVector3.h>
40 #include <TAxis.h>
41 
42 // ifstream constructor.
43 #include <fstream>
44 
45 const char* record_strings[] = { "Event", "Begin_run", "Pause", "Resume", "End_run" };
46 
47 using namespace std;
48 
49 namespace Belle2 {
56 //-----------------------------------------------------------------
57 // Register the Module
58 //-----------------------------------------------------------------
59  REG_MODULE(arichBtest)
60 
61 //-----------------------------------------------------------------
62 // Implementation
63 //-----------------------------------------------------------------
64 
65 
66  arichBtestModule::arichBtestModule() : Module(), m_end(0), m_events(0), m_file(NULL), m_timestart(0), m_mwpc(NULL)
67  {
68  //Set module properties
69  setDescription("Module for the ARICH Beamtest data analysis. It creates track form the MWPC hits and reads the HAPD hits");
70 
71  //Parameter definition
72  addParam("outputFileName", m_outFile, "Output Root Filename", string("output.root"));
73  vector<string> defaultList;
74  addParam("runList", m_runList, "Data Filenames.", defaultList);
75  vector<int> defaultMask;
76  addParam("mwpcTrackMask", m_MwpcTrackMask, "Create track from MWPC layers", defaultMask);
77  m_fp = NULL;
78  addParam("beamMomentum", m_beamMomentum, "Momentum of the beam [GeV]", 0.0);
79 
80  for (int i = 0; i < 32; i++) m_tdc[i] = 0;
81 
82  }
83 
84  TNtuple* m_tuple;
85  TH1F* hapd[6];
86  TH1F* mwpc_tdc[4][5];
87  TH1F* mwpc_diff[4][2];
88  TH1F* mwpc_sum[4][2];
89  TH1F* mwpc_sum_cut[4][2];
90  TH1F* mwpc_residuals[4][2];
91  TH2F* mwpc_xy[4];
92  TH2F* mwpc_residualsz[4][2];
93  //TGraph* m_hapdmap;
94  //TGraph* m_el2pos;
95 
96  void arichBtestModule::initialize()
97  {
98  B2INFO("arichBtestModule::initialize()");
99  StoreArray<ARICHAeroHit> aeroHits;
100  StoreArray<ARICHDigit> digits;
101  aeroHits.registerInDataStore();
102  digits.registerInDataStore();
103 
104  m_file = new TFile(m_outFile.c_str(), "RECREATE");
105  m_tuple = new TNtuple("nt", "Btest data hits", "x:y:xrec:yrec:m:c:gx:gy");
106  char hapdName[256];
107  for (int i = 0; i < 6; i++) {
108  sprintf(hapdName, "hapd%d", i);
109  hapd[i] = new TH1F(hapdName, hapdName, 145, -0.5, 144.5);
110  }
111  char name[256];
112 
113  for (int i = 0; i < 4; i++) {
114  for (int k = 0; k < 2; k++) {
115  sprintf(name, "mwpc%d_a%d_", i, k);
116  mwpc_diff[i][k] = new TH1F(strcat(name, "diff"), name, 300, -150, 150);
117  sprintf(name, "mwpc%d_a%d_", i, k);
118  mwpc_sum[i][k] = new TH1F(strcat(name, "sum"), name, 200, -0.5, 199.5);
119  sprintf(name, "mwpc%d_a%d_", i, k);
120  mwpc_sum_cut[i][k] = new TH1F(strcat(name, "sum_cut"), name, 200, -0.5, 199.5);
121  sprintf(name, "mwpc%d_a%d_", i, k);
122  mwpc_residuals[i][k] = new TH1F(strcat(name, "resd"), name, 200, -100, 100);
123  sprintf(name, "mwpc%d_a%d_", i, k);
124  mwpc_residualsz[i][k] = new TH2F(strcat(name, "resd_z"), name, 200, -25, 25, 400, -1000, 1000);
125  }
126  for (int k = 0; k < 5; k++) {
127  sprintf(name, "mwpc%d_a%d_", i, k);
128  mwpc_tdc[i][k] = new TH1F(strcat(name, "tdc"), name, 1024, -1024, 1024);
129  }
130  sprintf(name, "mwpc%d_", i);
131  mwpc_xy[i] = new TH2F(strcat(name, "xy"), name, 120, -30, 30, 120, -30, 30);
132  }
133  //m_hapdmap = new TGraph("arich/modules/arichBtest/geometry/hapd.map");
134  //m_el2pos = new TGraph("arich/modules/arichBtest/geometry/hapdchmap_v0.dat");
135 
136  /*
137  arich::ARICHGeometryPar* _arichgp = arich::ARICHGeometryPar::Instance();
138 
139  ofstream dout;
140  dout.open ("ChannelCenterGlob.txt");
141  for (int i=0;i<6;i++){
142  for (int k=0;k<144;k++){
143  TVector3 r = _arichgp->getChannelCenterGlob(i + 1, k);
144  dout << r.x() << " " << r.y() << endl;
145  }
146  }
147  dout.close();
148  */
149  time(&m_timestart);
150  }
151 
152  void arichBtestModule::beginRun()
153  {
154 
155  B2INFO("arichBtestModule::beginRun()");
156 
157  StoreObjPtr<EventMetaData> eventMetaDataPtr;
158  B2INFO("arichBtestModule::eventMetaDataPtr run:" << eventMetaDataPtr->getRun());
159  B2INFO("arichBtestModule::eventMetaDataPtr exp:" << eventMetaDataPtr->getExperiment());
160 
161  ARICHBtestGeometryPar* _arichbtgp = ARICHBtestGeometryPar::Instance();
162  m_mwpc = _arichbtgp->getMwpc();
163 
164  static int first = 1;
165  if (first) {
166  m_runCurrent = m_runList.begin();
167  first = 0;
168  m_mwpc[0].Print();
169  m_mwpc[1].Print();
170  m_mwpc[2].Print();
171  m_mwpc[3].Print();
172 
173  }
174  m_end = 1;
175 
176  if (m_runCurrent < m_runList.end()) {
177  if (m_fp) {
178  m_eveList.push_back(m_events);
179  gzclose(m_fp);
180  }
181  m_fp = gzopen((*m_runCurrent).c_str(), "rb");
182  if (m_fp == NULL) {
183  B2INFO("Cannot open " << *m_runCurrent);
184  m_end = 1;
185  } else {
186  B2INFO("File opened " << *m_runCurrent);
187  m_end = 0;
188  }
189  }
190  m_events = 0;
191  }
192 
193  int arichBtestModule::skipdata(gzFile fp)
194  {
195  unsigned int u;
196  gzread(fp, &u, sizeof(unsigned int));
197  gzseek(fp, u * sizeof(unsigned int), SEEK_CUR);
198  return u + 1;
199  }
200 
201  void arichBtestModule::readmwpc(unsigned int* dbuf, unsigned int len, int print1290)
202  {
203 
204  const int unsigned MAXV1290 = 32;
205 
206  unsigned int edge;
207  unsigned int tdcch;
208  unsigned int tdcl;
209  unsigned int nhits;
210 
211  for (int i = 0; i < 32; i++) m_tdc[i] = 0XFFFF;
212  //TDC V1290
213  //for (int i=0;i<len;i++) printf("V1290 %d(%d) 0x%08x \n",i,len, dbuf[i],n);
214  for (unsigned int i = 0; i < len; i++) {
215  int rid = (dbuf[i] >> 27) & 0x1F;
216  switch (rid) {
217  case 0x18: // Filler
218  if (print1290) printf("Filler 0x%08x %u.data\n", dbuf[i], i);
219  break;
220  case 0x8: // Global Header
221  if (print1290) printf("Global Header 0x%08x %u.data\n", dbuf[i], i);
222  break;
223  case 0x10: // Global Trailer --- Last word of data
224  nhits = ((dbuf[i] >> 5) & 0xFFFF); //Word Count = bit 31...21< Word Count: 20 ... 5 > 4...0
225  if (print1290) printf("Global Trailer 0x%08x %u.data STATUS=0x%03x nhits=%u\n", dbuf[i], i, (dbuf[i] >> 24) & 0x7, nhits);
226 
227  if (nhits != len) {
228  if (print1290) printf("V1290 nhits!=len %u %u\n", nhits, len);
229  };
230  break;
231  case 0x11: // Global Trigger TimeTag
232  if (print1290) printf("Global Trigger TimeTag 0x%08x %u.data\n", dbuf[i], i);
233  break;
234  case 0x1: // TDC header
235  if (print1290) printf("TDC header 0x%08x %u.data evid=%u wc=%u\n", dbuf[i], i, (dbuf[i] >> 12) & 0xFFF, dbuf[i] & 0xFFF);
236 
237  break;
238  case 0x3: // TDC trailer
239  if (print1290) printf("TDC trailer 0x%08x %u.data\n", dbuf[i], i);
240 
241  break;
242  case 0x4: // TDC Error
243  if (print1290) printf("TDC Error 0x%08x %u.data ERR=0x%x\n", dbuf[i], i, dbuf[i] & 0x3FFF);
244 
245  break;
246  case 0x0 : // TDC data
247 
248  edge = (dbuf[i] >> 26) & 0x1;
249  tdcch = (dbuf[i] >> 21) & 0x1F;
250  tdcl = dbuf[i] & 0x1FFFFF;
251  if (tdcch < MAXV1290) {
252  //if (tdcch>15) gV1290[tdcch]->Fill(tdcl/40);
253  //if (tdcch>15) if (tdcl< trg[tdcch-16] && tdcl>16000 && tdcl<20000 ) trg[tdcch-16]=tdcl;
254  if (print1290)
255  printf("V1290 0x%08x %u. [ch=%2u] edge=%u data=%u \n", dbuf[i], i, tdcch, edge, tdcl);
256  m_tdc[tdcch] = tdcl / 40;
257  }
258 
259  break;
260 
261  default:
262  if (print1290) printf("Unknown 0x%08x %u.data\n", dbuf[i], i);
263 
264  break;
265 
266  }
267  }
268 
269  }
270 
271  int arichBtestModule::readhapd(unsigned int len, unsigned int* data)
272  {
273 
274  ARICHGeometryPar* _arichgp = ARICHGeometryPar::Instance();
275  ARICHBtestGeometryPar* _arichbtgp = ARICHBtestGeometryPar::Instance();
276  //-----------------------------------------------------
277 
278  int bmask = 0xF;
279 
280  for (int module = 0; module < 6; module++) {
281  int istart = 19 * module;
282  int istop = 19 * module + 18;
283  for (int i = istart; i < istop; i++) {
284  for (int j = 0; j < 8; j++) {
285  unsigned int kdata = data[i] & (bmask << (j * 4));
286  if (kdata) {
287  int channelID = j + 8 * (i - istart);
288 
289  if (_arichgp->isActive(module, channelID)) {
290 
291  hapd[module]->Fill(channelID);
292 
293  StoreArray<ARICHDigit> arichDigits;
294  double globalTime = 0;
295 
296  double rposx = 0, rposy = 0;
297  pair<int, int> eposhapd(_arichbtgp->GetHapdElectronicMap(module * 144 + channelID));
298  int channel = eposhapd.second;
299 
300  if ((channel < 108 && channel > 71) || channel < 36) channel = 108 - (int(channel / 6) * 2 + 1) * 6 +
301  channel;
302  else channel = 144 - (int((channel - 36) / 6) * 2 + 1) * 6 + channel - 36;
303  TVector2 loc = _arichgp->getChannelCenterLoc(channel);
304  if (abs(loc.X()) > 2.3 || abs(loc.Y()) > 2.3) continue;
305 
306  arichDigits.appendNew(module + 1, channel, globalTime);
307 
308  TVector3 rechit = _arichgp->getChannelCenterGlob(module + 1, channel);
309  pair<double, double> poshapd(_arichbtgp->GetHapdChannelPosition(module * 144 + channelID));
310  m_tuple ->Fill(-poshapd.first, poshapd.second, rechit.x(), rechit.y(), module, channelID, rposx, rposy);
311  }
312  }
313  }
314  }
315  }
316 
317  return len;
318  }
319 
320 
321  int arichBtestModule::getTrack(int mask, TVector3& r, TVector3& dir)
322  {
323  int retval = 0;
324  //const int trgch = 13;
325  const double t0 = 0;// m_tdc[trgch];
326  for (int i = 0; i < 4; i++) {
327  ARICHTracking* w = &m_mwpc[i];
328 
329  for (int k = 0; k < 4; k++) mwpc_tdc[i][k]->Fill(m_tdc[w->tdc[k]] - t0);
330  mwpc_tdc[i][4]->Fill(m_tdc[w->atdc] - t0);
331 
332 
333  for (int k = 0; k < 2; k++) { // axis x,y
334  w->status[k] = 1;
335  w->diff[k] = m_tdc[w->tdc[2 * k]] - m_tdc[w->tdc[2 * k + 1]];
336  w->sum[k] = m_tdc[w->tdc[2 * k + 1]] + m_tdc[w->tdc[2 * k]] - 2 * m_tdc[w->atdc];
337  mwpc_sum[i][k]->Fill(w->sum[k]);
338  if (w->sum[k] < w->cutll[k] || w->sum[k] > w->cutul[k]) continue;
339  mwpc_sum_cut[i][k]->Fill(w->sum[k]);
340  w->status[k] = 0;
341  w->reco[k] = w->diff[k] * w->slp[k] + w->offset[k];
342  w->reco[k] += w->pos[k];
343 
344  mwpc_diff[i][k]->Fill(w->diff[k]);
345  }
346 
347  w->reco[2] = w->pos[2]; // update z axis
348  if (!w->status[0] && !w->status[1]) mwpc_xy[i]->Fill(w->reco[0], w->reco[1]);
349 
350  if (mask & (1 << i)) {
351  if (!w->status[0] && !w->status[1]) {
352  // add point to a fitter
353  } else {
354  retval = 1;
355  }
356  }
357  }
358 
359  // replace by fitter
360  int ii[4] = {0, 1, 0, 0};
361  int ncnt = 0;
362  for (int i = 0; i < 4; i++) {
363  if (mask & (1 << i)) {
364  //B2INFO(ncnt << " " << i << " Mask " << mask);
365  ii[ncnt++] = i;
366  }
367  }
368  int i0 = ii[0];
369  int i1 = ii[1];
370 
371  r.SetXYZ(m_mwpc[i0].reco[1], m_mwpc[i0].reco[0], m_mwpc[i0].reco[2]);
372  dir.SetXYZ(m_mwpc[i1].reco[1] - m_mwpc[i0].reco[1],
373  m_mwpc[i1].reco[0] - m_mwpc[i0].reco[0],
374  m_mwpc[i1].reco[2] - m_mwpc[i0].reco[2]);
375 
376  dir = dir.Unit();
377 
378 // end replace by fitter
379  if (dir.z() != 0) {
380  for (int i = 0; i < 4; i++) {
381  ARICHTracking* w = &m_mwpc[i];
382  double l = (w->reco[2] - r.z()) / dir.z() ;
383  TVector3 rext = r + dir * l;
384  if (!w->status[0]) mwpc_residuals[i][0]->Fill(w->reco[0] - rext.y());
385  if (!w->status[1]) mwpc_residuals[i][1]->Fill(w->reco[1] - rext.x());
386 
387  TAxis* axis = mwpc_residualsz[i][1]->GetYaxis();
388  for (int k = 0; k < axis->GetNbins(); k++) {
389  double ll = (w->reco[2] + axis->GetBinCenter(k + 1) - r.z()) / dir.z();
390  TVector3 rextt = r + dir * ll;
391  mwpc_residualsz[i][0]->Fill(w->reco[0] - rextt.y(), axis->GetBinCenter(k + 1));
392  mwpc_residualsz[i][1]->Fill(w->reco[1] - rextt.x(), axis->GetBinCenter(k + 1));
393 
394  }
395  }
396  }
397 
398  return retval;
399  }
400 
401  int arichBtestModule::readdata(gzFile fp, int rec_id, int)
402  {
403 
404  unsigned int len, data[10000];
405  gzread(fp, &len, sizeof(unsigned int));
406  //if (print) printf( "[%3d] %d: ", len, rec_id );
407  gzread(fp, data, sizeof(unsigned int)*len);
408 
409  TVector3 r;
410  TVector3 dir;
411  if (rec_id == 1) {
412  readmwpc(data, len);
413  int retval = getTrack(*(m_MwpcTrackMask.begin()), r, dir);
414  //dir = TVector3(0,0,1);
415 
416  if (!retval) {
417  // global transf, add track to datastore
418  StoreArray<ARICHAeroHit> arichAeroHits;
419 
420  int particleId = 11;// geant4
421  dir *= m_beamMomentum * Unit::GeV;
422  r *= Unit::mm /*/ CLHEP::mm*/;
423  static ARICHBtestGeometryPar* _arichbtgp = ARICHBtestGeometryPar::Instance();
424  static TVector3 dr = _arichbtgp->getTrackingShift();
425 
426  r += dr;
427 
428  //----------------------------------------
429  // Track rotation
430  //
431  TRotation rot = _arichbtgp->getFrameRotation();
432  TVector3 rc = _arichbtgp->getRotationCenter();
433 
434  TVector3 rrel = rc - rot * rc;
435  r = rot * r + rrel;
436  dir = rot * dir;
437  r.SetX(-r.X()); dir.SetX(-dir.X());
438 
439  //
440  // end track rotation
441  //----------------------------------------
442  r[1] = -r.y();
443  dir[1] = -dir.y();
444  B2DEBUG(50, "-----------> " << rc.x() << " " << rc.y() << " " << rc.z() << "::::" << rrel.x() << " " << rrel.y() << " " <<
445  rrel.z() << " ----> R " << r.x() << " " << r.y() << " " << r.z() << " ----> S " << dir.x() << " " << dir.y() << " " <<
446  dir.z());
447 
448  // Add new ARIHCAeroHit to datastore
449  arichAeroHits.appendNew(particleId, r, dir);
450 
451  }
452  }
453 
454  if (rec_id == 2) {
455  readhapd(len, data);
456  }
457 
458  for (unsigned int j = 0; j < len; j++) {
459  //if( j%8==0 && j!= 0 ) printf( " " );
460  //printf( " %08x", data[j] );
461  //if( j%8==7 || j==len-1 ) putchar('\n');
462  }
463  return len + 1;
464  }
465 
466 
467 
468 
469  void arichBtestModule::event()
470  {
471 
472 
473 
474 
475 
476  if (!m_fp) return;
477  unsigned int hdr[2];
478  EventRec rec;
479  BeginRec beginrec;
480  BeginRec endrec;
481  static char msg[1024];
482  int type;
483 
484  const int sint = sizeof(unsigned int);
485  do {
486  if (gzread(m_fp, hdr, sint * 2) != 2 * sint && m_end) {
487  B2INFO("[" << m_events << "] End of file: " << *m_runCurrent);
488  ++m_runCurrent;
489  if (m_runCurrent == m_runList.end()) {
490  StoreObjPtr<EventMetaData> eventMetaDataPtr;
491  eventMetaDataPtr->setEndOfData();
492  return;
493  }
494  beginRun();
495  }
496 
497  } while (m_end && m_runCurrent != m_runList.end());
498  if (m_events % 1000 == 0) {
499  time_t m_time;
500  time(&m_time);
501  B2INFO("neve= [" << m_events << "] in " << (double)(m_time - m_timestart) / 60. << " min (" << int(
502  m_time - m_timestart) << "s) from " << *m_runCurrent);
503 
504  }
505  m_events++;
506 
507  type = hdr[0];
508  // len = hdr[1];
509 
510  gzseek(m_fp, 2 * sizeof(unsigned int), SEEK_CUR);
511  switch (type) {
512  case BEGIN_RECORD_TYPE: {
513  gzread(m_fp, &beginrec, sizeof(beginrec));
514  time_t t = beginrec.time;
515  sprintf(msg, "BeginRec run %u time %s", beginrec.runno, ctime(&t));
516  B2INFO(msg);
517  break;
518  }
519  case END_RECORD_TYPE: {
520  gzread(m_fp, &endrec, sizeof(endrec));
521  time_t t = endrec.time;
522  sprintf(msg, "EndRec run %u time %s", endrec.runno, ctime(&t));
523  B2INFO(msg);
524  break;
525  }
526  case EVENT_RECORD_TYPE: {
527  gzread(m_fp, &rec, sizeof(rec));
528  int print = !(rec.evtno % 10000);
529  time_t t = rec.time;
530  if (print) {
531  sprintf(msg, "EventRec run %u evt %u mstime %u, time %s", rec.runno, rec.evtno, rec.mstime, ctime(&t));
532  B2INFO(msg);
533  }
534  /* if you just want to jump to the end */
535  int pos = gztell(m_fp);
536  /* try to read inside data */
537  int buflen = rec.len - sizeof(rec) / sizeof(int);
538 
539  //if (print) printf ("%d: buflen\n",buflen);
540  int n[5] = { 0 };
541  int i(0), j(0);
542  while (i < buflen && j < 5) {
543  n[j] = readdata(m_fp, j, print);
544  // n[j] = skipdata( m_fp );
545  i += n[j];
546  j++;
547  }
548  if (gzseek(m_fp, pos + buflen * sizeof(unsigned int), SEEK_SET) < 0) B2ERROR("gzseek returns -1 ");
549  break;
550  }
551  default:
552  B2ERROR("IO error unknown record type " << type);
553  break;
554  }
555  if (gzeof(m_fp)) m_end = 1;
556  }
557 
558  void arichBtestModule::endRun()
559  {
560  B2INFO(" arichBtestModule: End Run !!!");
561 
562  m_file->Write();
563  m_file->Close();
564 
565  if (m_fp) {
566  gzclose(m_fp);
567  m_eveList.push_back(m_events);
568  }
569  }
570 
571  void arichBtestModule::terminate()
572  {
573  int j = 1;
574  BOOST_FOREACH(const string & fname, m_runList) {
575  B2INFO(m_eveList[j] << " events processed from file " << fname);
576  j++;
577  }
578  for (int i = 0; i < 4; i++) {
579  //ARICHTracking* w = &m_mwpc[i];
580  B2INFO(i << " a1=" << m_mwpc[i].tdc[0] << " a2=" << m_mwpc[i].tdc[1] << " a3=" << m_mwpc[i].tdc[2] << " a2=" << m_mwpc[i].tdc[3]
581  << " A=" << m_mwpc[i].atdc);
582  }
583 
584  }
585 
587 }
The Class for ARICH Beamtest Geometry Parameters.
The Class for ARICH Geometry Parameters.
Beamtest ARICH Geometry Tracking Class.
Base class for Modules.
Definition: Module.h:72
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
The UserTutorial module.
TH2F * mwpc_xy[4]
calculated x-y track positions
TH1F * mwpc_residuals[4][2]
residuals from mwpcs
TH1F * mwpc_sum_cut[4][2]
tdc sum from mwpcs, with sum cut applied
TH2F * mwpc_residualsz[4][2]
z-residuals from mwpcs
TH1F * mwpc_diff[4][2]
tdc difference from mwpcs
TH1F * mwpc_sum[4][2]
tdc sum from mwpcs
TH1F * hapd[6]
histogram of hits for each hapd
TNtuple * m_tuple
ntuple containing hapd hits
TH1F * mwpc_tdc[4][5]
tdc information from mwpcs
TVector3 getRotationCenter()
Get the rotation center of the Aerogel RICH frame.
std::pair< double, double > GetHapdChannelPosition(int)
Get the position of the HAPD channel.
std::pair< int, int > GetHapdElectronicMap(int)
Get the mapping of the electronic channel to the HAPD module nr and the channel number.
TVector3 getChannelCenterGlob(int modID, int chanID)
get center of chanID channel of modID detector module (in global coordinates)
TVector2 getChannelCenterLoc(int chID)
get center position of chID channel (in detector module local coordinates)
TRotation getFrameRotation()
Get the rotation matrix of the Aerogel RICH frame.
ARICHTracking * getMwpc()
Get the pointer of the tracking MWPCs.
bool isActive(int module, int channel)
check the activity of the channel
TVector3 getTrackingShift()
Get the tracking shift.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
Begin record structure for the beamtest data.
unsigned int runno
Run number.
unsigned int time
Timestamp of the run.
Event record structure for the beamtest data.