]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFanalyzeSDigitsV2.C
Fixed (ininfluential) bug on TestHit method
[u/mrichter/AliRoot.git] / TOF / AliTOFanalyzeSDigitsV2.C
1 Int_t AliTOFanalyzeSDigitsV2(TString headersFile, Int_t ndump=15, Int_t iEvNum=0)
2 {
3   // V2: to be used with transient fSDigits data member
4   //
5   // Analyzes the TOF sdigits and fills QA-histograms 
6   // report problems to pierella@bo.infn.it
7   // iEvNum=0 means all events in the file
8
9   // Dynamically link some shared libs
10   if (gClassTable->GetID("AliRun") < 0) {
11     gROOT->LoadMacro("loadlibs.C");
12     loadlibs();
13   } else {
14     delete gAlice;
15     gAlice = 0;
16   }
17   
18   Int_t rc=0;
19   
20   
21   TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(headersFile.Data());
22   if(file){
23     cout<<"headerFile already open \n";
24   }
25   else {
26     if(!file)file=TFile::Open(headersFile.Data());
27   }
28   
29   // Get AliRun object from file
30   if (!gAlice) {
31     gAlice = (AliRun*)file->Get("gAlice");
32     if (gAlice) printf("AliRun object found on file\n");
33   }
34   
35   
36   if (iEvNum == 0) iEvNum = (Int_t) gAlice->TreeE()->GetEntries();
37   
38   AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
39
40   if (!tof) {
41     cout << "<AliTOFanalyzeSDigitsV2> No TOF detector found" << endl;
42     rc = 2;
43     return rc;
44   }
45
46   // adc and tdc
47   TH1F *htdc  = new TH1F("htdc","TDC [bin]",500,0.,15000.);
48   TH1F *hadc   = new TH1F("hadc","ADC [bin]",100,0., 3000.);
49
50   // number of slot for each sdigits
51   TH1F *hndigits  = new TH1F("hndigits","fNDigits",5,0.,5.);
52
53   // TOF sdigit volumes
54   TH1F *hsector  = new TH1F("hsector","Sector",20,0.,20.);
55   TH1F *hplate   = new TH1F("hplate","Plate ", 6,0., 6.);
56   TH1F *hstrip   = new TH1F("hstrip","Strip ",25,0.,25.);
57   TH1F *hpadz    = new TH1F("hpadz","Pad along z ",3,0.,3.);
58   TH1F *hpadx    = new TH1F("hpadx","Pad along x",50,0.,50.);
59   // ADC-TDC correlation
60   TH2F *h2tdcVSadc = new TH2F("h2tdcVSadc","TDC [bin] VS ADC [bin]",100,0.,3000.,500,0.,15000.);
61
62   cout << "First " << ndump << " SDigits found in TOF TreeS branch have:" << endl;
63
64
65   for (Int_t ievent = 0; ievent < iEvNum; ievent++) {
66     
67     char tname[100]; sprintf(tname,"TreeS%d",ievent);
68     
69     TTree *sdigitstree=(TTree*)file->Get(tname);                   
70     
71     TClonesArray * fSDigits= new TClonesArray("AliTOFSDigit",  1000); 
72     sdigitstree->GetBranch("TOF")->SetAddress(&fSDigits);           
73     
74     Int_t nEntries = sdigitstree->GetEntries();                                
75     //cout << nEntries << endl;
76     
77     // Loop through all entries in the tree
78     Int_t nbytes;
79     for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
80       
81       // Import the tree
82       nbytes += sdigitstree->GetEvent(iEntry);
83       
84       // Get the number of sdigits
85       Int_t ndig = fSDigits->GetEntriesFast();
86       cout << "------------------<AliTOFanalyzeSDigitsV2>------------------" << endl;
87       cout << "Found " << ndig << " TOF SDigits for event " << ievent << endl;
88       cout << "------------------------------------------------------------" << endl;
89
90       // start loop on sdigits
91       for (Int_t k = 0; k < ndig; k++) {
92         
93         // Get the information for this digit
94         AliTOFSDigit *tofsdigit = (AliTOFSDigit *) fSDigits->UncheckedAt(k);
95         Float_t firstTDC=tofsdigit->GetTdc(0);
96         Float_t firstADC=tofsdigit->GetAdc(0);
97         Int_t nslot=tofsdigit->GetNDigits();
98         hndigits->Fill(nslot);
99         htdc->Fill(firstTDC);
100         hadc->Fill(firstADC);
101         // TOF sdigit volumes
102         Int_t sector    = tofsdigit->GetSector(); // range [1-18]
103         Int_t plate     = tofsdigit->GetPlate();  // range [1- 5]
104         Int_t strip     = tofsdigit->GetStrip();  // range [1-20]
105         Int_t padz      = tofsdigit->GetPadz();   // range [1- 2]
106         Int_t padx      = tofsdigit->GetPadx();   // range [1-48]
107         // it is QA, then I perform QA!
108         Bool_t isSDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
109         
110         if (isSDigitBad) {
111           cout << "<AliTOFanalyzeSDigitsV2>  strange sdigit found" << endl;
112           rc = 3;
113           return rc;
114         }
115         
116         if(k<ndump){
117           cout << k << "-th | " << "Sector " << sector << " | Plate " << plate << " | Strip " << strip << " | PadZ " << padz << " | PadX " << padx << endl;
118           cout << k << "-th | ADC " << firstADC << " [bin] | TDC " << firstTDC << " [bin]" << endl;
119           cout << "----------------------------------------------------"<< endl;
120         }
121         
122         // filling sdigit volume histos
123         hsector->Fill(sector);
124         hplate->Fill(plate);
125         hstrip->Fill(strip);
126         hpadx->Fill(padx);
127         hpadz->Fill(padz);
128         h2tdcVSadc->Fill(firstADC,firstTDC);
129         
130         //cout << "firstTDC " << firstTDC << " firstADC " << firstADC << endl;
131       } // end loop on sdigits
132       
133     } // end loop on entries
134     
135   } // end loop on events 
136   
137   TFile *fout = new TFile("TOF_sdigitsQA.root","RECREATE");
138   htdc->Write();
139   hadc->Write();
140   hndigits->Write();
141   h2tdcVSadc->Write();
142   hsector->Write();
143   hplate->Write();
144   hstrip->Write();
145   hpadz->Write();
146   hpadx->Write();
147   fout->Close(); 
148
149
150   return rc;
151   
152
153 }