]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGGA/GammaConv/AliAnalysisTaskResolution.cxx
changed Resolution, Material, and PhotonQA task to be able to run on the grid
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAnalysisTaskResolution.cxx
CommitLineData
72395bd9 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Authors: Friederike Bock *
5 * Version 1.0 *
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////////////////////////////////////////////////
17//---------------------------------------------
18// QA Task for V0 Reader V1
19//---------------------------------------------
20////////////////////////////////////////////////
21
22#include "AliAnalysisTaskResolution.h"
23#include "TChain.h"
24#include "AliAnalysisManager.h"
25#include "TParticle.h"
26#include "TVectorF.h"
27#include "AliPIDResponse.h"
28#include "AliESDtrackCuts.h"
29#include "TFile.h"
30
31class iostream;
32
33using namespace std;
34
35ClassImp(AliAnalysisTaskResolution)
36
37//________________________________________________________________________
38AliAnalysisTaskResolution::AliAnalysisTaskResolution() : AliAnalysisTaskSE(),
39 fV0Reader(NULL),
40 fConversionGammas(NULL),
41 fConversionCuts(NULL),
42 fTreeEvent(NULL),
43 fTreeResolution(NULL),
44 fPrimVtxZ(0.),
45 fNContrVtx(0),
46 fNESDtracksEta09(0),
47 fNESDtracksEta0914(0),
48 fNESDtracksEta14(0),
49// fGammaRecCoords(NULL),
50// fGammaMCCoords(NULL),
51 fChi2ndf(0),
52 fIsHeavyIon(kFALSE),
53 fOutputList(NULL),
54 fEventList(NULL),
55 fResolutionList(NULL),
56 fESDEvent(NULL),
57 fMCEvent(NULL)
58{
59
60
61}
62
63
64//________________________________________________________________________
65AliAnalysisTaskResolution::AliAnalysisTaskResolution(const char *name) : AliAnalysisTaskSE(name),
66 fV0Reader(NULL),
67 fConversionGammas(NULL),
68 fConversionCuts(NULL),
69 fTreeEvent(NULL),
70 fTreeResolution(NULL),
71 fPrimVtxZ(0.),
72 fNContrVtx(0),
73 fNESDtracksEta09(0),
74 fNESDtracksEta0914(0),
75 fNESDtracksEta14(0),
76// fGammaRecCoords(NULL),
77// fGammaMCCoords(NULL),
78 fChi2ndf(0),
79 fIsHeavyIon(kFALSE),
80 fOutputList(NULL),
81 fEventList(NULL),
82 fResolutionList(NULL),
83 fESDEvent(NULL),
84 fMCEvent(NULL)
85{
86 // Default constructor
87
88 DefineInput(0, TChain::Class());
89 DefineOutput(1, TList::Class());
90}
91
92//________________________________________________________________________
93AliAnalysisTaskResolution::~AliAnalysisTaskResolution()
94{
95 // default deconstructor
96
97}
98//________________________________________________________________________
99void AliAnalysisTaskResolution::UserCreateOutputObjects()
100{
101 // Create User Output Objects
102
103 if(fOutputList != NULL){
104 delete fOutputList;
105 fOutputList = NULL;
106 }
107 if(fOutputList == NULL){
108 fOutputList = new TList();
109 fOutputList->SetOwner(kTRUE);
110 }
111
112 fEventList = new TList();
113 fEventList->SetName("EventList");
114 fEventList->SetOwner(kTRUE);
115 fOutputList->Add(fEventList);
116
117 fTreeEvent = new TTree("Event","Event");
118 fTreeEvent->Branch("primVtxZ",&fPrimVtxZ,"fPrimVtxZ/F");
119 fTreeEvent->Branch("nContrVtx",&fNContrVtx,"fNContrVtx/I");
120 fTreeEvent->Branch("nGoodTracksEta09",&fNESDtracksEta09,"fNESDtracksEta09/I");
121 fTreeEvent->Branch("nGoodTracksEta0914",&fNESDtracksEta0914,"fNESDtracksEta0914/I");
122 fTreeEvent->Branch("nGoodTracksEta14",&fNESDtracksEta14,"fNESDtracksEta14/I");
123 fEventList->Add(fTreeEvent);
124
125 fResolutionList = new TList();
126 fResolutionList->SetName("ResolutionList");
127 fResolutionList->SetOwner(kTRUE);
128 fOutputList->Add(fResolutionList);
129
130// fGammaRecCoords = new Float_t[5];
131// fGammaMCCoords = new Float_t[5];
132
133 fTreeResolution = new TTree("Resolution","Resolution");
134 fTreeResolution->Branch("RecCoords",fGammaRecCoords,"fGammaRecCoords[5]/F");
135 fTreeResolution->Branch("MCCoords",fGammaMCCoords,"fGammaMCCoords[5]/F");
136 fTreeResolution->Branch("chi2ndf",&fChi2ndf,"fChi2ndf/F");
137 fResolutionList->Add(fTreeResolution);
138 // V0 Reader Cuts
139 TString cutnumber = fConversionCuts->GetCutNumber();
140
141 PostData(1, fOutputList);
142}
143
144//________________________________________________________________________
145void AliAnalysisTaskResolution::UserExec(Option_t *){
146
147 fV0Reader=(AliV0ReaderV1*)AliAnalysisManager::GetAnalysisManager()->GetTask("V0ReaderV1");
148
149 Int_t eventQuality = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetEventQuality();
150 if(eventQuality != 0){// Event Not Accepted
151 return;
152 }
153 fESDEvent = (AliESDEvent*) InputEvent();
154 if (fESDEvent==NULL) return;
155 if(fIsHeavyIon && !fConversionCuts->IsCentralitySelected(fESDEvent)) return;
156 fNESDtracksEta09 = CountTracks09(); // Estimate Event Multiplicity
157 fNESDtracksEta0914 = CountTracks0914(); // Estimate Event Multiplicity
158 fNESDtracksEta14 = fNESDtracksEta09 + fNESDtracksEta0914;
159 if(fESDEvent){
160 if(fESDEvent->GetPrimaryVertexTracks()->GetNContributors()>0) {
161 fNContrVtx = fESDEvent->GetPrimaryVertexTracks()->GetNContributors();
162 } else if(fESDEvent->GetPrimaryVertexTracks()->GetNContributors()<1) {
163 if(fESDEvent->GetPrimaryVertexSPD()->GetNContributors()>0) {
164 fNContrVtx = fESDEvent->GetPrimaryVertexSPD()->GetNContributors();
165 } else if(fESDEvent->GetPrimaryVertexSPD()->GetNContributors()<1) {
166 fNContrVtx = 0;
167 }
168 }
169 }
170 fPrimVtxZ = fESDEvent->GetPrimaryVertex()->GetZ();
171
172 if (fTreeEvent){
173 fTreeEvent->Fill();
174 }
175
176 fConversionGammas=fV0Reader->GetReconstructedGammas();
177 if(MCEvent()){
178 fMCEvent = MCEvent();
179 }
180 ProcessPhotons();
181 PostData(1, fOutputList);
182}
183
184
185///________________________________________________________________________
186void AliAnalysisTaskResolution::ProcessPhotons(){
187
188 // Fill Histograms for QA and MC
189 for(Int_t firstGammaIndex=0;firstGammaIndex<fConversionGammas->GetEntriesFast();firstGammaIndex++){
190 AliAODConversionPhoton *gamma=dynamic_cast<AliAODConversionPhoton*>(fConversionGammas->At(firstGammaIndex));
191 if (gamma ==NULL) continue;
192 if(!fConversionCuts->PhotonIsSelected(gamma,fESDEvent)) continue;
193 fGammaRecCoords[0] = gamma->GetPhotonPt();
194 fGammaRecCoords[1] = gamma->GetPhotonPhi();
195 fGammaRecCoords[2] = gamma->GetPhotonEta();
196 fGammaRecCoords[3] = gamma->GetConversionRadius();
197 fGammaRecCoords[4] = gamma->GetConversionZ();
198 fChi2ndf = gamma->GetChi2perNDF();
199 if(fMCEvent){
200// cout << "generating MC stack"<< endl;
201 AliStack *fMCStack = fMCEvent->Stack();
202 if (!fMCStack) continue;
203 TParticle *posDaughter = gamma->GetPositiveMCDaughter(fMCStack);
204 TParticle *negDaughter = gamma->GetNegativeMCDaughter(fMCStack);
205// cout << "generate Daughters: "<<posDaughter << "\t" << negDaughter << endl;
206
207 if(posDaughter == NULL || negDaughter == NULL){
208 continue;
209 } else if(posDaughter->GetMother(0) != negDaughter->GetMother(0) || (posDaughter->GetMother(0) == negDaughter->GetMother(0) && posDaughter->GetMother(0) ==-1)){
210 continue;
211 } else {
212// cout << "same mother" << endl;
213 Int_t pdgCodePos;
214 if (posDaughter->GetPdgCode()) pdgCodePos = posDaughter->GetPdgCode(); else continue;
215 Int_t pdgCodeNeg;
216 if (negDaughter->GetPdgCode()) pdgCodeNeg = negDaughter->GetPdgCode(); else continue;
217// cout << "PDG codes daughters: " << pdgCodePos << "\t" << pdgCodeNeg << endl;
218 Int_t pdgCode;
219 if (gamma->GetMCParticle(fMCStack)->GetPdgCode()) pdgCode = gamma->GetMCParticle(fMCStack)->GetPdgCode(); else continue;
220// cout << "PDG code: " << pdgCode << endl;
221 if(TMath::Abs(pdgCodePos)!=11 || TMath::Abs(pdgCodeNeg)!=11)
222 continue;
223 else if ( !(pdgCodeNeg==pdgCodePos)){
224 TParticle *truePhotonCanditate = gamma->GetMCParticle(fMCStack);
225 if(pdgCode == 111)
226 continue;
227 else if (pdgCode == 221)
228 continue;
229 else if (!(negDaughter->GetUniqueID() != 5 || posDaughter->GetUniqueID() !=5)){
230 if(pdgCode == 22){
231 fGammaMCCoords[0] = truePhotonCanditate->Pt();
232 fGammaMCCoords[1] = gamma->GetNegativeMCDaughter(fMCStack)->Phi();
233 fGammaMCCoords[2] = gamma->GetNegativeMCDaughter(fMCStack)->Eta();
234 fGammaMCCoords[3] = gamma->GetNegativeMCDaughter(fMCStack)->R();
235 fGammaMCCoords[4] = gamma->GetNegativeMCDaughter(fMCStack)->Vz();
236
237 if (fTreeResolution){
238 fTreeResolution->Fill();
239 }
240 }
241 } else continue; //garbage
242 } else continue; //garbage
243 }
244 }
245 }
246}
247
248//________________________________________________________________________
249Int_t AliAnalysisTaskResolution::CountTracks09(){
250 Int_t fNumberOfESDTracks = 0;
251 if(fInputEvent->IsA()==AliESDEvent::Class()){
252 // Using standard function for setting Cuts
253 Bool_t selectPrimaries=kTRUE;
254 AliESDtrackCuts *EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
255 EsdTrackCuts->SetMaxDCAToVertexZ(2);
256 EsdTrackCuts->SetEtaRange(-0.9, 0.9);
257 EsdTrackCuts->SetPtRange(0.15);
258
259 for(Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++){
260 AliESDtrack* curTrack = (AliESDtrack*) fInputEvent->GetTrack(iTracks);
261 if(!curTrack) continue;
262 // if(fMCEvent && ((AliConversionCuts*)fCutArray->At(fiCut))->GetSignalRejection() != 0){
263 // if(!((AliConversionCuts*)fCutArray->At(fiCut))->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack)) continue;
264 // }
265 if(EsdTrackCuts->AcceptTrack(curTrack) ) fNumberOfESDTracks++;
266 }
267 delete EsdTrackCuts;
268 EsdTrackCuts=0x0;
269 }
270 else if(fInputEvent->IsA()==AliAODEvent::Class()){
271 for(Int_t iTracks = 0; iTracks<fInputEvent->GetNumberOfTracks(); iTracks++){
272 AliAODTrack* curTrack = (AliAODTrack*) fInputEvent->GetTrack(iTracks);
273 if(!curTrack->IsPrimaryCandidate()) continue;
274 if(abs(curTrack->Eta())>0.9) continue;
275 if(curTrack->Pt()<0.15) continue;
276 if(abs(curTrack->ZAtDCA())>2) continue;
277 fNumberOfESDTracks++;
278 }
279 }
280
281 return fNumberOfESDTracks;
282}
283
284Int_t AliAnalysisTaskResolution::CountTracks0914(){
285
286 // Using standard function for setting Cuts ; We use TPCOnlyTracks for outer eta region
287 //Bool_t selectPrimaries=kTRUE;
288 // EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
289 Int_t fNumberOfESDTracks = 0;
290 if(fInputEvent->IsA()==AliESDEvent::Class()){
291 // Using standard function for setting Cuts
292 AliESDtrackCuts *EsdTrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
293 EsdTrackCuts->SetMaxDCAToVertexZ(5);
294 EsdTrackCuts->SetEtaRange(0.9, 1.4);
295 EsdTrackCuts->SetPtRange(0.15);
296
297 for(Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++){
298 AliESDtrack* curTrack = (AliESDtrack*) fInputEvent->GetTrack(iTracks);
299 if(!curTrack) continue;
300 // if(fMCEvent && ((AliConversionCuts*)fCutArray->At(fiCut))->GetSignalRejection() != 0){
301 // if(!((AliConversionCuts*)fCutArray->At(fiCut))->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack)) continue;
302 // }
303 if(EsdTrackCuts->AcceptTrack(curTrack) ) fNumberOfESDTracks++;
304 }
305 EsdTrackCuts->SetEtaRange(-1.4, -0.9);
306 for(Int_t iTracks = 0; iTracks < fESDEvent->GetNumberOfTracks(); iTracks++){
307 AliESDtrack* curTrack = fESDEvent->GetTrack(iTracks);
308 if(!curTrack) continue;
309 if(EsdTrackCuts->AcceptTrack(curTrack) ) fNumberOfESDTracks++;
310 }
311 delete EsdTrackCuts;
312 EsdTrackCuts=0x0;
313 }
314 else if(fInputEvent->IsA()==AliAODEvent::Class()){
315 for(Int_t iTracks = 0; iTracks<fInputEvent->GetNumberOfTracks(); iTracks++){
316 AliAODTrack* curTrack = (AliAODTrack*) fInputEvent->GetTrack(iTracks);
317// if(!curTrack->IsPrimaryCandidate()) continue;
318 if(abs(curTrack->Eta())<0.9 || abs(curTrack->Eta())>1.4 ) continue;
319 if(curTrack->Pt()<0.15) continue;
320 if(abs(curTrack->ZAtDCA())>5) continue;
321 fNumberOfESDTracks++;
322 }
323 }
324
325 return fNumberOfESDTracks;
326}
327
328//________________________________________________________________________
329void AliAnalysisTaskResolution::Terminate(Option_t *)
330{
331// if (fStreamMaterial){
332// fStreamMaterial->GetFile()->Write();
333// }
334// if (fStreamResolution){
335// fStreamResolution->GetFile()->Write();
336// }
337}