]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFanalyzeDigits.C
Correcting some trivial warnings on Alpha
[u/mrichter/AliRoot.git] / TOF / AliTOFanalyzeDigits.C
1 Int_t AliTOFanalyzeDigits(TString headersFile, 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
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
39   AliTOFdigit *tofdigit;
40
41   AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
42
43   if (!tof) {
44     cout << "<AliTOFanalyzeDigits> No TOF detector found" << endl;
45     rc = 2;
46     return rc;
47   }
48
49   // adc and tdc
50   TH1F *htdc  = new TH1F("htdc","TDC [bin]",500,0.,15000.);
51   TH1F *hadc   = new TH1F("hadc","ADC [bin]",100,0., 3000.);
52
53   // TOF digit 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 << " Digits found in TOF TreeD branch have:" << endl;
63
64   for (Int_t ievent = 0; ievent < iEvNum; ievent++) {
65
66     gAlice->GetEvent(ievent) ;
67     if(gAlice->TreeD()==0) {
68       cout << "<AliTOFanalyzeDigits> No  TreeD found" << endl;
69       rc = 4;
70       return rc;
71     }
72     
73     
74     
75     Int_t ndig, k;
76     gAlice->ResetDigits();
77     gAlice->TreeD()->GetEvent(ievent);
78     TClonesArray * TOFdigits   = tof->Digits();
79     
80     ndig=TOFdigits->GetEntries();
81     
82     cout << "<AliTOFanalyzeDigits> found " << ndig
83          << " TOF digits for event " << ievent << endl;
84     
85     for (k=0; k<ndig; k++) {
86       tofdigit= (AliTOFdigit*) TOFdigits->UncheckedAt(k);
87       Float_t tdc=tofdigit->GetTdc();
88       Float_t adc=tofdigit->GetAdc();
89       htdc->Fill(tdc);
90       hadc->Fill(adc);
91       // TOF digit volumes
92       Int_t sector    = tofdigit->GetSector(); // range [1-18]
93       Int_t plate     = tofdigit->GetPlate();  // range [1- 5]
94       Int_t strip     = tofdigit->GetStrip();  // range [1-20]
95       Int_t padz      = tofdigit->GetPadz();   // range [1- 2]
96       Int_t padx      = tofdigit->GetPadx();   // range [1-48]
97       // it is QA, then I perform QA!
98       Bool_t isDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
99
100       if (isDigitBad) {
101         cout << "<AliTOFanalyzeDigits>  strange digit found" << endl;
102         rc = 3;
103         return rc;
104       }
105       
106       if(k<ndump){
107         cout << k << "-th | " << "Sector " << sector << " | Plate " << plate << " | Strip " << strip << " | PadZ " << padz << " | PadX " << padx << endl;
108         cout << k << "-th | ADC " << adc << " [bin] | TDC " << tdc << " [bin]" << endl;
109         cout << "----------------------------------------------------"<< endl;
110       }
111
112       // filling digit volume histos
113       hsector->Fill(sector);
114       hplate->Fill(plate);
115       hstrip->Fill(strip);
116       hpadx->Fill(padx);
117       hpadz->Fill(padz);
118       h2tdcVSadc->Fill(adc,tdc);
119
120     }
121   
122   } // end loop on events
123
124   TFile *fout = new TFile("TOF_digitsQA.root","RECREATE");
125   htdc->Write();
126   hadc->Write();
127   h2tdcVSadc->Write();
128   hsector->Write();
129   hplate->Write();
130   hstrip->Write();
131   hpadz->Write();
132   hpadx->Write();
133   fout->Close(); 
134
135
136   return rc;
137   
138 }
139