]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFanalyzeSDigitsV2.C
updating of the TOF volume numbering (ALICE-INT-2003-038)
[u/mrichter/AliRoot.git] / TOF / AliTOFanalyzeSDigitsV2.C
1 Int_t AliTOFanalyzeSDigitsV2(Int_t ndump=15, Int_t iEvNum=0)
2 {
3   //
4   // Analyzes the TOF sdigits and fills QA-histograms 
5   // report problems to pierella@bo.infn.it
6   // iEvNum=0 means all events in the file
7   //
8   // Updated to the new I/O by: A. De Caro, C. Zampolli
9
10   Int_t rc = 0;
11
12   // adc and tdc
13   TH1F *htdc     = new TH1F("htdc","TDC [bin]",5000,0.,150000.);
14   TH1F *hadc     = new TH1F("hadc","ADC [bin]",100,0., 3000.);
15   // TOF sdigit volumes
16   TH1F *hsector  = new TH1F("hsector","Sector",18,0.,18.);
17   TH1F *hplate   = new TH1F("hplate","Plate ", 5,0., 5.);
18   TH1F *hstrip   = new TH1F("hstrip","Strip ",20,0.,20.);
19   TH1F *hpadz    = new TH1F("hpadz","Pad along z ",2,0.,2.);
20   TH1F *hpadx    = new TH1F("hpadx","Pad along x",48,0.,48.);
21   // ADC-TDC correlation
22   TH2F *h2tdcVSadc = new TH2F("h2tdcVSadc","TDC [bin] VS ADC [bin]",500,0.,150000.,100,0.,3000.);
23   
24   // Dynamically link some shared libs
25   if (gClassTable->GetID("AliRun") < 0) {
26     gROOT->LoadMacro("loadlibs.C");
27     loadlibs();
28   }
29
30   if (gAlice) 
31     {
32       delete gAlice->GetRunLoader();
33       delete gAlice;
34       gAlice = 0x0;
35    }
36   
37   AliRunLoader *rl = AliRunLoader::Open("galice.root",AliConfig::fgkDefaultEventFolderName,"read");
38   if (!rl)
39     {
40       cerr<<"Can't load RunLoader from file"<<"!\n";
41       rc = 1;
42       return rc;
43     }
44
45   rl->LoadgAlice();
46   gAlice = rl->GetAliRun();
47
48   if (!gAlice)
49     {
50       cerr << "<AliTOFanalyzeSDigits> AliRun object not found on file\n ";
51       rc = 2;
52       return rc;
53     }
54
55   AliLoader *tofl = rl->GetLoader("TOFLoader");
56   AliTOF *tof = (AliTOF *) rl->GetAliRun()->GetDetector("TOF");
57
58   if (tof==0x0 || tofl==0x0)
59     {
60       cerr << "<AliTOFanalyzeSDigits> no TOF detector found" << endl;
61       rc = 3;
62       return rc;
63     }
64   
65   cout << "First " << ndump << " SDigits found in TOF TreeS branch have:" << endl;
66
67   if (iEvNum == 0)
68     {
69       rl->LoadHeader();
70       TTree *TE = rl->TreeE();
71       iEvNum = (Int_t)TE->GetEntries();
72     }
73
74   AliTOFSDigit *tofsdigit;
75   
76   for (Int_t ievent = 0; ievent < iEvNum; ievent++) {
77     printf ("Processing event %d \n", ievent);
78     rl->GetEvent(ievent);
79     
80     // Get the pointer SDigit tree
81     tofl->LoadSDigits();
82     TTree *TS=tofl->TreeS();
83     tof->SetTreeAddress();
84     
85     if(!TS)
86       {
87         cout << "<AliTOFanalyzeSDigits> No TreeS found" << endl;
88         rc = 4;
89         return rc;
90       }
91     
92     TClonesArray * TOFsdigits = new TClonesArray("AliTOFSDigit",1000);
93     TOFsdigits = tof->SDigits();
94     TOFsdigits = TS->GetBranch("TOF")->SetAddress(&TOFsdigits); 
95
96     Int_t nEntries = TS->GetEntries(); 
97
98     for (Int_t iEntry = 0; iEntry < nEntries; iEntry ++) 
99       {
100         tof->ResetDigits();
101         TS->GetEvent(iEntry);
102         Int_t ndig = TOFsdigits->GetEntriesFast();
103         cout << "<AliTOFanalyzeSDigits> found " << ndig
104              << " TOF sdigits for event " << ievent << endl;
105         
106         for (Int_t k=0; k<ndig; k++) {
107           tofsdigit = (AliTOFSDigit*) TOFsdigits->UncheckedAt(k);
108           Float_t firstTDC = tofsdigit->GetTdc(0);
109           Float_t firstADC = tofsdigit->GetAdc(0);
110           htdc->Fill(firstTDC);
111           hadc->Fill(firstADC);
112           // TOF sdigit volumes
113           Int_t sector = tofsdigit->GetSector(); // range [0-17]
114           Int_t plate  = tofsdigit->GetPlate();  // range [0- 4]
115           Int_t strip  = tofsdigit->GetStrip();  // range [0-19]
116           Int_t padz   = tofsdigit->GetPadz();   // range [0- 1]
117           Int_t padx   = tofsdigit->GetPadx();   // range [0-47]
118           // it is QA, then I perform QA!
119           Bool_t isSDigitBad = (sector<0 || sector>17 || plate<0 || plate >4 || padz<0 || padz>1 || padx<0 || padx>47);
120           
121
122           if (isSDigitBad) {
123             cout << "<AliTOFanalyzeSDigits>  strange sdigit found" << endl;
124             rc = 4;
125             return rc;
126           }
127           
128           if(k<ndump){
129             cout << k << "-th | " << "Sector " << sector << " | Plate " << plate << " | Strip " << strip << " | PadZ " << padz << " | PadX " << padx << endl;
130             cout << k << "-th | ADC " << firstADC << " [bin] | TDC " << firstTDC << " [bin]" << endl;
131             cout << "----------------------------------------------------"<< endl;
132           }
133           
134           // filling sdigit volume histos
135           hsector->Fill(sector);
136           hplate->Fill(plate);
137           hstrip->Fill(strip);
138           hpadx->Fill(padx);
139           hpadz->Fill(padz);
140           h2tdcVSadc->Fill(firstTDC,firstADC);
141           
142         }
143       }
144     
145     tofl->UnloadSDigits();
146   
147   } // end loop on events
148   
149   rl->UnloadHeader();
150   rl->UnloadgAlice();
151
152   TFile *fout = new TFile("TOF_sdigitsQA.root","RECREATE");
153   htdc->Write();
154   hadc->Write();
155   h2tdcVSadc->Write();
156   hsector->Write();
157   hplate->Write();
158   hstrip->Write();
159   hpadz->Write();
160   hpadx->Write();
161   fout->Close();
162
163   delete htdc;
164   delete hadc;
165   delete h2tdcVSadc;
166   delete hsector;
167   delete hplate;
168   delete hstrip;   
169   delete hpadz;
170   delete hpadx;
171
172   if (gAlice) 
173     {
174       delete gAlice->GetRunLoader();
175       delete gAlice;
176       gAlice = 0x0;
177     }
178
179   return rc;
180   
181 }
182