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