]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ESDCheck/AliHMPIDQATask.cxx
fixed compiler problem
[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 **************************************************************************/
0b28fd57 15
16/* $Id$ */
17
1dfe075f 18//_________________________________________________________________________
19// An analysis task to check the HMPID data in simulated data
20//
21//*-- Annalisa Mastroserio
22//////////////////////////////////////////////////////////////////////////////
23
0b28fd57 24#include <TCanvas.h>
1dfe075f 25#include <TChain.h>
0b28fd57 26#include <TF1.h>
27#include <TFile.h>
1dfe075f 28#include <TH1F.h>
29#include <TH2F.h>
1dfe075f 30#include <TLegend.h>
0b28fd57 31#include <TROOT.h>
1dfe075f 32#include <TVector3.h>
1dfe075f 33
34#include "AliHMPIDQATask.h"
35#include "AliESD.h"
36#include "AliLog.h"
37#include "AliPID.h"
38
39//______________________________________________________________________________
40AliHMPIDQATask::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//______________________________________________________________________________
61AliHMPIDQATask::~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//______________________________________________________________________________
75void AliHMPIDQATask::Init(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 if (!fESD) {
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 if (!fESD)
94 fChain->SetBranchAddress("ESD", &fESD) ;
95 }
96 // The output objects will be written to
97 TDirectory * cdir = gDirectory ;
98 // Open a file for output #0
99 char outputName[1024] ;
100 sprintf(outputName, "%s.root", GetName() ) ;
101 OpenFile(0, outputName , "RECREATE") ;
102 if (cdir)
103 cdir->cd() ;
104
105 // create histograms
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//______________________________________________________________________________
141void 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//______________________________________________________________________________
177void AliHMPIDQATask::Terminate(Option_t *)
178{
179 // Processing when the event loop is ended
180
181 Float_t n = 1.292 ; //mean freon ref idx
182 TF1 * hHMPIDpPi = new TF1("RiPiTheo", "acos(sqrt(x*x+[0]*[0])/(x*[1]))", 1.2, 7) ;
183 hHMPIDpPi->SetLineWidth(1) ;
184 hHMPIDpPi->SetParameter(1,n) ;
185
186 AliPID ppp ;
187 hHMPIDpPi->SetLineColor(kRed);
188 hHMPIDpPi->SetParameter(0,AliPID::ParticleMass(AliPID::kPion)); //mass
189
190 TF1 * hHMPIDK = static_cast<TF1*>(hHMPIDpPi->Clone()) ;
191 hHMPIDK ->SetLineColor(kGreen) ;
192 hHMPIDK ->SetParameter(0, AliPID::ParticleMass(AliPID::kKaon)) ;
193
194 TF1 * hHMPIDP=static_cast<TF1*>(hHMPIDpPi->Clone()) ;
195 hHMPIDP ->SetLineColor(kBlue) ;
196 hHMPIDP ->SetParameter(0,AliPID::ParticleMass(AliPID::kProton)) ;
197
198 TCanvas * cHMPID = new TCanvas("cHMPID","HMPID ESD Test") ;
199 cHMPID->SetFillColor(10) ;
200 cHMPID->SetHighLightColor(10) ;
201 cHMPID->Divide(3,2) ;
202
203 cHMPID->cd(1);
204 fhHMPIDCkovP->Draw() ;
205 hHMPIDpPi->Draw("same") ;
206 hHMPIDK->Draw("same") ;
207 hHMPIDP->Draw("same") ;
208
209 cHMPID->cd(2) ;
210 fhHMPIDMipXY->Draw() ;
211
212 cHMPID->cd(3) ;
213 fhHMPIDProb[0]->Draw() ;
214 fhHMPIDProb[1]->Draw("same") ;
215
216 cHMPID->cd(4) ;
217 fhHMPIDSigP ->Draw() ;
218
219 cHMPID->cd(5) ;
220 fhHMPIDDifXY->Draw() ;
221
222 cHMPID->cd(6) ;
223 fhHMPIDProb[2]->Draw() ;
224 fhHMPIDProb[3]->Draw("same") ;
225 fhHMPIDProb[4]->Draw("same") ;
226
227 cHMPID->Print("HMPID.eps");
228
229 char line[1024] ;
230 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
231 gROOT->ProcessLine(line);
232 sprintf(line, ".!rm -fR *.eps");
233 gROOT->ProcessLine(line);
234
235 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
236}