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