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