]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFanalyzeSDigits.C
added method to print cpass0/x calibration status
[u/mrichter/AliRoot.git] / TOF / AliTOFanalyzeSDigits.C
1 Int_t AliTOFanalyzeSDigits(TString headersFile, 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   // Dynamically link some shared libs
9   if (gClassTable->GetID("AliRun") < 0) {
10     gROOT->LoadMacro("loadlibs.C");
11     loadlibs();
12   } else {
13     delete gAlice;
14     gAlice = 0;
15   }
16   
17   Int_t rc=0;
18
19   
20   TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(headersFile.Data());
21   if(file){
22     cout<<"headerFile already open \n";
23   }
24   else {
25     if(!file)file=TFile::Open(headersFile.Data());
26   }
27   
28   // Get AliRun object from file
29   if (!gAlice) {
30     gAlice = (AliRun*)file->Get("gAlice");
31     if (gAlice) printf("AliRun object found on file\n");
32   }
33   
34
35   if (iEvNum == 0) iEvNum = (Int_t) gAlice->TreeE()->GetEntries();
36
37
38   AliTOFSDigit *tofsdigit;
39
40   AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
41
42   if (!tof) {
43     cout << "<AliTOFanalyzeSDigits> No TOF detector found" << endl;
44     rc = 2;
45     return rc;
46   }
47
48   // adc and tdc
49   TH1F *htdc  = new TH1F("htdc","TDC [bin]",5000,0.,150000.);
50   TH1F *hadc   = new TH1F("hadc","ADC [bin]",100,0., 3000.);
51
52   // TOF sdigit volumes
53   TH1F *hsector  = new TH1F("hsector","Sector",20,0.,20.);
54   TH1F *hplate   = new TH1F("hplate","Plate ", 6,0., 6.);
55   TH1F *hstrip   = new TH1F("hstrip","Strip ",25,0.,25.);
56   TH1F *hpadz    = new TH1F("hpadz","Pad along z ",3,0.,3.);
57   TH1F *hpadx    = new TH1F("hpadx","Pad along x",50,0.,50.);
58   // ADC-TDC correlation
59   TH2F *h2tdcVSadc = new TH2F("h2tdcVSadc","TDC [bin] VS ADC [bin]",500,0.,150000.,100,0.,3000.);
60
61   cout << "First " << ndump << " SDigits found in TOF TreeS branch have:" << endl;
62
63   for (Int_t ievent = 0; ievent < iEvNum; ievent++) {
64
65     gAlice->GetEvent(ievent) ;
66     if(gAlice->TreeS()==0) {
67       cout << "<AliTOFanalyzeSDigits> No  TreeS found" << endl;
68       rc = 4;
69       return rc;
70     }
71     
72     
73     
74     Int_t ndig, k;
75     gAlice->ResetDigits();
76     gAlice->TreeS()->GetEvent(ievent);
77     TClonesArray * TOFdigits   = tof->SDigits();
78     
79     ndig=TOFdigits->GetEntries();
80     
81     cout << "<AliTOFanalyzeSDigits> found " << ndig
82          << " TOF sdigits for event " << ievent << endl;
83     
84     for (k=0; k<ndig; k++) {
85       tofsdigit= (AliTOFSDigit*) TOFdigits->UncheckedAt(k);
86       Float_t firstTDC=tofsdigit->GetTdc(0);
87       Float_t firstADC=tofsdigit->GetAdc(0);
88       htdc->Fill(firstTDC);
89       hadc->Fill(firstADC);
90       // TOF sdigit volumes
91       Int_t sector    = tofsdigit->GetSector(); // range [1-18]
92       Int_t plate     = tofsdigit->GetPlate();  // range [1- 5]
93       Int_t strip     = tofsdigit->GetStrip();  // range [1-20]
94       Int_t padz      = tofsdigit->GetPadz();   // range [1- 2]
95       Int_t padx      = tofsdigit->GetPadx();   // range [1-48]
96       // it is QA, then I perform QA!
97       Bool_t isSDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
98
99       if (isSDigitBad) {
100         cout << "<AliTOFanalyzeSDigits>  strange sdigit found" << endl;
101         rc = 3;
102         return rc;
103       }
104       
105       if(k<ndump){
106         cout << k << "-th | " << "Sector " << sector << " | Plate " << plate << " | Strip " << strip << " | PadZ " << padz << " | PadX " << padx << endl;
107         cout << k << "-th | ADC " << firstADC << " [bin] | TDC " << firstTDC << " [bin]" << endl;
108         cout << "----------------------------------------------------"<< endl;
109       }
110
111       // filling sdigit volume histos
112       hsector->Fill(sector);
113       hplate->Fill(plate);
114       hstrip->Fill(strip);
115       hpadx->Fill(padx);
116       hpadz->Fill(padz);
117       h2tdcVSadc->Fill(firstTDC,firstADC);
118
119       //cout << "firstTDC " << firstTDC << " firstADC " << firstADC << endl;
120     }
121   
122   } // end loop on events
123
124   TFile *fout = new TFile("TOF_sdigitsQA.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