]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ESDCheck/AliHMPIDQATask.cxx
Adding includes now needed by ROOT
[u/mrichter/AliRoot.git] / ESDCheck / AliHMPIDQATask.cxx
CommitLineData
1dfe075f 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// An analysis task to check the HMPID data in simulated data
17//
18//*-- Annalisa Mastroserio
19//////////////////////////////////////////////////////////////////////////////
20
21#include <TChain.h>
22#include <TH1F.h>
23#include <TH2F.h>
24#include <TF1.h>
25#include <TCanvas.h>
26#include <TLegend.h>
27#include <TVector3.h>
28#include <TFile.h>
29
30#include "AliHMPIDQATask.h"
31#include "AliESD.h"
32#include "AliLog.h"
33#include "AliPID.h"
34
35//______________________________________________________________________________
36AliHMPIDQATask::AliHMPIDQATask(const char *name) :
37 AliAnalysisTask(name,""),
38 fChain(0),
39 fESD(0),
40 fhHMPIDCkovP(0),
41 fhHMPIDMipXY(0),
42 fhHMPIDDifXY(0),
43 fhHMPIDSigP(0)
44{
45 // Constructor.
46 // Input slot #0 works with an Ntuple
47 DefineInput(0, TChain::Class());
48 // Output slot #0 writes into a TH1 container
49 DefineOutput(0, TObjArray::Class()) ;
50
51 Int_t i ;
52 for(i = 0 ; i < 5 ; i++)
53 fhHMPIDProb[i]=0;
54}
55
56//______________________________________________________________________________
57AliHMPIDQATask::~AliHMPIDQATask()
58{
59 // dtor
60 fOutputContainer->Clear() ;
61 delete fOutputContainer ;
62
63 delete fhHMPIDCkovP ;
64 delete fhHMPIDMipXY ;
65 delete fhHMPIDDifXY ;
66 delete fhHMPIDSigP ;
67 delete [] fhHMPIDProb ;
68}
69
70//______________________________________________________________________________
71void AliHMPIDQATask::Init(const Option_t*)
72{
73 // Initialisation of branch container and histograms
74
75 AliInfo(Form("*** Initialization of %s", GetName())) ;
76
77 // Get input data
78 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
79 if (!fChain) {
80 AliError(Form("Input 0 for %s not found\n", GetName()));
81 return ;
82 }
83
84 if (!fESD) {
85 // One should first check if the branch address was taken by some other task
86 char ** address = (char **)GetBranchAddress(0, "ESD") ;
87 if (address)
88 fESD = (AliESD *)(*address) ;
89 if (!fESD)
90 fChain->SetBranchAddress("ESD", &fESD) ;
91 }
92 // The output objects will be written to
93 TDirectory * cdir = gDirectory ;
94 // Open a file for output #0
95 char outputName[1024] ;
96 sprintf(outputName, "%s.root", GetName() ) ;
97 OpenFile(0, outputName , "RECREATE") ;
98 if (cdir)
99 cdir->cd() ;
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//______________________________________________________________________________
137void 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//______________________________________________________________________________
173void AliHMPIDQATask::Terminate(Option_t *)
174{
175 // Processing when the event loop is ended
176
177 Float_t n = 1.292 ; //mean freon ref idx
178 TF1 * hHMPIDpPi = new TF1("RiPiTheo", "acos(sqrt(x*x+[0]*[0])/(x*[1]))", 1.2, 7) ;
179 hHMPIDpPi->SetLineWidth(1) ;
180 hHMPIDpPi->SetParameter(1,n) ;
181
182 AliPID ppp ;
183 hHMPIDpPi->SetLineColor(kRed);
184 hHMPIDpPi->SetParameter(0,AliPID::ParticleMass(AliPID::kPion)); //mass
185
186 TF1 * hHMPIDK = static_cast<TF1*>(hHMPIDpPi->Clone()) ;
187 hHMPIDK ->SetLineColor(kGreen) ;
188 hHMPIDK ->SetParameter(0, AliPID::ParticleMass(AliPID::kKaon)) ;
189
190 TF1 * hHMPIDP=static_cast<TF1*>(hHMPIDpPi->Clone()) ;
191 hHMPIDP ->SetLineColor(kBlue) ;
192 hHMPIDP ->SetParameter(0,AliPID::ParticleMass(AliPID::kProton)) ;
193
194 TCanvas * cHMPID = new TCanvas("cHMPID","HMPID ESD Test") ;
195 cHMPID->SetFillColor(10) ;
196 cHMPID->SetHighLightColor(10) ;
197 cHMPID->Divide(3,2) ;
198
199 cHMPID->cd(1);
200 fhHMPIDCkovP->Draw() ;
201 hHMPIDpPi->Draw("same") ;
202 hHMPIDK->Draw("same") ;
203 hHMPIDP->Draw("same") ;
204
205 cHMPID->cd(2) ;
206 fhHMPIDMipXY->Draw() ;
207
208 cHMPID->cd(3) ;
209 fhHMPIDProb[0]->Draw() ;
210 fhHMPIDProb[1]->Draw("same") ;
211
212 cHMPID->cd(4) ;
213 fhHMPIDSigP ->Draw() ;
214
215 cHMPID->cd(5) ;
216 fhHMPIDDifXY->Draw() ;
217
218 cHMPID->cd(6) ;
219 fhHMPIDProb[2]->Draw() ;
220 fhHMPIDProb[3]->Draw("same") ;
221 fhHMPIDProb[4]->Draw("same") ;
222
223 cHMPID->Print("HMPID.eps");
224
225 char line[1024] ;
226 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
227 gROOT->ProcessLine(line);
228 sprintf(line, ".!rm -fR *.eps");
229 gROOT->ProcessLine(line);
230
231 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
232}