]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDRaw2Digits.C
Coding conventions
[u/mrichter/AliRoot.git] / TRD / AliTRDRaw2Digits.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TH1.h>
3 #include <TCanvas.h>
4 #include "AliTRDrawData.h"
5 #include "AliTRDdigitsManager.h"
6 #include "AliTRDdigit.h"
7 #include "AliTRDgeometryFull.h"
8 #include "AliTRDparameter.h"
9 #include "AliTRDmatrix.h"
10 #include "AliRawReaderFile.h"
11 #include "AliRunLoader.h"
12 #endif
13
14 void AliTRDRaw2Digits(Int_t iEvent = 0, Int_t iDet = 0)
15 {
16
17   AliTRDrawData *raw = new AliTRDrawData();
18   raw->SetDebug(2);
19   AliRawReaderFile rawReader(iEvent);
20
21   AliTRDdigitsManager *digitsManagerRaw = raw->Raw2Digits(&rawReader);
22
23   // The geometry object
24   AliTRDgeometryFull *geo = new AliTRDgeometryFull();
25
26   // The parameter object
27   AliTRDparameter    *par = new AliTRDparameter("TRDparameter"
28                                                ,"TRD parameter class");
29
30   // Print the event and detector number
31   cout << " iEvent = " << iEvent << endl;
32   cout << " iDet = " << iDet << endl;
33
34   // Define the detector matrix for one chamber
35   const Int_t iSec = geo->GetSector(iDet);
36   const Int_t iCha = geo->GetChamber(iDet);
37   const Int_t iPla = geo->GetPlane(iDet);
38   Int_t  rowMax = par->GetRowMax(iPla,iCha,iSec);
39   Int_t  colMax = par->GetColMax(iPla);
40   Int_t timeMax = par->GetTimeMax();
41   Int_t timeTotal = par->GetTimeTotal();
42   cout << "Geometry: rowMax = "  <<  rowMax
43                 << " colMax = "  <<  colMax
44                 << " timeMax = " << timeMax
45                 << " timeTotal = " << timeTotal << endl;
46   AliTRDmatrix *matrix = new AliTRDmatrix(rowMax,colMax,timeTotal,iSec,iCha,iPla);
47
48   AliRunLoader* rl = AliRunLoader::Open("galice.root");
49   AliLoader* loader = rl->GetLoader("TRDLoader");
50   rl->GetEvent(iEvent);
51   loader->LoadDigits();
52   
53   // Create the digits manager
54   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
55   digitsManager->SetDebug(1);
56
57   // Read the digits from the file
58   digitsManager->ReadDigits(loader->TreeD());
59
60   TH1F *DigDiff = new TH1F("DigDiff","DigDiff",100,-10,+10);
61   Int_t DigAmpRaw, DigAmp;
62
63   // Loop through the detector pixel
64   for (Int_t time = 0; time < timeTotal; time++) {
65     for (Int_t  col = 0;  col <  colMax;  col++) {
66       for (Int_t  row = 0;  row <  rowMax;  row++) {
67
68         AliTRDdigit* digit = digitsManagerRaw->GetDigit(row,col,time,iDet);
69         
70         matrix->SetSignal(row,col,time,digit->GetAmp());
71
72         DigAmpRaw = digit->GetAmp();
73
74         digit = digitsManager->GetDigit(row,col,time,iDet);
75
76         DigAmp = digit->GetAmp();
77
78         DigDiff->Fill((Float_t)DigAmp-(Float_t)DigAmpRaw);
79
80         delete digit;
81
82       }
83     }
84   }
85
86   // Display the detector matrix
87   matrix->Draw();
88   matrix->ProjRow();
89   matrix->ProjCol();
90   matrix->ProjTime();
91
92   TCanvas *c1 = new TCanvas("c1","Canvas 1",10,10,600,500);
93   DigDiff->Draw();
94
95   printf("Number of digits expected = %d \n",timeTotal*rowMax*colMax);
96
97 }
98
99