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