]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliHMPIDQATask.cxx
Pass event number as argument of AliVEventHandler::BeginEvent
[u/mrichter/AliRoot.git] / ESDCheck / AliHMPIDQATask.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // An analysis task to check the HMPID data in simulated data
20 //
21 //*-- Annalisa Mastroserio
22 //////////////////////////////////////////////////////////////////////////////
23
24 #include <TCanvas.h>
25 #include <TChain.h>
26 #include <TF1.h>
27 #include <TFile.h> 
28 #include <TH1F.h>
29 #include <TH2F.h>
30 #include <TLegend.h> 
31 #include <TROOT.h>
32 #include <TVector3.h> 
33 #include <TString.h> 
34
35 #include "AliHMPIDQATask.h" 
36 #include "AliESD.h" 
37 #include "AliLog.h"
38 #include "AliPID.h"
39
40 //______________________________________________________________________________
41 AliHMPIDQATask::AliHMPIDQATask(const char *name) : 
42   AliAnalysisTask(name,""),  
43   fChain(0),
44   fESD(0), 
45   fhHMPIDCkovP(0),
46   fhHMPIDMipXY(0),
47   fhHMPIDDifXY(0),
48   fhHMPIDSigP(0)
49 {
50   // Constructor.
51   // Input slot #0 works with an Ntuple
52   DefineInput(0, TChain::Class());
53   // Output slot #0 writes into a TH1 container
54   DefineOutput(0,  TObjArray::Class()) ; 
55
56   Int_t i ; 
57   for(i = 0 ; i < 5 ; i++) 
58     fhHMPIDProb[i]=0;
59 }
60
61 //______________________________________________________________________________
62 AliHMPIDQATask::~AliHMPIDQATask()
63 {
64   // dtor
65   fOutputContainer->Clear() ; 
66   delete fOutputContainer ; 
67   
68   delete fhHMPIDCkovP ;  
69   delete fhHMPIDMipXY ;  
70   delete fhHMPIDDifXY ;  
71   delete fhHMPIDSigP ;   
72   delete [] fhHMPIDProb ;
73 }
74
75 //______________________________________________________________________________
76 void AliHMPIDQATask::ConnectInputData(const Option_t*)
77 {
78   // Initialisation of branch container and histograms 
79     
80   AliInfo(Form("*** Initialization of %s", GetName())) ; 
81   
82   // Get input data
83   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
84   if (!fChain) {
85     AliError(Form("Input 0 for %s not found\n", GetName()));
86     return ;
87   }
88   
89   // One should first check if the branch address was taken by some other task
90   char ** address = (char **)GetBranchAddress(0, "ESD");
91   if (address) {
92     fESD = (AliESD*)(*address);
93   } else {
94     fESD = new AliESD();
95     SetBranchAddress(0, "ESD", &fESD);
96   }
97 }
98
99 //________________________________________________________________________
100 void AliHMPIDQATask::CreateOutputObjects()
101 {  
102   // create histograms 
103
104   OpenFile(0) ; 
105
106   fhHMPIDCkovP    = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]", 150,   0,  7  ,100, -3, 1); 
107   fhHMPIDSigP     = new TH2F("SigP"  ,"#sigma_{#theta_c}"          , 150,   0,  7  ,100, 0, 1e20);
108   fhHMPIDMipXY    = new TH2F("MipXY" ,"mip position"               , 260,   0,130  ,252,0,126); 
109   fhHMPIDDifXY    = new TH2F("DifXY" ,"diff"                       , 260, -10, 10  ,252,-10,10); 
110   
111   fhHMPIDProb[0] = new TH1F("PidE" ,"PID: e yellow #mu magenta"  ,100,0,1); 
112   fhHMPIDProb[0]->SetLineColor(kYellow);
113   fhHMPIDProb[1] = new TH1F("PidMu","pid of #mu"                 ,100,0,1); 
114   fhHMPIDProb[1]->SetLineColor(kMagenta);
115   fhHMPIDProb[2] = new TH1F("PidPi","PID: #pi red K green p blue",100,0,1); 
116   fhHMPIDProb[2]->SetLineColor(kRed);
117   fhHMPIDProb[3] = new TH1F("PidK" ,"pid of K"                   ,100,0,1); 
118   fhHMPIDProb[3]->SetLineColor(kGreen);
119   fhHMPIDProb[4] = new TH1F("PidP" ,"pid of p"                   ,100,0,1); 
120   fhHMPIDProb[4]->SetLineColor(kBlue);
121  
122
123   
124   // create output container
125   
126   fOutputContainer = new TObjArray(9) ; 
127   fOutputContainer->SetName(GetName()) ; 
128
129   fOutputContainer->AddAt(fhHMPIDCkovP,      0) ; 
130   fOutputContainer->AddAt(fhHMPIDSigP,       1) ; 
131   fOutputContainer->AddAt(fhHMPIDMipXY,      2) ; 
132   fOutputContainer->AddAt(fhHMPIDDifXY,      3) ; 
133   fOutputContainer->AddAt(fhHMPIDProb[0],    4) ; 
134   fOutputContainer->AddAt(fhHMPIDProb[1],    5) ; 
135   fOutputContainer->AddAt(fhHMPIDProb[2],    6) ; 
136   fOutputContainer->AddAt(fhHMPIDProb[3],    7) ; 
137   fOutputContainer->AddAt(fhHMPIDProb[4],    8) ; 
138 }
139
140 //______________________________________________________________________________
141 void AliHMPIDQATask::Exec(Option_t *) 
142 {
143   // Processing of one event
144     
145   Long64_t entry = fChain->GetReadEntry() ;
146   
147   if (!fESD) {
148     AliError("fESD is not connected to the input!") ; 
149     return ; 
150   }
151   
152   if ( !((entry-1)%100) ) 
153     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
154   
155   // ************************  HMPID *************************************
156   Int_t iTrk ; 
157   for(iTrk = 0 ; iTrk < fESD->GetNumberOfTracks() ; iTrk++){
158     AliESDtrack *pTrk = fESD->GetTrack(iTrk) ;
159
160     fhHMPIDCkovP->Fill( pTrk->GetP(), pTrk->GetHMPIDsignal() ) ; 
161     fhHMPIDSigP ->Fill( pTrk->GetP(), TMath::Sqrt(pTrk->GetHMPIDchi2()) ) ;
162      
163 //     Float_t xm,ym; Int_t q,np;  pTrk->GetHMPIDmip(xm,ym,q,np);  fMipXY->Fill(xm,ym); //mip info
164 //     Float_t xd,yd,th,ph;        pTrk->GetHMPIDtrk(xd,yd,th,ph); fDifXY->Fill(xd,yd); //track info 
165      
166     Double_t pid[5] ;  
167     pTrk->GetHMPIDpid(pid) ; 
168     Int_t i ; 
169     for(i = 0 ; i < 5 ; i++) 
170       fhHMPIDProb[i]->Fill(pid[i]) ;
171   }//tracks loop 
172        
173   PostData(0, fOutputContainer);
174 }
175
176 //______________________________________________________________________________
177 void AliHMPIDQATask::Terminate(Option_t *)
178 {
179   // Processing when the event loop is ended
180   fOutputContainer = (TObjArray*)GetOutputData(0);
181   fhHMPIDCkovP   = (TH2F*)fOutputContainer->At(0);
182   fhHMPIDSigP    = (TH2F*)fOutputContainer->At(1);
183   fhHMPIDMipXY   = (TH2F*)fOutputContainer->At(2);
184   fhHMPIDDifXY   = (TH2F*)fOutputContainer->At(3);
185   fhHMPIDProb[0] = (TH1F*)fOutputContainer->At(4);
186   fhHMPIDProb[1] = (TH1F*)fOutputContainer->At(5);
187   fhHMPIDProb[2] = (TH1F*)fOutputContainer->At(6);
188   fhHMPIDProb[3] = (TH1F*)fOutputContainer->At(7);
189   fhHMPIDProb[4] = (TH1F*)fOutputContainer->At(8);
190
191   Bool_t problem = kFALSE ; 
192   AliInfo(Form(" *** %s Report:", GetName())) ; 
193
194   Float_t n = 1.292 ; //mean freon ref idx 
195   TF1 * hHMPIDpPi = new TF1("RiPiTheo", "acos(sqrt(x*x+[0]*[0])/(x*[1]))", 1.2, 7) ; 
196   hHMPIDpPi->SetLineWidth(1) ; 
197   hHMPIDpPi->SetParameter(1,n) ; 
198
199   AliPID ppp ;                 
200   hHMPIDpPi->SetLineColor(kRed);   
201   hHMPIDpPi->SetParameter(0,AliPID::ParticleMass(AliPID::kPion));    //mass
202
203   TF1 * hHMPIDK = static_cast<TF1*>(hHMPIDpPi->Clone()) ; 
204   hHMPIDK ->SetLineColor(kGreen) ; 
205   hHMPIDK ->SetParameter(0, AliPID::ParticleMass(AliPID::kKaon)) ; 
206
207   TF1 * hHMPIDP=static_cast<TF1*>(hHMPIDpPi->Clone()) ; 
208   hHMPIDP ->SetLineColor(kBlue) ;  
209   hHMPIDP ->SetParameter(0,AliPID::ParticleMass(AliPID::kProton)) ; 
210
211   TCanvas * cHMPID = new TCanvas("cHMPID","HMPID ESD Test") ;
212   cHMPID->SetFillColor(10) ; 
213   cHMPID->SetHighLightColor(10) ; 
214   cHMPID->Divide(3,2) ;
215
216   cHMPID->cd(1); 
217   fhHMPIDCkovP->Draw() ; 
218   hHMPIDpPi->Draw("same") ; 
219   hHMPIDK->Draw("same") ; 
220   hHMPIDP->Draw("same") ;   
221
222   cHMPID->cd(2) ; 
223   fhHMPIDMipXY->Draw() ;   
224
225   cHMPID->cd(3) ; 
226   fhHMPIDProb[0]->Draw() ; 
227   fhHMPIDProb[1]->Draw("same") ;
228   
229   cHMPID->cd(4) ; 
230   fhHMPIDSigP ->Draw() ;                                                          
231
232   cHMPID->cd(5) ; 
233   fhHMPIDDifXY->Draw() ;   
234
235   cHMPID->cd(6) ; 
236   fhHMPIDProb[2]->Draw() ; 
237   fhHMPIDProb[3]->Draw("same") ; 
238   fhHMPIDProb[4]->Draw("same") ; 
239
240   cHMPID->Print("HMPID.eps");
241   
242   char line[1024] ; 
243   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
244   gROOT->ProcessLine(line);
245   sprintf(line, ".!rm -fR *.eps"); 
246   gROOT->ProcessLine(line);
247   
248   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
249
250   TString report ; 
251   if(problem)
252     report="Problems found, please check!!!";  
253   else 
254     report="OK";
255
256   AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ; 
257 }