]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFanalyzeDigits.C
Added ConnectionMapCIterator typedef
[u/mrichter/AliRoot.git] / TOF / AliTOFanalyzeDigits.C
1 Int_t AliTOFanalyzeDigits(Int_t ndump=15, Int_t iEvNum=0)
2 {
3   //
4   // Analyzes the TOF digits and fills QA-histograms 
5   // report problems to pierella@bo.infn.it
6   // iEvNum=0 means all events in the file
7   // Author: F. Pierella (Bologna University)
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]",500,0.,15000.);
14   TH1F *hadc   = new TH1F("hadc","ADC [bin]",100,0., 3000.);
15   // TOF digit volumes
16   TH1F *hsector  = new TH1F("hsector","Sector",20,0.,20.);
17   TH1F *hplate   = new TH1F("hplate","Plate ", 6,0., 6.);
18   TH1F *hstrip   = new TH1F("hstrip","Strip ",25,0.,25.);
19   TH1F *hpadz    = new TH1F("hpadz","Pad along z ",3,0.,3.);
20   TH1F *hpadx    = new TH1F("hpadx","Pad along x",50,0.,50.);
21   // ADC-TDC correlation
22   TH2F *h2tdcVSadc = new TH2F("h2tdcVSadc","TDC [bin] VS ADC [bin]",100,0.,3000.,500,0.,15000.);
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 << "<AliTOFanalyzeDigits> 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 << "<AliTOFanalyzeDigits> No TOF detector found" << endl;
61       rc = 3;
62       return rc;
63     }
64   
65   cout << "First " << ndump << " Digits found in TOF TreeD 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   AliTOFdigit *tofdigit;  
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 Digit tree
81     tofl->LoadDigits();
82     TTree *TD=tofl->TreeD();
83     tof->SetTreeAddress();
84
85     if(!TD)
86       {
87         cout << "<AliTOFanalyzeDigits> No  TreeD found" << endl;
88         rc = 3;
89         return rc;
90       }
91
92     TClonesArray * TOFdigits = new TClonesArray("AliTOFdigit",1000);
93     TOFdigits = tof->Digits();
94     TOFdigits = TD->GetBranch("TOF")->SetAddress(&TOFdigits); 
95
96     Int_t nEntries = TD->GetEntries(); 
97
98     for (Int_t iEntry = 0; iEntry < nEntries; iEntry ++) 
99       {
100         tof->ResetDigits();
101         TD->GetEvent(iEntry);
102         Int_t ndig = TOFdigits->GetEntriesFast();
103         cout << "<AliTOFanalyzeDigits> found " << ndig
104              << " TOF digits for event " << ievent << endl;
105         
106         for (Int_t k=0; k<ndig; k++) { 
107           tofdigit = (AliTOFdigit*) TOFdigits->UncheckedAt(k);
108           Float_t tdc = tofdigit->GetTdc();
109           Float_t adc = tofdigit->GetAdc();
110           htdc->Fill(tdc);
111           hadc->Fill(adc);
112           // TOF digit volumes
113           Int_t sector = tofdigit->GetSector(); // range [1-18]
114           Int_t plate  = tofdigit->GetPlate();  // range [1- 5]
115           Int_t strip  = tofdigit->GetStrip();  // range [1-20]
116           Int_t padz   = tofdigit->GetPadz();   // range [1- 2]
117           Int_t padx   = tofdigit->GetPadx();   // range [1-48]
118           // it is QA, then I perform QA!
119           Bool_t isDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
120           
121           if (isDigitBad)
122             {
123               cout << "<AliTOFanalyzeDigits>  strange digit 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 " << adc << " [bin] | TDC " << tdc << " [bin]" << endl;
131             cout << "----------------------------------------------------"<< endl;
132           }
133           
134           // filling digit 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(adc,tdc);
141
142         }
143       }
144     
145     tofl->UnloadDigits();
146     rl->UnloadHeader();
147     rl->UnloadgAlice();
148   
149   } // end loop on events
150   
151   TFile *fout = new TFile("TOF_digitsQA.root","RECREATE");
152   htdc->Write();
153   hadc->Write();
154   h2tdcVSadc->Write();
155   hsector->Write();
156   hplate->Write();
157   hstrip->Write();
158   hpadz->Write();
159   hpadx->Write();
160   fout->Close(); 
161
162   delete htdc;
163   delete hadc;
164   delete h2tdcVSadc;
165   delete hsector;
166   delete hplate;
167   delete hstrip;   
168   delete hpadz;
169   delete hpadx;
170
171   if (gAlice) 
172     {
173       delete gAlice->GetRunLoader();
174       delete gAlice;
175       gAlice = 0x0;
176     }
177
178   return rc;
179   
180 }
181