1 /**************************************************************************
\r
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
\r
4 * Author: The ALICE Off-line Project. *
\r
5 * Contributors are mentioned in the code where appropriate. *
\r
7 * Permission to use, copy, modify and distribute this software and its *
\r
8 * documentation strictly for non-commercial purposes is hereby granted *
\r
9 * without fee, provided that the above copyright notice appears in all *
\r
10 * copies and that both the copyright notice and this permission notice *
\r
11 * appear in the supporting documentation. The authors make no claims *
\r
12 * about the suitability of this software for any purpose. It is *
\r
13 * provided "as is" without express or implied warranty. *
\r
14 **************************************************************************/
\r
16 //-----------------------------------------------------------------
\r
17 // AliAnalysisTaskSpectraAllChAOD class
\r
18 //-----------------------------------------------------------------
\r
22 #include "TLegend.h"
\r
25 #include "THnSparse.h"
\r
26 #include "TCanvas.h"
\r
27 #include "AliAnalysisTask.h"
\r
28 #include "AliAODTrack.h"
\r
29 #include "AliAODMCParticle.h"
\r
30 #include "AliVParticle.h"
\r
31 #include "AliAODEvent.h"
\r
32 #include "AliAODInputHandler.h"
\r
33 #include "AliAnalysisTaskSpectraAllChAOD.h"
\r
34 #include "AliAnalysisTaskESDfilter.h"
\r
35 #include "AliAnalysisDataContainer.h"
\r
36 #include "AliSpectraAODTrackCuts.h"
\r
37 #include "AliSpectraAODEventCuts.h"
\r
38 #include "AliCentrality.h"
\r
40 #include "AliVEvent.h"
\r
41 #include "AliStack.h"
\r
42 #include <TMCProcess.h>
\r
46 using namespace std;
\r
48 ClassImp(AliAnalysisTaskSpectraAllChAOD)
\r
50 //________________________________________________________________________
\r
51 AliAnalysisTaskSpectraAllChAOD::AliAnalysisTaskSpectraAllChAOD(const char *name) : AliAnalysisTaskSE(name), fAOD(0), fTrackCuts(0), fEventCuts(0), fIsMC(0), fOutput(0)
\r
53 // Default constructor
\r
55 DefineInput(0, TChain::Class());
\r
56 DefineOutput(1, TList::Class());
\r
57 DefineOutput(2, AliSpectraAODEventCuts::Class());
\r
58 DefineOutput(3, AliSpectraAODTrackCuts::Class());
\r
61 //________________________________________________________________________
\r
62 //________________________________________________________________________
\r
63 void AliAnalysisTaskSpectraAllChAOD::UserCreateOutputObjects()
\r
65 // create output objects
\r
66 fOutput=new TList();
\r
67 fOutput->SetOwner();
\r
68 fOutput->SetName("fOutput");
\r
70 if (!fTrackCuts) AliFatal("Track Cuts should be set in the steering macro");
\r
71 if (!fEventCuts) AliFatal("Event Cuts should be set in the steering macro");
\r
73 //dimensions of standard THnSparse
\r
75 const Double_t ptBins[] = {0.05,0.1,0.12,0.14,0.16,0.18,0.20,0.25,0.30,0.35,0.4,0.45,0.5,0.55,0.6,0.65,0.7,0.75,0.8,0.85,0.9,0.95,1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0,2.1,2.2,2.3,2.4,2.5,2.6,2.7,2.8,2.9,3.0,3.2,3.4,3.6,3.8,4.0,4.2,4.4,4.6,4.8,5.0,5.5,6.0,6.5,7.0,7.5,8.0,8.5,9.0,9.5,10.,11.,12.,13.,14.,15.};
\r
78 Int_t binsHistReal[nvar] = { 100, 100, nptBins, 20};
\r
79 Double_t xminHistReal[nvar] = { 0., 0, 0., 0.};
\r
80 Double_t xmaxHistReal[nvar] = { 10., 10., 5., 100.};
\r
82 THnSparseF* NSparseHist = new THnSparseF("NSparseHist","NSparseHist",nvar,binsHistReal,xminHistReal,xmaxHistReal);
\r
83 NSparseHist->GetAxis(0)->SetTitle("Q vec VZERO-C");
\r
84 NSparseHist->GetAxis(1)->SetTitle("Q vec VZERO-A");
\r
85 NSparseHist->GetAxis(2)->SetTitle("#it{p}_{T}");
\r
86 NSparseHist->SetBinEdges(2,ptBins);
\r
87 NSparseHist->GetAxis(3)->SetTitle(Form("%s cent",fEventCuts->GetCentralityMethod().Data()));
\r
88 fOutput->Add(NSparseHist);
\r
90 TH2F* hGen = new TH2F("hGen","hGen",nptBins,ptBins,20,0.,100.);
\r
91 hGen->GetXaxis()->SetTitle("#it{p}_{T,Gen}");
\r
92 hGen->GetYaxis()->SetTitle(Form("%s cent",fEventCuts->GetCentralityMethod().Data()));
\r
95 TH2F* hRec = new TH2F("hRec","hRec",nptBins,ptBins,20,0.,100.);
\r
96 hRec->GetXaxis()->SetTitle("#it{p}_{T,Rec}");
\r
97 hRec->GetYaxis()->SetTitle(Form("%s cent",fEventCuts->GetCentralityMethod().Data()));
\r
100 PostData(1, fOutput );
\r
101 PostData(2, fEventCuts);
\r
102 PostData(3, fTrackCuts);
\r
105 //________________________________________________________________________
\r
106 void AliAnalysisTaskSpectraAllChAOD::UserExec(Option_t *)
\r
109 fAOD = dynamic_cast<AliAODEvent*>(fInputEvent);
\r
111 AliWarning("ERROR: AliAODEvent not available \n");
\r
115 if (strcmp(fAOD->ClassName(), "AliAODEvent"))
\r
117 AliFatal("Not processing AODs");
\r
120 if(!fEventCuts->IsSelected(fAOD,fTrackCuts))return;//event selection
\r
122 // First do MC to fill up the MC particle array, such that we can use it later
\r
123 TClonesArray *arrayMC = 0;
\r
126 arrayMC = (TClonesArray*) fAOD->GetList()->FindObject(AliAODMCParticle::StdBranchName());
\r
128 AliFatal("Error: MC particles branch not found!\n");
\r
130 Int_t nMC = arrayMC->GetEntries();
\r
131 for (Int_t iMC = 0; iMC < nMC; iMC++)
\r
133 AliAODMCParticle *partMC = (AliAODMCParticle*) arrayMC->At(iMC);
\r
134 if(!partMC->Charge()) continue;//Skip neutrals
\r
135 if(!partMC->IsPhysicalPrimary()) continue; //skip secondaries
\r
136 if(partMC->Eta() < fTrackCuts->GetEtaMin() || partMC->Eta() > fTrackCuts->GetEtaMax())continue;//ETA CUT ON GENERATED!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
138 ((TH2F*)fOutput->FindObject("hGen"))->Fill(partMC->Pt(),fEventCuts->GetCent());
\r
139 //Printf("a particle");
\r
144 //main loop on tracks
\r
145 for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) {
\r
146 AliAODTrack* track = fAOD->GetTrack(iTracks);
\r
147 if (!fTrackCuts->IsSelected(track,kTRUE)) continue;
\r
149 ((TH2F*)fOutput->FindObject("hRec"))->Fill(track->Pt(),fEventCuts->GetCent());
\r
152 var[0]=fEventCuts->GetqV0C();
\r
153 var[1]=fEventCuts->GetqV0A();
\r
154 var[2]=track->Pt();
\r
155 var[3]=fEventCuts->GetCent();
\r
156 ((THnSparseF*)fOutput->FindObject("NSparseHist"))->Fill(var);
\r
158 //Printf("a track");
\r
160 } // end loop on tracks
\r
162 PostData(1, fOutput );
\r
163 PostData(2, fEventCuts);
\r
164 PostData(3, fTrackCuts);
\r
167 //_________________________________________________________________
\r
168 void AliAnalysisTaskSpectraAllChAOD::Terminate(Option_t *)
\r