1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 /* AliAnalysisTaskEx01.cxx
20 * Template task producing a P_t spectrum and pseudorapidity distribution.
21 * Includes explanations of physics and primary track selections
23 * Instructions for adding histograms can be found below, starting with NEW HISTO
25 * Based on tutorial example from offline pages
26 * Edited by Arvinder Palaha
28 #include "AliAnalysisTaskEx01.h"
30 #include "Riostream.h"
38 #include "AliAnalysisTaskSE.h"
39 #include "AliAnalysisManager.h"
41 #include "AliESDtrackCuts.h"
42 #include "AliESDEvent.h"
43 #include "AliESDInputHandler.h"
44 #include "AliAODEvent.h"
45 #include "AliMCEvent.h"
47 ClassImp(AliAnalysisTaskEx01)
49 //________________________________________________________________________
50 AliAnalysisTaskEx01::AliAnalysisTaskEx01() // All data members should be initialised here
55 fHistEta(0) // The last in the above list should not have a comma after it
57 // Dummy constructor ALWAYS needed for I/O.
60 //________________________________________________________________________
61 AliAnalysisTaskEx01::AliAnalysisTaskEx01(const char *name) // All data members should be initialised here
62 :AliAnalysisTaskSE(name),
66 fHistEta(0) // The last in the above list should not have a comma after it
69 // Define input and output slots here (never in the dummy constructor)
70 // Input slot #0 works with a TChain - it is connected to the default input container
71 // Output slot #1 writes into a TH1 container
72 DefineOutput(1, TList::Class()); // for output list
75 //________________________________________________________________________
76 AliAnalysisTaskEx01::~AliAnalysisTaskEx01()
78 // Destructor. Clean-up the output list, but not the histograms that are put inside
79 // (the list is owner and will clean-up these histograms). Protect in PROOF case.
80 if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
83 if (fTrackCuts) delete fTrackCuts;
86 //________________________________________________________________________
87 void AliAnalysisTaskEx01::UserCreateOutputObjects()
90 // Called once (on the worker node)
92 fOutput = new TList();
93 fOutput->SetOwner(); // IMPORTANT!
95 fTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kTRUE);
96 // === Primary Track Selection ===
98 // The definition of a primary track is taken from the ALICE Twiki
99 // page https://twiki.cern.ch/twiki/bin/view/ALICE/SelectionOfPrimaryTracksForPpDataAnalysis
100 // using the following parameters for a standard dN/dPt analysis:
101 // track quality cuts:
102 // esdTrackCuts->SetMinNClustersTPC(70);
103 // esdTrackCuts->SetMaxChi2PerClusterTPC(4);
104 // esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
105 // esdTrackCuts->SetRequireTPCRefit(kTRUE);
106 // esdTrackCuts->SetRequireITSRefit(kTRUE);
107 // esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
108 // AliESDtrackCuts::kAny);
110 // if(selPrimaries) {
111 // // 7*(0.0026+0.0050/pt^1.01)
112 // esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0182+0.0350/pt^1.01");
114 // esdTrackCuts->SetMaxDCAToVertexZ(2);
115 // esdTrackCuts->SetDCAToVertex2D(kFALSE);
116 // esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
118 // The Primary Track Selection is implemented here by creating an
119 // AliESDtrackCuts object, with kTRUE argument to choose primary tracks.
121 // By default, it is set to the above conditions which are suitable for
122 // a standard inclusive dN/dPt analysis. For others, such as identified
123 // dN/dPt or strangeness as well as others, follow the above link for
124 // the specific changes to include in the selection.
126 // To change cuts after selecting some default set, one can use
127 // esdtrackcuts->SetMinNClustersTPC(70) for example
131 Float_t ptlow = 0.1, ptup = 3.1;
132 fHistPt = new TH1F("fHistPt", "P_{T} distribution for reconstructed", ptbins, ptlow, ptup);
133 fHistPt->GetXaxis()->SetTitle("P_{T} (GeV/c)");
134 fHistPt->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");
135 fHistPt->SetMarkerStyle(kFullCircle);
138 Float_t etalow = -2.0, etaup = 2.0;
139 fHistEta = new TH1F("fHistEta","#eta distribution for reconstructed",etabins, etalow, etaup);
140 fHistEta->GetXaxis()->SetTitle("#eta");
141 fHistEta->GetYaxis()->SetTitle("counts");
143 // NEW HISTO should be defined here, with a sensible name,
145 fOutput->Add(fHistPt);
146 fOutput->Add(fHistEta);
147 // NEW HISTO added to fOutput here
148 PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram
151 //________________________________________________________________________
152 void AliAnalysisTaskEx01::UserExec(Option_t *)
155 // Called for each event
158 // Create pointer to reconstructed event
159 AliVEvent *event = InputEvent();
160 if (!event) { Printf("ERROR: Could not retrieve event"); return; }
162 // If the task accesses MC info, this can be done as in the commented block below:
164 // Create pointer to reconstructed event
165 AliMCEvent *mcEvent = MCEvent();
166 if (!mcEvent) { Printf("ERROR: Could not retrieve MC event"); return; }
167 Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
169 // set up a stack for use in check for primary/stable particles
170 AliStack* stack = mcEvent->Stack();
171 if( !stack ) { Printf( "Stack not available"); return; }
174 // create pointer to event
175 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(event);
177 AliError("Cannot get the ESD event");
180 // AliESDHeader* esdheader = (AliESDHeader*)esd->GetHeader();
182 // === Physics Selection Task ===
184 // To perform a physics selection here, a bitwise operation is used against
185 // the UInt_t mask which is extracted in the following way:
187 // UInt_t mask = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
189 // This can be tested to produce the following
191 // Bool_t bMinBias = (mask == AliVEvent::kMB) ? 1 : 0; // check if minimum bias trigger class fired
192 // Bool_t bHighMult = (mask == AliVEvent::kHighMult) ? 1 : 0; // check if high multiplicity trigger class fired
194 // For more complicated trigger selections, one can directly test both
195 // trigger classes and fired trigger inputs for a particular event, for e.g.
197 // Bool_t bCSH1 = (esd->IsTriggerClassFired("CSH1-B-NOPF-ALLNOTRD")) ? 1 : 0;
198 // Bool_t b0SH1 = (esdheader->IsTriggerInputFired("0SH1")) ? 1 : 0;
200 // These booleans can then be used to fill different histograms for specific
201 // conditions, or summed to make one cut for all events that fill histos.
203 // Do some fast cuts first
204 // check for good reconstructed vertex
205 if(!(esd->GetPrimaryVertex()->GetStatus())) return;
206 // if vertex is from spd vertexZ, require more stringent cut
207 if (esd->GetPrimaryVertex()->IsFromVertexerZ()) {
208 if (esd->GetPrimaryVertex()->GetDispersion()>0.02 || esd->GetPrimaryVertex()->GetZRes()>0.25 ) return; // bad vertex from VertexerZ
211 // Track loop for reconstructed event
212 Int_t ntracks = esd->GetNumberOfTracks();
213 for(Int_t i = 0; i < ntracks; i++) {
214 AliESDtrack* esdtrack = esd->GetTrack(i); // pointer to reconstructed to track
216 AliError(Form("ERROR: Could not retrieve esdtrack %d",i));
220 // Some MC checks, if MC is used
221 //if(esdtrack->GetLabel() < 0) continue; // get rid of "ghost" tracks
223 // ... and the thorough checking of ESD cuts after.
224 // if this is not a primary track, skip to the next one
225 if(!fTrackCuts->AcceptTrack(esdtrack)) continue;
227 fHistPt->Fill(esdtrack->Pt());
228 fHistEta->Fill(esdtrack->Eta());
230 // NEW HISTO should be filled before this point, as PostData puts the
231 // information for this iteration of the UserExec in the container
232 PostData(1, fOutput);
236 //________________________________________________________________________
237 void AliAnalysisTaskEx01::Terminate(Option_t *)
239 // Draw result to screen, or perform fitting, normalizations
240 // Called once at the end of the query
242 fOutput = dynamic_cast<TList*> (GetOutputData(1));
243 if(!fOutput) { Printf("ERROR: could not retrieve TList fOutput"); return; }
245 fHistPt = dynamic_cast<TH1F*> (fOutput->FindObject("fHistPt"));
246 if (!fHistPt) { Printf("ERROR: could not retrieve fHistPt"); return;}
247 fHistEta = dynamic_cast<TH1F*> (fOutput->FindObject("fHistEta"));
248 if (!fHistEta) { Printf("ERROR: could not retrieve fHistEta"); return;}
250 // Get the physics selection histograms with the selection statistics
251 //AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
252 //AliESDInputHandler *inputH = dynamic_cast<AliESDInputHandler*>(mgr->GetInputEventHandler());
253 //TH2F *histStat = (TH2F*)inputH->GetStatistics();
256 // NEW HISTO should be retrieved from the TList container in the above way,
257 // so it is available to draw on a canvas such as below
259 TCanvas *c = new TCanvas("AliAnalysisTaskEx01","P_{T} & #eta",10,10,1020,510);
262 fHistPt->DrawCopy("E");
264 fHistEta->DrawCopy("E");