]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliHMPIDQATask.cxx
Remove also QA files before re-running the reconstruction
[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   fOutputContainer(0),
46   fhHMPIDCkovP(0),
47   fhHMPIDMipXY(0),
48   fhHMPIDDifXY(0),
49   fhHMPIDSigP(0)
50 {
51   // Constructor.
52   // Input slot #0 works with an Ntuple
53   DefineInput(0, TChain::Class());
54   // Output slot #0 writes into a TH1 container
55   DefineOutput(0,  TObjArray::Class()) ; 
56
57   Int_t i ; 
58   for(i = 0 ; i < 5 ; i++) 
59     fhHMPIDProb[i]=0;
60 }
61
62 //______________________________________________________________________________
63 AliHMPIDQATask::~AliHMPIDQATask()
64 {
65   // dtor
66   fOutputContainer->Clear() ; 
67   delete fOutputContainer ; 
68   
69   delete fhHMPIDCkovP ;  
70   delete fhHMPIDMipXY ;  
71   delete fhHMPIDDifXY ;  
72   delete fhHMPIDSigP ;   
73   delete [] fhHMPIDProb ;
74 }
75
76 //______________________________________________________________________________
77 void AliHMPIDQATask::ConnectInputData(const Option_t*)
78 {
79   // Initialisation of branch container and histograms 
80     
81   AliInfo(Form("*** Initialization of %s", GetName())) ; 
82   
83   // Get input data
84   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
85   if (!fChain) {
86     AliError(Form("Input 0 for %s not found\n", GetName()));
87     return ;
88   }
89   
90   // One should first check if the branch address was taken by some other task
91   char ** address = (char **)GetBranchAddress(0, "ESD");
92   if (address) {
93     fESD = (AliESD*)(*address);
94   } else {
95     fESD = new AliESD();
96     SetBranchAddress(0, "ESD", &fESD);
97   }
98 }
99
100 //________________________________________________________________________
101 void AliHMPIDQATask::CreateOutputObjects()
102 {  
103   // create histograms 
104
105   OpenFile(0) ; 
106
107   fhHMPIDCkovP    = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]", 150,   0,  7  ,100, -3, 1); 
108   fhHMPIDSigP     = new TH2F("SigP"  ,"#sigma_{#theta_c}"          , 150,   0,  7  ,100, 0, 1e20);
109   fhHMPIDMipXY    = new TH2F("MipXY" ,"mip position"               , 260,   0,130  ,252,0,126); 
110   fhHMPIDDifXY    = new TH2F("DifXY" ,"diff"                       , 260, -10, 10  ,252,-10,10); 
111   
112   fhHMPIDProb[0] = new TH1F("PidE" ,"PID: e yellow #mu magenta"  ,100,0,1); 
113   fhHMPIDProb[0]->SetLineColor(kYellow);
114   fhHMPIDProb[1] = new TH1F("PidMu","pid of #mu"                 ,100,0,1); 
115   fhHMPIDProb[1]->SetLineColor(kMagenta);
116   fhHMPIDProb[2] = new TH1F("PidPi","PID: #pi red K green p blue",100,0,1); 
117   fhHMPIDProb[2]->SetLineColor(kRed);
118   fhHMPIDProb[3] = new TH1F("PidK" ,"pid of K"                   ,100,0,1); 
119   fhHMPIDProb[3]->SetLineColor(kGreen);
120   fhHMPIDProb[4] = new TH1F("PidP" ,"pid of p"                   ,100,0,1); 
121   fhHMPIDProb[4]->SetLineColor(kBlue);
122  
123
124   
125   // create output container
126   
127   fOutputContainer = new TObjArray(9) ; 
128   fOutputContainer->SetName(GetName()) ; 
129
130   fOutputContainer->AddAt(fhHMPIDCkovP,      0) ; 
131   fOutputContainer->AddAt(fhHMPIDSigP,       1) ; 
132   fOutputContainer->AddAt(fhHMPIDMipXY,      2) ; 
133   fOutputContainer->AddAt(fhHMPIDDifXY,      3) ; 
134   fOutputContainer->AddAt(fhHMPIDProb[0],    4) ; 
135   fOutputContainer->AddAt(fhHMPIDProb[1],    5) ; 
136   fOutputContainer->AddAt(fhHMPIDProb[2],    6) ; 
137   fOutputContainer->AddAt(fhHMPIDProb[3],    7) ; 
138   fOutputContainer->AddAt(fhHMPIDProb[4],    8) ; 
139 }
140
141 //______________________________________________________________________________
142 void AliHMPIDQATask::Exec(Option_t *) 
143 {
144   // Processing of one event
145     
146   Long64_t entry = fChain->GetReadEntry() ;
147   
148   if (!fESD) {
149     AliError("fESD is not connected to the input!") ; 
150     return ; 
151   }
152   
153   if ( !((entry-1)%100) ) 
154     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
155   
156   // ************************  HMPID *************************************
157   Int_t iTrk ; 
158   for(iTrk = 0 ; iTrk < fESD->GetNumberOfTracks() ; iTrk++){
159     AliESDtrack *pTrk = fESD->GetTrack(iTrk) ;
160
161     fhHMPIDCkovP->Fill( pTrk->GetP(), pTrk->GetHMPIDsignal() ) ; 
162     fhHMPIDSigP ->Fill( pTrk->GetP(), TMath::Sqrt(pTrk->GetHMPIDchi2()) ) ;
163      
164 //     Float_t xm,ym; Int_t q,np;  pTrk->GetHMPIDmip(xm,ym,q,np);  fMipXY->Fill(xm,ym); //mip info
165 //     Float_t xd,yd,th,ph;        pTrk->GetHMPIDtrk(xd,yd,th,ph); fDifXY->Fill(xd,yd); //track info 
166      
167     Double_t pid[5] ;  
168     pTrk->GetHMPIDpid(pid) ; 
169     Int_t i ; 
170     for(i = 0 ; i < 5 ; i++) 
171       fhHMPIDProb[i]->Fill(pid[i]) ;
172   }//tracks loop 
173        
174   PostData(0, fOutputContainer);
175 }
176
177 //______________________________________________________________________________
178 void AliHMPIDQATask::Terminate(Option_t *)
179 {
180   // Processing when the event loop is ended
181   fOutputContainer = (TObjArray*)GetOutputData(0);
182   fhHMPIDCkovP   = (TH2F*)fOutputContainer->At(0);
183   fhHMPIDSigP    = (TH2F*)fOutputContainer->At(1);
184   fhHMPIDMipXY   = (TH2F*)fOutputContainer->At(2);
185   fhHMPIDDifXY   = (TH2F*)fOutputContainer->At(3);
186   fhHMPIDProb[0] = (TH1F*)fOutputContainer->At(4);
187   fhHMPIDProb[1] = (TH1F*)fOutputContainer->At(5);
188   fhHMPIDProb[2] = (TH1F*)fOutputContainer->At(6);
189   fhHMPIDProb[3] = (TH1F*)fOutputContainer->At(7);
190   fhHMPIDProb[4] = (TH1F*)fOutputContainer->At(8);
191
192   Bool_t problem = kFALSE ; 
193   AliInfo(Form(" *** %s Report:", GetName())) ; 
194
195   Float_t n = 1.292 ; //mean freon ref idx 
196   TF1 * hHMPIDpPi = new TF1("RiPiTheo", "acos(sqrt(x*x+[0]*[0])/(x*[1]))", 1.2, 7) ; 
197   hHMPIDpPi->SetLineWidth(1) ; 
198   hHMPIDpPi->SetParameter(1,n) ; 
199
200   AliPID ppp ;                 
201   hHMPIDpPi->SetLineColor(kRed);   
202   hHMPIDpPi->SetParameter(0,AliPID::ParticleMass(AliPID::kPion));    //mass
203
204   TF1 * hHMPIDK = static_cast<TF1*>(hHMPIDpPi->Clone()) ; 
205   hHMPIDK ->SetLineColor(kGreen) ; 
206   hHMPIDK ->SetParameter(0, AliPID::ParticleMass(AliPID::kKaon)) ; 
207
208   TF1 * hHMPIDP=static_cast<TF1*>(hHMPIDpPi->Clone()) ; 
209   hHMPIDP ->SetLineColor(kBlue) ;  
210   hHMPIDP ->SetParameter(0,AliPID::ParticleMass(AliPID::kProton)) ; 
211
212   TCanvas * cHMPID = new TCanvas("cHMPID","HMPID ESD Test") ;
213   cHMPID->SetFillColor(10) ; 
214   cHMPID->SetHighLightColor(10) ; 
215   cHMPID->Divide(3,2) ;
216
217   cHMPID->cd(1); 
218   fhHMPIDCkovP->Draw() ; 
219   hHMPIDpPi->Draw("same") ; 
220   hHMPIDK->Draw("same") ; 
221   hHMPIDP->Draw("same") ;   
222
223   cHMPID->cd(2) ; 
224   fhHMPIDMipXY->Draw() ;   
225
226   cHMPID->cd(3) ; 
227   fhHMPIDProb[0]->Draw() ; 
228   fhHMPIDProb[1]->Draw("same") ;
229   
230   cHMPID->cd(4) ; 
231   fhHMPIDSigP ->Draw() ;                                                          
232
233   cHMPID->cd(5) ; 
234   fhHMPIDDifXY->Draw() ;   
235
236   cHMPID->cd(6) ; 
237   fhHMPIDProb[2]->Draw() ; 
238   fhHMPIDProb[3]->Draw("same") ; 
239   fhHMPIDProb[4]->Draw("same") ; 
240
241   cHMPID->Print("HMPID.eps");
242   
243   char line[1024] ; 
244   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
245   gROOT->ProcessLine(line);
246   sprintf(line, ".!rm -fR *.eps"); 
247   gROOT->ProcessLine(line);
248   
249   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
250
251   TString report ; 
252   if(problem)
253     report="Problems found, please check!!!";  
254   else 
255     report="OK";
256
257   AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ; 
258 }