1 Int_t AliTRDanalyzeDigits()
9 const Int_t kNpla = AliTRDgeometry::Nplan();
12 cout << "<AliTRDanalyzeDigits> No AliRun object found" << endl;
16 Int_t nPart = gAlice->GetEvent(0);
18 // Get the pointer to the TRD detector
19 AliTRD *trd = (AliTRD *) gAlice->GetDetector("TRD");
21 cout << "<AliTRDanalyzeDigits> No TRD detector found" << endl;
26 // Get the digitizer object
27 TFile *file = (TFile *) gROOT->GetListOfFiles()->FindObject("TRD_test.root");
28 AliTRDdigitizer *digitizer = (AliTRDdigitizer *) file->Get("digitizer");
30 cout << "<AliTRDanalyzeDigits> No digitizer object found" << endl;
35 // Define the histograms
36 Int_t adcRange = ((Int_t) digitizer->GetADCoutRange());
37 TH1F *hAmpAll = new TH1F("hAmpAll" ,"Amplitude of the digits (all)"
38 ,adcRange+1,-0.5,((Float_t) adcRange)+0.5);
39 TH1F *hAmpEl = new TH1F("hAmpEl" ,"Amplitude of the digits (electrons)"
40 ,adcRange+1,-0.5,((Float_t) adcRange)+0.5);
41 TH1F *hAmpPi = new TH1F("hAmpPi" ,"Amplitude of the digits (pions)"
42 ,adcRange+1,-0.5,((Float_t) adcRange)+0.5);
43 TH1F *hAmpNoise = new TH1F("hAmpNoise","Amplitude of the digits (noise)"
46 // Get the pointer to the geometry object
49 geo = trd->GetGeometry();
52 cout << "<AliTRDanalyzeDigits> No TRD geometry found" << endl;
57 // Create the digits manager
58 AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
60 digitsManager->Open("TRD_test.root");
62 // Read the digits from the file
63 if (!(digitsManager->ReadDigits())) {
64 cout << "<AliTRDanalyzeDigits> Cannot read the digits" << endl;
71 Int_t countDigits = 0;
72 Int_t iSec = trd->GetSensSector();
73 Int_t iCha = trd->GetSensChamber();
74 Int_t timeMax = geo->GetTimeTotal();
76 TProfile *hAmpTimeEl = new TProfile("hAmpTimeEl","Amplitude of the digits (electrons)"
77 ,timeMax,-0.5,((Double_t) timeMax)-0.5);
78 TProfile *hAmpTimePi = new TProfile("hAmpTimePi","Amplitude of the digits (pions)"
79 ,timeMax,-0.5,((Double_t) timeMax)-0.5);
81 // Loop over all planes
82 for (Int_t iPla = 0; iPla < kNpla; iPla++) {
84 Int_t iDet = geo->GetDetector(iPla,iCha,iSec);
85 Int_t rowMax = geo->GetRowMax(iPla,iCha,iSec);
86 Int_t colMax = geo->GetColMax(iPla);
89 matrix = new AliTRDmatrix(rowMax,colMax,timeMax,iSec,iCha,iPla);
92 // Loop through the detector pixel
93 for (Int_t time = 0; time < timeMax; time++) {
94 for (Int_t col = 0; col < colMax; col++) {
95 for (Int_t row = 0; row < rowMax; row++) {
97 AliTRDdigit *digit = digitsManager->GetDigit(row,col,time,iDet);
98 Int_t amp = digit->GetAmp();
99 Int_t track0 = digitsManager->GetTrack(0,row,col,time,iDet);
100 Int_t track1 = digitsManager->GetTrack(1,row,col,time,iDet);
101 TParticle *particle = 0;
103 particle = gAlice->Particle(track0);
109 matrix->SetSignal(row,col,time,amp);
118 hAmpNoise->Fill(amp);
122 if ((particle) && (particle->GetPdgCode() == 11) && (track1 < 0)) {
124 hAmpTimeEl->Fill(time,amp);
128 if ((particle) && (particle->GetPdgCode() == -211) && (track1 < 0)) {
130 hAmpTimePi->Fill(time,amp);
141 cout << "<AliTRDanalyzeDigits> Found " << countDigits << " digits in total" << endl;
143 // Display the detector matrix
149 TCanvas *cDigits = new TCanvas("cDigits","AliTRDanalyzeDigits",50,50,600,800);
150 cDigits->Divide(2,3);
153 hAmpAll->SetXTitle("Amplitude (ADC-channels)");
154 hAmpAll->SetYTitle("Entries");
158 hAmpNoise->SetXTitle("Amplitude (ADC-channels)");
159 hAmpNoise->SetYTitle("Entries");
163 hAmpEl->SetXTitle("Amplitude (ADC-channels)");
164 hAmpEl->SetYTitle("Entries");
168 hAmpPi->SetXTitle("Amplitude (ADC-channels)");
169 hAmpPi->SetYTitle("Entries");
172 hAmpTimeEl->SetXTitle("Timebin number");
173 hAmpTimeEl->SetYTitle("Mean amplitude");
174 hAmpTimeEl->Draw("HIST");
176 hAmpTimePi->SetXTitle("Timebin number");
177 hAmpTimePi->SetYTitle("Mean amplitude");
178 hAmpTimePi->Draw("HIST");
180 TFile *fileOut = new TFile("digits_test.root","RECREATE");