]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Test beam analysis class by Sylwester
authorcblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 5 Nov 2007 17:14:50 +0000 (17:14 +0000)
committercblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 5 Nov 2007 17:14:50 +0000 (17:14 +0000)
TRD/AliTRDtestBeam.cxx [new file with mode: 0644]
TRD/AliTRDtestBeam.h [new file with mode: 0644]
TRD/Macros/AliTRDanalyzeTestBeam.C [new file with mode: 0644]
TRD/TRDbaseLinkDef.h
TRD/libTRDbase.pkg

diff --git a/TRD/AliTRDtestBeam.cxx b/TRD/AliTRDtestBeam.cxx
new file mode 100644 (file)
index 0000000..443b7c4
--- /dev/null
@@ -0,0 +1,329 @@
+
+#include "AliTRDtestBeam.h"
+
+#include "AliTRDRawStreamTB.h"
+#include "AliRawReaderMemory.h"
+
+#include <iostream>
+#include <fstream>
+
+#include <sys/types.h>
+#include <sys/stat.h>
+#include <fcntl.h>
+#include <unistd.h>
+
+//#include <>
+
+ClassImp(AliTRDtestBeam)
+
+const Long_t AliTRDtestBeam::file_head_size = 544; // ?
+const Long_t AliTRDtestBeam::event_head_size = 68; //?
+const Long_t AliTRDtestBeam::ldc_head_size = 68; //?
+const Long_t AliTRDtestBeam::equip_head_size = 28; //
+const Int_t AliTRDtestBeam::vme_in =1; //VME event in
+const Int_t AliTRDtestBeam::sim_in =1; //Si-strips in
+
+//typedef char byte;
+
+//offsets in bytes
+const Int_t AliTRDtestBeam::pos_run = 20; //run nr. (in file and event header)
+const Int_t AliTRDtestBeam::pos_length = 0; //event/equip. length
+const Int_t AliTRDtestBeam::pos_eqid = 8;  //equipment id.
+const Int_t AliTRDtestBeam::pos_sioff = 12;  //Si data size offset (3 extra words!!!)      
+     
+using namespace std;
+
+//____________________________________________________________________________ 
+AliTRDtestBeam::AliTRDtestBeam() :
+  fDataStream(0),
+  fHeaderIsRead(0),
+  fEventCount(0),
+  fLimit(4), 
+  fCurrent(0),
+  fDdlOff(0),
+  fSiOff(0),
+  fQdcOff(0),
+  fDdlSize(0),
+  fFileHeader(0),
+  fEventHeader(0),
+  fEventData(0),
+  fNSi1(0),
+  fNSi2(0),
+  fCher(0),
+  fPb(0)
+{
+
+}
+//____________________________________________________________________________ 
+AliTRDtestBeam::AliTRDtestBeam(const char *filename) :
+  fDataStream(0),
+  fHeaderIsRead(0),
+  fEventCount(0),
+  fLimit(4), 
+  fCurrent(0),
+  fDdlOff(0),
+  fSiOff(0),
+  fQdcOff(0),
+  fDdlSize(0),
+  fFileHeader(0),
+  fEventHeader(0),
+  fEventData(0),
+  fNSi1(0),
+  fNSi2(0),
+  fCher(0),
+  fPb(0)
+{
+
+  fDataStream = new ifstream(filename, ifstream::in | ifstream::binary );
+  cout << fDataStream->is_open() << endl;
+  //fHeaderIsRead = kTRUE;
+  fHeaderIsRead = kTRUE;
+
+  fFileHeader = new Char_t[file_head_size];
+  fEventHeader = new Char_t[event_head_size];
+  fEventData = new Char_t[fLimit];
+}
+//____________________________________________________________________________ 
+
+Int_t AliTRDtestBeam::NextEvent() {
+  
+  Long_t data_size=0,ldc_off; //,ldc_id,ldc2_id;
+  Long_t ldc_size,eq_id; //,ev_l2;
+  Long_t event_nr,ev_l1;
+  Long_t word;
+  
+  if ( !fHeaderIsRead ) {
+    fDataStream->read(fFileHeader, file_head_size);
+    if(fDataStream->fail()) {
+      cerr << "Error reading file header! " << endl;   
+      return false;
+    }
+    cout  << " Run nr.  " << Int(pos_run, fFileHeader) << endl;    
+    fHeaderIsRead=kTRUE;
+  }
+
+  fDataStream->read(fEventHeader, event_head_size);
+  if(fDataStream->fail()) {
+    cerr << "End of file, Event " << fEventCount  << endl;     
+    return false;
+  }
+  
+  data_size = Int(pos_length, fEventHeader)-event_head_size; //?
+  event_nr = Int((4+pos_run), fEventHeader); //ev.nr.
+  //cout << " Event " << event_nr <<" size "<< data_size <<endl;
+  
+    if (event_nr <= fEventCount-1) { //watch-out ...event counter starts at 1?
+      cout << fEventCount << " End of file?, Event " << fEventCount << endl;   
+      return false;
+    }
+    //cout <<  "Run " << Int(pos_run, header)<< " , Event " <<event_nr <<endl;
+    
+    // enough space for data?
+    if (fLimit < data_size) {
+      delete[] fEventData;
+      fEventData = new Char_t[data_size];
+      fLimit = data_size;
+    }
+    
+    fDataStream->read(fEventData, data_size);
+    
+    if(fDataStream->fail()) {
+      cerr << "End of file, Event " << fEventCount; // << endl;        
+       return false;
+    }
+    
+    //cout  << " ...IDs (size) : ";
+    
+    ldc_off=0; // size of data from one DDL link
+    
+    for ( size_t k = 0; k < 2; k++ ) { // 2 LDCs (DDL & VME)
+      
+      ldc_size = Int(ldc_off+pos_length, fEventData); //
+      //ldc_size1=(ldc_size-ldc_head_size);
+      eq_id = Int(ldc_off+ldc_head_size+pos_eqid, fEventData);     
+      //cout  << eq_id <<" ("<<ldc_size<<") ";     
+      
+      ev_l1 = Int((4+ldc_off+pos_run), fEventData); //ev.nr.
+      if ( ev_l1 != event_nr ){
+       //cerr << "Eq_id " <<eq_id<<" event nr. mismatch? " << event_nr <<" / "<< ev_l1 <<" ...LDC data size (header:68) " <<ldc_size<<endl;
+      }
+      
+      if (eq_id == 1024) {  //DDL data
+       fDdlOff = ldc_off; //+ldc_head_size+equip_head_size + 32;
+       fDdlSize = ldc_size;
+      }
+      
+      if (eq_id == 550) {  //Si-strip data (+QDC)
+       //cout << "550" << endl;
+       fSiOff=ldc_off+ldc_head_size+equip_head_size+pos_sioff;
+       word = Int(fSiOff, fEventData);
+       Short_t LenSi1 = (word >> 16) & 0xffff;
+       Short_t LenSi2 = word & 0xffff;
+       fQdcOff=fSiOff+4*(LenSi1+LenSi2+1)+equip_head_size+4; 
+      } 
+      else if (eq_id == 1182) {  //QDC first...
+       //cout << "1182" << endl;
+       fQdcOff=ldc_off+ldc_head_size+equip_head_size+pos_sioff;
+       fSiOff=fQdcOff+equip_head_size+4;
+      }
+      
+      ldc_off=ldc_size;
+      
+    }
+    //cout << endl;
+
+    //cout << "DDL = " << fDdlOff << endl;
+    // cout << "Si  = " << fSiOff << endl;
+    //cout << "QDC = " << fQdcOff << endl;
+    
+    DecodeSi();
+
+    fEventCount++; //event counter
+    return true;
+}
+//____________________________________________________________________________
+Int_t AliTRDtestBeam::DecodeSi() {
+  
+  if (fSiOff < 0) return 0;
+  
+  // cout << "decoding Si data" << endl;
+
+  Long_t word;
+  
+  word=Int(fSiOff, fEventData);
+  fNSi1 = (word >> 16) & 0xffff;
+  fNSi2 = word & 0xffff;
+  
+  Int_t cSi=fSiOff; //   
+  for (int i = 0; i < fNSi1; i++) {
+    fSi1Address[i] =  ( Int(cSi, fEventData) >> 12 ) & 0x7ff;
+    fSi1Charge[i] = Int(cSi, fEventData)  & 0xfff;
+    cSi+=4;
+  }
+    
+  for (int i = 0; i < fNSi2; i++) {  //1,for Date!
+    fSi2Address[i] =  ( Int(cSi, fEventData) >> 12 ) & 0x7ff;
+    fSi2Charge[i] = Int(cSi, fEventData)  & 0xfff;
+    cSi+=4;
+  }  
+  
+  // reconstruction
+
+  int LenSiX = 640;
+
+  int qmaxX; int amaxX;
+  int qmaxY; int amaxY;
+  
+  qmaxX = 5;
+  qmaxY = 5;
+  amaxX = -1;
+  amaxY = -1+LenSiX;
+  for( int i = 0; i < GetNSi1(); i++ ) {
+    if (fSi1Address[i] == 0) continue; // noise
+   
+    if (fSi1Address[i] < LenSiX ) {
+      if( fSi1Charge[i] > qmaxX ) {
+       qmaxX = fSi1Charge[i];
+       amaxX = fSi1Address[i];
+      }
+    } else  {
+      if( fSi1Charge[i] > qmaxY ) {
+       qmaxY = fSi1Charge[i];
+       amaxY = fSi1Address[i];
+      }
+    }
+  }
+  
+  fX[0] = (float)(amaxX*0.05);  // [mm]
+  fY[0] = (float)((amaxY-LenSiX)*0.05);
+  fQx[0] = (float)qmaxX;
+  fQy[0] = (float)qmaxY;
+  
+  // 
+  qmaxX = 5;
+  qmaxY = 5;
+  amaxX = -1;
+  amaxY = -1+LenSiX;
+
+  for( int i = 0; i < GetNSi2(); i++ ) {
+    
+    if (fSi2Address[i] == 1279) continue; // noise
+    if (fSi2Address[i] == 0) continue;    // noise
+    
+    if(fSi2Address[i] < LenSiX) {
+      if( fSi2Charge[i] > qmaxX ) {
+       qmaxX = fSi2Charge[i];
+       amaxX = fSi2Address[i];
+      }
+    } else {
+      if( fSi2Charge[i] > qmaxY ) {
+       //if (fSi2Charge[i] > 50) cout << fSi2Charge[i] << " " << i << " " <<  fSi2Address[i] << endl;
+       qmaxY = fSi2Charge[i];
+       amaxY = fSi2Address[i];
+      }
+    }
+  }
+  
+  fX[1] = (float)(amaxX*0.05);  // [mm]
+  fY[1] = (float)((amaxY-LenSiX)*0.05);
+  fQx[1] = (float)qmaxX;
+  fQy[1] = (float)qmaxY;
+  
+  if (fQdcOff < 0) return 0;
+  word=Int(fQdcOff, fEventData);
+  fPb   = (Double_t)((word >> 16) & 0xFFF);
+  fCher = (Double_t)((word ) & 0xFFF);
+
+  //cout << fCher << " " << fPb << endl;
+  return 1;
+}
+//____________________________________________________________________________ 
+/**/
+AliTRDRawStreamTB *AliTRDtestBeam::GetTRDrawStream() {
+  
+  // needs AliTRDRawStreamTB  
+  //cout << "Chamber reader:" << (Int_t)(fEventData+fDdlOff) << " " << fDdlSize << endl;
+  //int ifout = open("dump.dat", O_WRONLY | O_TRUNC | O_CREAT);
+  //write(ifout, (void*)(fEventData+fDdlOff+16), fDdlSize);
+  //close(ifout);
+
+  AliRawReaderMemory *reader = new AliRawReaderMemory((UChar_t*)(fEventData+fDdlOff), (UInt_t)fDdlSize);
+  reader->SetEquipmentID(1024);
+  reader->ReadHeader();
+  AliTRDRawStreamTB::RawBufferMissAligned(kTRUE);
+  AliTRDRawStreamTB::SupressWarnings(kTRUE);
+  AliTRDRawStreamTB *tb = new AliTRDRawStreamTB(reader); 
+  tb->Init();
+  return tb;
+  /*
+    return 
+
+    AliRawReaderMemory *rmem = data->GetRawReader();
+    rmem->ReadHeader();
+    
+    AliTRDRawStreamTB tb(rmem);
+    tb.Init();
+    AliTRDRawStreamTB::SupressWarnings(kTRUE);
+    
+  */
+}
+/**/
+//____________________________________________________________________________ 
+
+Int_t AliTRDtestBeam::Int(Int_t i, Char_t *start) {
+  
+  bool swap = kFALSE;
+
+  if(swap) {
+    char *q=(char*)(start+i); 
+    char p[] = {q[3], q[2], q[1], q[0]};
+    return *((int*) p);
+  } else return *((int*)(start+i));
+}
+
+//____________________________________________________________________________ 
diff --git a/TRD/AliTRDtestBeam.h b/TRD/AliTRDtestBeam.h
new file mode 100644 (file)
index 0000000..8f17b38
--- /dev/null
@@ -0,0 +1,116 @@
+#ifndef AliTRDtestBeam_h
+#define AliTRDtestBeam_h
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/*
+The class to read the test beam 2007 data
+*/
+
+#define MAX_SI 2000
+
+#include "TObject.h"
+
+class AliTRDrawStreamV2;
+class AliTRDRawStreamTB;
+using namespace std;
+
+class AliTRDtestBeam: public TObject {
+
+public:
+  AliTRDtestBeam();       // ctor
+  AliTRDtestBeam(const char *filename); // constructor
+  //AliTRDtestBeam(const AliTRDtestBeam &tb);  
+  //AliTRDtestBeam& operator = (const AliTRDtestBeam& tb) {return *this}
+  virtual ~AliTRDtestBeam() {;} // dtor
+
+  Int_t NextEvent();
+  
+  AliTRDRawStreamTB *GetTRDrawStream(); // needs RawStreamTB
+
+  // silicon
+  Short_t GetNSi1() {return fNSi1;}
+  Short_t GetNSi2() {return fNSi2;}
+  
+  Int_t GetSi1Address(Int_t i) {return (i<fNSi1)? fSi1Address[i] : -1;};
+  Int_t GetSi2Address(Int_t i) {return (i<fNSi2)? fSi2Address[i] : -1;};
+  
+  Int_t GetSi1Charge(Int_t i) {return (i<fNSi1)? fSi1Charge[i] : -1;};
+  Int_t GetSi2Charge(Int_t i) {return (i<fNSi2)? fSi1Charge[i] : -1;};
+  
+  Double_t GetX(Int_t n) {return (n<2)? fX[n] : -1;}
+  Double_t GetY(Int_t n) {return (n<2)? fY[n] : -1;}
+  Double_t GetQx(Int_t n) {return (n<2)? fQx[n] : -1;}
+  Double_t GetQy(Int_t n) {return (n<2)? fQy[n] : -1;}
+
+  // calo
+  Double_t GetCher() {return fCher;}
+  Double_t GetPb() {return fPb;}
+
+protected:
+  
+  ifstream *fDataStream;
+  
+  Bool_t fHeaderIsRead;
+  Int_t fEventCount;
+  Int_t fLimit; // = 4
+  Int_t fCurrent;
+
+  Int_t fDdlOff;
+  Int_t fSiOff;
+  Int_t fQdcOff;
+  Int_t fDdlSize;
+
+  Char_t *fFileHeader;
+  Char_t *fEventHeader;
+  Char_t *fEventData;
+
+  // silicon data
+  
+  Short_t fNSi1;
+  Short_t fNSi2;
+  
+  Int_t fSi1Address[MAX_SI];
+  Int_t fSi2Address[MAX_SI];
+  
+  Int_t fSi1Charge[MAX_SI];
+  Int_t fSi2Charge[MAX_SI];
+  
+  // reconstructed Silicon data 
+
+  Double_t fX[2];
+  Double_t fY[2];
+  Double_t fQx[2];
+  Double_t fQy[2];
+
+  // cherenkov glass
+  Double_t fCher;
+  Double_t fPb;
+  
+
+  // data reading
+  
+  Int_t Int(Int_t i, Char_t *start);
+  Int_t DecodeSi();
+
+  //
+  static const Long_t file_head_size; //= 544; // ?
+  static const Long_t event_head_size; // = 68; //?
+  static const Long_t ldc_head_size; // = 68; //?
+  static const Long_t equip_head_size; // = 28; //
+  static const Int_t vme_in; //=1; //VME event in
+  static const Int_t sim_in; //=1; //Si-strips in
+  
+  //typedef char byte;
+  
+  //offsets in bytes
+  static const Int_t pos_run; // = 20; //run nr. (in file and event header)
+  static const Int_t pos_length; // = 0; //event/equip. length
+  static const Int_t pos_eqid; // = 8;  //equipment id.
+  static const Int_t pos_sioff; // = 12;  //Si data size offset (3 extra words!!!)
+
+  ClassDef(AliTRDtestBeam,1)  // description 
+};
+
+#endif // AliTRDQADatamaker_H
+
diff --git a/TRD/Macros/AliTRDanalyzeTestBeam.C b/TRD/Macros/AliTRDanalyzeTestBeam.C
new file mode 100644 (file)
index 0000000..25cfead
--- /dev/null
@@ -0,0 +1,269 @@
+
+void AliTRDanalyzeTestBeam(Int_t run, Int_t begin, Int_t end) {
+  
+  gROOT->SetStyle("Plain");
+  //gStyle->SetPadTopMargin(0.02);
+  //gStyle->SetPadRightMargin(0.02);
+  //gStyle->SetPadLeftMargin(0.1);
+
+  gStyle->SetPalette(1);
+  gStyle->SetOptStat(0);
+  gStyle->SetOptDate();
+  
+  TGaxis::SetMaxDigits(3);
+
+  TStopwatch st;
+  st.Start();
+
+  const Int_t N = 640;
+
+  // declare histograms
+  const int nbins = 100;
+  const int nstart = 0;
+
+  TH1D *mSi1L = new TH1D("mSi1L", ";number of Si1 fired pads", 50, nstart-0.5, nstart+nbins-0.5);
+  TH1D *mSi2L = new TH1D("mS21L", ";number of Si2 fired pads", 50, nstart-0.5, nstart+nbins-0.5);
+  TProfile *mSi1ChP = new TProfile("mSi1ChP",";Si1 pad number;signal", 1280, -0.5, 1279.5, 0, 200, "s");
+  TProfile *mSi2ChP = new TProfile("mS21ChP",";Si2 pad number;signal", 1280, -0.5, 1279.5, 0, 200, "s");
+  
+  TH1D *mSi1N = new TH1D("mSi1N", "Noise Dist Si1;ADC", 100, 0, 50);
+  TH1D *mSi2N = new TH1D("mSi2N", "Noise Dist Si2;ADC", 100, 0, 50);
+  
+  TH1D *mSiCh[4];
+  mSiCh[0] = new TH1D("mSi1ChX", ";Si1X pad amplitude (ADC)", 250, -0.5, 249.5);
+  mSiCh[1] = new TH1D("mSi1ChY", ";Si1Y pad amplitude (ADC)", 250, -0.5, 249.5);
+  mSiCh[2] = new TH1D("mSi2ChX", ";Si2X pad amplitude (ADC)", 250, -0.5, 249.5);
+  mSiCh[3] = new TH1D("mSi2ChY", ";Si2Y pad amplitude (ADC)", 250, -0.5, 249.5);
+  
+  TH1D *mSiFullCh[4];
+  mSiFullCh[0] = new TH1D("mSi1fChX", "Si1X;max amplitude (ADC)", 300, -0.5, 299.5);
+  mSiFullCh[1] = new TH1D("mSi1fChY", "Si1Y;max amplitude (ADC)", 300, -0.5, 299.5);
+  mSiFullCh[2] = new TH1D("mSi2fChX", "Si2X;max amplitude (ADC)", 300, -0.5, 299.5);
+  mSiFullCh[3] = new TH1D("mSi2fChY", "Si2Y;max amplitude (ADC)", 300, -0.5, 299.5);  
+
+  TH2D *mPos[2];
+  mPos[0] = new TH2D("posSi1", ";x Si1 (mm);y Si1 (mm)", 128, 0, 32, 128, 0, 32);
+  mPos[1] = new TH2D("posSi2", ";x Si2 (mm);y Si2 (mm)", 128, 0, 32, 128, 0, 32);
+
+  TH2D *mCor[2];
+  mCor[0] = new TH2D("corX", ";Si1 X (mm);Si2 X (mm)", 128, 0, 32, 128, 0, 32);
+  mCor[1] = new TH2D("corY", ";Si1 Y (mm);Si2 Y (mm)", 128, 0, 32, 128, 0, 32);
+
+  TH2D *mChCor[2];
+  mChCor[0] = new TH2D("ChCorSi1", ";Si1 amp X;Si1 amp Y", 100, 0, 200, 100, 0, 200);
+  mChCor[1] = new TH2D("ChCorSi2", ";Si2 amp X;Si2 amp Y", 100, 0, 200, 100, 0, 200);
+  
+  gStyle->SetOptStat(11);
+  TH1D *mPb   = new TH1D("mPb", ";Amp Pb (ADC)", 150, -0.5, 4499.5);
+  TH1D *mCher = new TH1D("mCher", ";Amp Cherenkov (ADC)", 150, -0.5, 4499.5);
+  TH2D *mPbCher = new TH2D("mPbCher", ";amp Cherenkov;amp Pb", 150, -0.5, 4499.5, 150, -0.5, 4599.5);
+  // gStyle->SetOptStat(0);
+
+  // needed by the AliTRDRawStreamTB
+  AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
+  AliTRDcalibDB *calib = AliTRDcalibDB::Instance();
+  calib->SetRun(0);
+
+  // TRD data monitoring
+  TH1D *mDet = new TH1D("mDet", ";chamber", 20, -0.5, 19.5);
+  TH1D *mROB = new TH1D("mROB", ";rob", 20, -0.5, 19.5);
+  TH1D *mTRDsig = new TH1D("mTRDsig", ";trdSignal", 100, -0.5, 299.5);
+  
+  //AliLog::SetClassDebugLevel("AliTRDRawStreamTB", 10);
+
+  int counter = 0;
+  // for(Int_t run = 365; run < 386; run++) {
+  //for(Int_t run = 369; run < 382; run++) {
+  //for(int run = 387; run < 389; run++) {
+
+  // if (run == 389) continue;
+    cout << run << endl;
+
+    for(Int_t fn=begin; fn<end; fn++) {
+      
+      // connect to data 
+      //const char *base="/Users/radomski/data/1GeV/";
+      const char *base = "./";
+      const char *filename = Form("%s/run%d_gdc_daq09.%03d.raw",base,run,fn);
+      if (gSystem->AccessPathName(filename)) continue;
+      cout << filename << endl;
+
+      AliTRDtestBeam *data = new AliTRDtestBeam(filename);
+      
+      // process data
+      while (data->NextEvent()) {
+      
+       if (!(counter%1000)) cout << "Event = " << counter << endl;
+       counter++;
+       
+       /*
+       AliTRDRawStreamTB *tb = data->GetTRDrawStream();
+       while(tb->Next()) {
+         mROB->Fill(tb->GetROB());
+         mDet->Fill(tb->GetDet());
+         int *sig = tb->GetSignals();
+         mTRDsig->Fill(sig[0]);
+         mTRDsig->Fill(sig[1]);
+         mTRDsig->Fill(sig[2]);
+       }
+       delete tb;
+       */
+
+       mCher->Fill(data->GetCher());
+       mPb->Fill(data->GetPb());
+       mPbCher->Fill(data->GetCher(), data->GetPb());
+
+       mSi1L->Fill(data->GetNSi1());
+       mSi2L->Fill(data->GetNSi2());
+       
+       for(int i=0; i<data->GetNSi1(); i++) {
+         Int_t q = data->GetSi1Charge(i);
+         Int_t a = data->GetSi1Address(i);
+         if (a == 0) continue; // noisy channels
+         mSi1ChP->Fill(a, q);
+         if (a < N) mSiCh[0]->Fill(q);
+         else mSiCh[1]->Fill(q);              
+       }
+       
+       for(int i=0; i<data->GetNSi2(); i++) {
+         Int_t q = data->GetSi2Charge(i);
+         Int_t a = data->GetSi2Address(i);
+         if (a == 0 || a == 1279) continue; // noisy channels
+         mSi2ChP->Fill(a, q);
+         if (a < N) mSiCh[2]->Fill(q);
+         else mSiCh[3]->Fill(q);              
+       }
+       
+       mSiFullCh[0]->Fill(data->GetQx(0));
+       mSiFullCh[1]->Fill(data->GetQy(0));
+       mSiFullCh[2]->Fill(data->GetQx(1));
+       mSiFullCh[3]->Fill(data->GetQy(1));
+       
+       for(int k=0; k<2; k++)
+         mChCor[k]->Fill(data->GetQx(k), data->GetQy(k));
+       
+       /*
+         if (data->GetQx(0) < 20) continue;
+         if (data->GetQx(1) < 20) continue;
+         if (data->GetQy(0) < 20) continue;
+         if (data->GetQy(1) < 20) continue;
+       */
+       
+       for(int k=0; k<2; k++)
+         mPos[k]->Fill(data->GetX(k), data->GetY(k));
+      
+       mCor[0]->Fill(data->GetX(0), data->GetX(1));
+       mCor[1]->Fill(data->GetY(0), data->GetY(1));
+      }
+      
+      delete data;
+    }
+    //}
+
+  // process histograms
+  for(int i=1; i<1281; i++) mSi1N->Fill(mSi1ChP->GetBinError(i));
+  for(int i=1; i<1281; i++) mSi2N->Fill(mSi2ChP->GetBinError(i));
+
+  // display
+  cout << "Number of Events = " << counter << endl;
+  
+  /**/
+  TCanvas *c = new TCanvas("siliconSignal", "silicon signal");
+  c->Divide(2,2, 0.01, 0.01);
+  c->cd(1);
+  mSi1L->Draw();
+  c->cd(2);
+  mSi2L->Draw();
+  c->cd(3);
+  mSi1ChP->Draw();
+  c->cd(4);
+  mSi2ChP->Draw();
+  /* */
+
+  /**/
+  c = new TCanvas();
+  c->Divide(2,2,0.01,0.01);
+  c->cd(1);
+  mSi1N->Draw();
+  c->cd(2);
+  mSi2N->Draw();
+  /**/
+
+  /**/
+  // pads
+  c = new TCanvas("siPads", "Silicon Pads");
+  c->Divide(2,2, 0.01, 0.01);
+  for(int i=0; i<4; i++) {
+    c->cd(i+1);
+    gPad->SetLogy();
+    mSiCh[i]->Draw();
+  }
+  
+  // clusters
+  c = new TCanvas("siCluster", "silicon clusters");
+  c->Divide(2,2, 0.01, 0.01);
+  for(int i=0; i<4; i++) {
+    c->cd(i+1);
+    gPad->SetLogy();
+    mSiFullCh[i]->Draw();
+  }
+
+  // position and correlation
+  c = new TCanvas("siPosition", "reconstructed position");
+  c->Divide(2,2, 0.01, 0.01);
+  for(int i=0; i<2; i++) {
+    c->cd(1+i);
+    mPos[i]->Draw("col");
+    c->cd(3+i);
+    mCor[i]->Draw("col");  
+  }
+  
+  c = new TCanvas("siCharge", "si charge correlation");
+  c->Divide(2,2, 0.01, 0.01);
+  c->cd(1);
+  gPad->SetLogz();
+  mChCor[0]->Draw("col");
+  c->cd(2);
+  gPad->SetLogz();
+  mChCor[1]->Draw("col");
+  /**/
+
+  new TCanvas();
+  gPad->SetLogy();
+  mCher->Draw();
+  
+  // electron sample
+  int bin = mCher->FindBin(500.);
+  double ef = (mCher->Integral(bin, 151)/ mCher->GetSum());
+  TLatex *l = new TLatex(2e3, 0.02*mCher->GetSum(), Form("Electron fraction = %.2f ", ef));
+  l->Draw();
+  
+  new TCanvas();
+  gPad->SetLogy();
+  mPb->Draw();
+
+  new TCanvas();
+  gPad->SetLogz();
+  mPbCher->Draw("colz");
+  
+  /*
+  c = new TCanvas();
+  c->Divide(2,2,0.01, 0.01);
+  c->cd(1);
+  gPad->SetLogy();
+  mTRDsig->Draw();
+
+  c->cd(2);
+  mROB->Draw();
+
+  c->cd(3);
+  mDet->Draw();
+  */
+  
+  st.Stop();
+  st.Print();
+}
+
+//Int_t SelectEvent() {
+//}
index defda61ccb157c09fee543b704fc754e8b629500..8289d85ed4fa0b19f7df46400940d954d761bab4 100644 (file)
@@ -34,6 +34,8 @@
 #pragma link C++ class  AliTRDRawStreamV2+;
 #pragma link C++ class  AliTRDRawStreamTB+;
 
+#pragma link C++ class  AliTRDtestBeam+;
+
 #pragma link C++ class  AliTRDCommonParam+;
 #pragma link C++ class  AliTRDfeeParam+;
 
index 425aa23cef67008d82074aa25a46adec31f0c234..f1ce07b121b24bc68be0bba103eaae02ab319d68 100644 (file)
@@ -15,6 +15,7 @@ SRCS= AliTRDarrayI.cxx \
       AliTRDRawStream.cxx \
       AliTRDRawStreamV2.cxx \
       AliTRDRawStreamTB.cxx \
+      AliTRDtestBeam.cxx \
       AliTRDCommonParam.cxx \
       AliTRDfeeParam.cxx \
       AliTRDcalibDB.cxx \