9 #include <arich/modules/arichBtest/arichBtestModule.h>
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>
23 #include <boost/algorithm/string.hpp>
24 #include <boost/foreach.hpp>
28 #include <arich/dataobjects/ARICHDigit.h>
29 #include <arich/dataobjects/ARICHAeroHit.h>
30 #include "arich/geometry/ARICHGeometryPar.h"
31 #include "arich/geometry/ARICHBtestGeometryPar.h"
33 #include "arich/modules/arichBtest/arichBtestData.h"
65 setDescription(
"Module for the ARICH Beamtest data analysis. It creates track form the MWPC hits and reads the HAPD hits");
68 addParam(
"outputFileName",
m_outFile,
"Output Root Filename", std::string(
"output.root"));
69 std::vector<std::string> defaultList;
71 std::vector<int> defaultMask;
76 for (
int i = 0; i < 32; i++)
m_tdc[i] = 0;
94 B2INFO(
"arichBtestModule::initialize()");
98 digits.registerInDataStore();
101 m_tuple =
new TNtuple(
"nt",
"Btest data hits",
"x:y:xrec:yrec:m:c:gx:gy");
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);
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);
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);
126 sprintf(name,
"mwpc%d_", i);
127 mwpc_xy[i] =
new TH2F(strcat(name,
"xy"), name, 120, -30, 30, 120, -30, 30);
151 B2INFO(
"arichBtestModule::beginRun()");
154 B2INFO(
"arichBtestModule::eventMetaDataPtr run:" << eventMetaDataPtr->getRun());
155 B2INFO(
"arichBtestModule::eventMetaDataPtr exp:" << eventMetaDataPtr->getExperiment());
160 static int first = 1;
177 m_fp = gzopen((*m_runCurrent).c_str(),
"rb");
192 gzread(fp, &u,
sizeof(
unsigned int));
193 gzseek(fp, u *
sizeof(
unsigned int), SEEK_CUR);
200 const int unsigned MAXV1290 = 32;
207 for (
int i = 0; i < 32; i++)
m_tdc[i] = 0XFFFF;
210 for (
unsigned int i = 0; i < len; i++) {
211 int rid = (dbuf[i] >> 27) & 0x1F;
214 if (print1290) printf(
"Filler 0x%08x %u.data\n", dbuf[i], i);
217 if (print1290) printf(
"Global Header 0x%08x %u.data\n", dbuf[i], i);
220 nhits = ((dbuf[i] >> 5) & 0xFFFF);
221 if (print1290) printf(
"Global Trailer 0x%08x %u.data STATUS=0x%03x nhits=%u\n", dbuf[i], i, (dbuf[i] >> 24) & 0x7, nhits);
224 if (print1290) printf(
"V1290 nhits!=len %u %u\n", nhits, len);
228 if (print1290) printf(
"Global Trigger TimeTag 0x%08x %u.data\n", dbuf[i], i);
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);
235 if (print1290) printf(
"TDC trailer 0x%08x %u.data\n", dbuf[i], i);
239 if (print1290) printf(
"TDC Error 0x%08x %u.data ERR=0x%x\n", dbuf[i], i, dbuf[i] & 0x3FFF);
244 edge = (dbuf[i] >> 26) & 0x1;
245 tdcch = (dbuf[i] >> 21) & 0x1F;
246 tdcl = dbuf[i] & 0x1FFFFF;
247 if (tdcch < MAXV1290) {
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;
258 if (print1290) printf(
"Unknown 0x%08x %u.data\n", dbuf[i], i);
274 unsigned int bmask = 0xF;
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));
283 int channelID = j + 8 * (i - istart);
285 if (_arichgp->
isActive(module, channelID)) {
287 hapd[module]->Fill(channelID);
290 double globalTime = 0;
292 double rposx = 0, rposy = 0;
294 int channel = eposhapd.second;
296 if ((channel < 108 && channel > 71) || channel < 36) channel = 108 - (int(channel / 6) * 2 + 1) * 6 +
298 else channel = 144 - (int((channel - 36) / 6) * 2 + 1) * 6 + channel - 36;
300 if (abs(loc.X()) > 2.3 || abs(loc.Y()) > 2.3)
continue;
302 arichDigits.
appendNew(module + 1, channel, globalTime);
306 m_tuple ->Fill(-poshapd.first, poshapd.second, rechit.X(), rechit.Y(), module, channelID, rposx, rposy);
322 for (
int i = 0; i < 4; i++) {
325 for (
int k = 0; k < 4; k++)
mwpc_tdc[i][k]->Fill(
m_tdc[w->tdc[k]] - t0);
329 for (
int k = 0; k < 2; k++) {
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];
334 if (w->sum[k] < w->cutll[k] || w->sum[k] > w->cutul[k])
continue;
337 w->reco[k] = w->diff[k] * w->slp[k] + w->offset[k];
338 w->reco[k] += w->pos[k];
343 w->reco[2] = w->pos[2];
344 if (!w->status[0] && !w->status[1])
mwpc_xy[i]->Fill(w->reco[0], w->reco[1]);
346 if (mask & (1 << i)) {
347 if (!w->status[0] && !w->status[1]) {
356 int ii[4] = {0, 1, 0, 0};
358 for (
int i = 0; i < 4; i++) {
359 if (mask & (1 << i)) {
376 for (
int i = 0; i < 4; 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());
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));
400 unsigned int len, data[10000];
401 gzread(fp, &len,
sizeof(
unsigned int));
403 gzread(fp, data,
sizeof(
unsigned int)*len);
430 TVector3 rrel = rc - rot * rc;
433 r.SetX(-r.X()); dir.SetX(-dir.X());
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() <<
" " <<
445 arichAeroHits.
appendNew(particleId, r, dir);
454 for (
unsigned int j = 0; j < len; j++) {
477 static char msg[1024];
480 const int sint =
sizeof(
unsigned int);
482 if (gzread(
m_fp, hdr, sint * 2) != 2 * sint &&
m_end) {
487 eventMetaDataPtr->setEndOfData();
497 B2INFO(
"neve= [" <<
m_events <<
"] in " << (
double)(m_time -
m_timestart) / 60. <<
" min (" <<
int(
506 gzseek(
m_fp, 2 *
sizeof(
unsigned int), SEEK_CUR);
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));
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));
522 case EVENT_RECORD_TYPE: {
523 gzread(
m_fp, &rec,
sizeof(rec));
524 int print = !(rec.evtno % 10000);
527 sprintf(msg,
"EventRec run %u evt %u mstime %u, time %s", rec.runno, rec.evtno, rec.mstime, ctime(&t));
531 int pos = gztell(
m_fp);
533 int buflen = rec.len -
sizeof(rec) /
sizeof(
int);
538 while (i < buflen && j < 5) {
544 if (gzseek(
m_fp, pos + buflen *
sizeof(
unsigned int), SEEK_SET) < 0) B2ERROR(
"gzseek returns -1 ");
548 B2ERROR(
"IO error unknown record type " << type);
556 B2INFO(
" arichBtestModule: End Run !!!");
570 BOOST_FOREACH(
const std::string & fname,
m_runList) {
571 B2INFO(
m_eveList[j] <<
" events processed from file " << fname);
574 for (
int i = 0; i < 4; 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);
The Class for ARICH Beamtest Geometry Parameters.
The Class for ARICH Geometry Parameters.
Beamtest ARICH Geometry Tracking Class.
void Print()
Debug printouts.
void setDescription(const std::string &description)
Sets the description of the module.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Accessor to arrays stored in the data store.
T * appendNew()
Construct a new T object at the end of the array.
Type-safe access to single objects in the data store.
static const double mm
[millimeters]
static const double GeV
Standard of [energy, momentum, mass].
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
arichBtestModule()
Constructor.
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 ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
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.