1 /*************************************************************************
2 * Copyright(c) 1998-2008, 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 **************************************************************************/
23 #include "TClonesArray.h"
24 #include "TParticle.h"
25 #include "TObjString.h"
27 #include "TDatabasePDG.h"
28 #include "TLorentzVector.h"
31 #include "AliAnalysisManager.h"
32 #include "AliInputEventHandler.h"
33 #include "AliESDEvent.h"
34 #include "AliAODEvent.h"
35 #include "AliMCEvent.h"
36 #include "AliAODVZERO.h"
37 #include "AliAODZDC.h"
38 #include "AliESDVZERO.h"
39 #include "AliESDZDC.h"
40 #include "AliPIDResponse.h"
41 #include "AliAODTrack.h"
42 #include "AliAODPid.h"
43 #include "AliAODVertex.h"
44 #include "AliESDVertex.h"
45 #include "AliMultiplicity.h"
46 #include "AliESDtrack.h"
47 #include "AliESDMuonTrack.h"
48 #include "AliAODMCParticle.h"
49 #include "AliMCParticle.h"
50 #include "AliCentrality.h"
51 #include "AliKFVertex.h"
52 #include "AliExternalTrackParam.h"
53 #include "AliTriggerAnalysis.h"
56 #include "AliAnalysisTaskUpcPsi2s.h"
58 ClassImp(AliAnalysisTaskUpcPsi2s);
63 //trees for UPC analysis,
64 // michal.broz@cern.ch
66 //_____________________________________________________________________________
67 AliAnalysisTaskUpcPsi2s::AliAnalysisTaskUpcPsi2s()
68 : AliAnalysisTaskSE(),fType(0),isMC(kFALSE),fRunTree(kTRUE),fRunHist(kTRUE),fPIDResponse(0),fJPsiTree(0),fPsi2sTree(0),
69 fRunNum(0),fPerNum(0),fOrbNum(0),fL0inputs(0),fL1inputs(0),
70 fTOFtrig1(0), fTOFtrig2(0),
71 fVtxContrib(0),fVtxChi2(0),fVtxNDF(0),
72 fBCrossNum(0),fNtracklets(0),fNLooseTracks(0),
73 fZDCAenergy(0),fZDCCenergy(0),fV0Adecision(0),fV0Cdecision(0),
74 fDataFilnam(0),fRecoPass(0),fEvtNum(0),
75 fJPsiAODTracks(0),fJPsiESDTracks(0),fPsi2sAODTracks(0),fPsi2sESDTracks(0),fGenPart(0),
76 fListTrig(0),fHistCcup4TriggersPerRun(0), fHistCcup7TriggersPerRun(0), fHistCcup2TriggersPerRun(0),
77 fHistZedTriggersPerRun(0),fHistCvlnTriggersPerRun(0), fHistMBTriggersPerRun(0),fHistCentralTriggersPerRun(0),fHistSemiCentralTriggersPerRun(0),
78 fListHist(0),fHistNeventsJPsi(0),fHistTPCsignalJPsi(0),fHistDiLeptonPtJPsi(0),fHistDiElectronMass(0),fHistDiMuonMass(0),
79 fHistNeventsPsi2s(0),fHistPsi2sMassVsPt(0),fHistPsi2sMassCoherent(0)
85 }//AliAnalysisTaskUpcPsi2s
88 //_____________________________________________________________________________
89 AliAnalysisTaskUpcPsi2s::AliAnalysisTaskUpcPsi2s(const char *name)
90 : AliAnalysisTaskSE(name),fType(0),isMC(kFALSE),fRunTree(kTRUE),fRunHist(kTRUE),fPIDResponse(0),fJPsiTree(0),fPsi2sTree(0),
91 fRunNum(0),fPerNum(0),fOrbNum(0),fL0inputs(0),fL1inputs(0),
92 fTOFtrig1(0), fTOFtrig2(0),
93 fVtxContrib(0),fVtxChi2(0),fVtxNDF(0),
94 fBCrossNum(0),fNtracklets(0),fNLooseTracks(0),
95 fZDCAenergy(0),fZDCCenergy(0),fV0Adecision(0),fV0Cdecision(0),
96 fDataFilnam(0),fRecoPass(0),fEvtNum(0),
97 fJPsiAODTracks(0),fJPsiESDTracks(0),fPsi2sAODTracks(0),fPsi2sESDTracks(0),fGenPart(0),
98 fListTrig(0),fHistCcup4TriggersPerRun(0), fHistCcup7TriggersPerRun(0), fHistCcup2TriggersPerRun(0),
99 fHistZedTriggersPerRun(0),fHistCvlnTriggersPerRun(0), fHistMBTriggersPerRun(0),fHistCentralTriggersPerRun(0),fHistSemiCentralTriggersPerRun(0),
100 fListHist(0),fHistNeventsJPsi(0),fHistTPCsignalJPsi(0),fHistDiLeptonPtJPsi(0),fHistDiElectronMass(0),fHistDiMuonMass(0),
101 fHistNeventsPsi2s(0),fHistPsi2sMassVsPt(0),fHistPsi2sMassCoherent(0)
106 if( strstr(name,"ESD") ) fType = 0;
107 if( strstr(name,"AOD") ) fType = 1;
111 DefineOutput(1, TTree::Class());
112 DefineOutput(2, TTree::Class());
113 DefineOutput(3, TList::Class());
114 DefineOutput(4, TList::Class());
116 }//AliAnalysisTaskUpcPsi2s
118 //_____________________________________________________________________________
119 void AliAnalysisTaskUpcPsi2s::Init()
122 for(Int_t i=0; i<ntrg; i++) fTrigger[i] = kFALSE;
123 for(Int_t i=0; i<4; i++) {
126 fPIDElectron[i] = -666;
128 fTriggerInputsMC[i] = kFALSE;
130 for(Int_t i=0; i<3; i++){
138 //_____________________________________________________________________________
139 AliAnalysisTaskUpcPsi2s::~AliAnalysisTaskUpcPsi2s()
159 }//~AliAnalysisTaskUpcPsi2s
162 //_____________________________________________________________________________
163 void AliAnalysisTaskUpcPsi2s::UserCreateOutputObjects()
166 AliAnalysisManager *man = AliAnalysisManager::GetAnalysisManager();
167 AliInputEventHandler *inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
168 fPIDResponse = inputHandler->GetPIDResponse();
171 fDataFilnam = new TObjString();
172 fDataFilnam->SetString("");
175 fJPsiAODTracks = new TClonesArray("AliAODTrack", 1000);
176 fJPsiESDTracks = new TClonesArray("AliESDtrack", 1000);
177 fPsi2sAODTracks = new TClonesArray("AliAODTrack", 1000);
178 fPsi2sESDTracks = new TClonesArray("AliESDtrack", 1000);
179 fGenPart = new TClonesArray("TParticle", 1000);
181 //output tree with JPsi candidate events
182 fJPsiTree = new TTree("fJPsiTree", "fJPsiTree");
183 fJPsiTree ->Branch("fRunNum", &fRunNum, "fRunNum/I");
184 fJPsiTree ->Branch("fPerNum", &fPerNum, "fPerNum/i");
185 fJPsiTree ->Branch("fOrbNum", &fOrbNum, "fOrbNum/i");
187 fJPsiTree ->Branch("fBCrossNum", &fBCrossNum, "fBCrossNum/s");
188 fJPsiTree ->Branch("fTrigger", &fTrigger[0], Form("fTrigger[%i]/O", ntrg));
189 fJPsiTree ->Branch("fL0inputs", &fL0inputs, "fL0inputs/i");
190 fJPsiTree ->Branch("fL1inputs", &fL1inputs, "fL1inputs/i");
191 fJPsiTree ->Branch("fNtracklets", &fNtracklets, "fNtracklets/s");
192 fJPsiTree ->Branch("fNLooseTracks", &fNLooseTracks, "fNLooseTracks/s");
193 fJPsiTree ->Branch("fVtxContrib", &fVtxContrib, "fVtxContrib/I");
195 fJPsiTree ->Branch("fTOFtrig1", &fTOFtrig1, "fTOFtrig1/O");
196 fJPsiTree ->Branch("fTOFtrig2", &fTOFtrig2, "fTOFtrig2/O");
197 fJPsiTree ->Branch("fTOFphi", &fTOFphi[0], "fTOFphi[4]/D");
199 fJPsiTree ->Branch("fPIDMuon", &fPIDMuon[0], "fPIDMuon[4]/D");
200 fJPsiTree ->Branch("fPIDElectron", &fPIDElectron[0], "fPIDElectron[4]/D");
201 fJPsiTree ->Branch("fPIDPion", &fPIDPion[0], "fPIDPion[4]/D");
203 fJPsiTree ->Branch("fVtxPos", &fVtxPos[0], "fVtxPos[3]/D");
204 fJPsiTree ->Branch("fVtxErr", &fVtxErr[0], "fVtxErr[3]/D");
205 fJPsiTree ->Branch("fVtxChi2", &fVtxChi2, "fVtxChi2/D");
206 fJPsiTree ->Branch("fVtxNDF", &fVtxNDF, "fVtxNDF/D");
208 fJPsiTree ->Branch("fKfVtxPos", &fKfVtxPos[0], "fKfVtxPos[3]/D");
210 fJPsiTree ->Branch("fZDCAenergy", &fZDCAenergy, "fZDCAenergy/D");
211 fJPsiTree ->Branch("fZDCCenergy", &fZDCCenergy, "fZDCCenergy/D");
212 fJPsiTree ->Branch("fV0Adecision", &fV0Adecision, "fV0Adecision/I");
213 fJPsiTree ->Branch("fV0Cdecision", &fV0Cdecision, "fV0Cdecision/I");
214 fJPsiTree ->Branch("fDataFilnam", &fDataFilnam);
215 fJPsiTree ->Branch("fRecoPass", &fRecoPass, "fRecoPass/S");
216 fJPsiTree ->Branch("fEvtNum", &fEvtNum, "fEvtNum/L");
218 fJPsiTree ->Branch("fJPsiESDTracks", &fJPsiESDTracks);
221 fJPsiTree ->Branch("fJPsiAODTracks", &fJPsiAODTracks);
224 fJPsiTree ->Branch("fGenPart", &fGenPart);
225 fJPsiTree ->Branch("fTriggerInputsMC", &fTriggerInputsMC[0], "fTriggerInputsMC[4]/O");
229 //output tree with Psi2s candidate events
230 fPsi2sTree = new TTree("fPsi2sTree", "fPsi2sTree");
231 fPsi2sTree ->Branch("fRunNum", &fRunNum, "fRunNum/I");
232 fPsi2sTree ->Branch("fPerNum", &fPerNum, "fPerNum/i");
233 fPsi2sTree ->Branch("fOrbNum", &fOrbNum, "fOrbNum/i");
235 fPsi2sTree ->Branch("fBCrossNum", &fBCrossNum, "fBCrossNum/s");
236 fPsi2sTree ->Branch("fTrigger", &fTrigger[0], Form("fTrigger[%i]/O", ntrg));
237 fPsi2sTree ->Branch("fL0inputs", &fL0inputs, "fL0inputs/i");
238 fPsi2sTree ->Branch("fL1inputs", &fL1inputs, "fL1inputs/i");
239 fPsi2sTree ->Branch("fNtracklets", &fNtracklets, "fNtracklets/s");
240 fPsi2sTree ->Branch("fNLooseTracks", &fNLooseTracks, "fNLooseTracks/s");
241 fPsi2sTree ->Branch("fVtxContrib", &fVtxContrib, "fVtxContrib/I");
243 fPsi2sTree ->Branch("fTOFtrig1", &fTOFtrig1, "fTOFtrig1/O");
244 fPsi2sTree ->Branch("fTOFtrig2", &fTOFtrig2, "fTOFtrig2/O");
245 fPsi2sTree ->Branch("fTOFphi", &fTOFphi[0], "fTOFphi[4]/D");
247 fPsi2sTree ->Branch("fPIDMuon", &fPIDMuon[0], "fPIDMuon[4]/D");
248 fPsi2sTree ->Branch("fPIDElectron", &fPIDElectron[0], "fPIDElectron[4]/D");
249 fPsi2sTree ->Branch("fPIDPion", &fPIDPion[0], "fPIDPion[4]/D");
251 fPsi2sTree ->Branch("fVtxPos", &fVtxPos[0], "fVtxPos[3]/D");
252 fPsi2sTree ->Branch("fVtxErr", &fVtxErr[0], "fVtxErr[3]/D");
253 fPsi2sTree ->Branch("fVtxChi2", &fVtxChi2, "fVtxChi2/D");
254 fPsi2sTree ->Branch("fVtxNDF", &fVtxNDF, "fVtxNDF/D");
256 fPsi2sTree ->Branch("fKfVtxPos", &fKfVtxPos[0], "fKfVtxPos[3]/D");
258 fPsi2sTree ->Branch("fZDCAenergy", &fZDCAenergy, "fZDCAenergy/D");
259 fPsi2sTree ->Branch("fZDCCenergy", &fZDCCenergy, "fZDCCenergy/D");
260 fPsi2sTree ->Branch("fV0Adecision", &fV0Adecision, "fV0Adecision/I");
261 fPsi2sTree ->Branch("fV0Cdecision", &fV0Cdecision, "fV0Cdecision/I");
262 fPsi2sTree ->Branch("fDataFilnam", &fDataFilnam);
263 fPsi2sTree ->Branch("fRecoPass", &fRecoPass, "fRecoPass/S");
264 fPsi2sTree ->Branch("fEvtNum", &fEvtNum, "fEvtNum/L");
266 fPsi2sTree ->Branch("fPsi2sESDTracks", &fPsi2sESDTracks);
269 fPsi2sTree ->Branch("fPsi2sAODTracks", &fPsi2sAODTracks);
272 fPsi2sTree ->Branch("fGenPart", &fGenPart);
273 fPsi2sTree ->Branch("fTriggerInputsMC", &fTriggerInputsMC[0], "fTriggerInputsMC[4]/O");
277 fListTrig = new TList();
278 fListTrig ->SetOwner();
280 fHistCcup4TriggersPerRun = new TH1D("fHistCcup4TriggersPerRun", "fHistCcup4TriggersPerRun", 33000, 167000.5, 200000.5);
281 fListTrig->Add(fHistCcup4TriggersPerRun);
283 fHistCcup7TriggersPerRun = new TH1D("fHistCcup7TriggersPerRun", "fHistCcup7TriggersPerRun", 33000, 167000.5, 200000.5);
284 fListTrig->Add(fHistCcup7TriggersPerRun);
286 fHistCcup2TriggersPerRun = new TH1D("fHistCcup2TriggersPerRun", "fHistCcup2TriggersPerRun", 33000, 167000.5, 200000.5);
287 fListTrig->Add(fHistCcup2TriggersPerRun);
289 fHistZedTriggersPerRun = new TH1D("fHistZedTriggersPerRun", "fHistZedTriggersPerRun", 33000, 167000.5, 200000.5);
290 fListTrig->Add(fHistZedTriggersPerRun);
292 fHistCvlnTriggersPerRun = new TH1D("fHistCvlnTriggersPerRun", "fHistCvlnTriggersPerRun", 33000, 167000.5, 200000.5);
293 fListTrig->Add(fHistCvlnTriggersPerRun);
295 fHistMBTriggersPerRun = new TH1D("fHistMBTriggersPerRun", "fHistMBTriggersPerRun", 33000, 167000.5, 200000.5);
296 fListTrig->Add(fHistMBTriggersPerRun);
298 fHistCentralTriggersPerRun = new TH1D("fHistCentralTriggersPerRun", "fHistCentralTriggersPerRun", 33000, 167000.5, 200000.5);
299 fListTrig->Add(fHistCentralTriggersPerRun);
301 fHistSemiCentralTriggersPerRun = new TH1D("fHistSemiCentralTriggersPerRun", "fHistSemiCentralTriggersPerRun", 33000, 167000.5, 200000.5);
302 fListTrig->Add(fHistSemiCentralTriggersPerRun);
304 fListHist = new TList();
305 fListHist ->SetOwner();
307 TString CutNameJPsi[13] = {"Analyzed","Triggered","Vertex cut","V0 decision","Neutron ZDC cut","Two good tracks",
308 "Like sign","Oposite sign","One p_{T}>1", "Both p_{T}>1","PID","Dimuom","Dielectron"};
309 fHistNeventsJPsi = new TH1D("fHistNeventsJPsi","fHistNeventsPsi2s",13,0.5,13.5);
310 for (Int_t i = 0; i<13; i++) fHistNeventsJPsi->GetXaxis()->SetBinLabel(i+1,CutNameJPsi[i].Data());
311 fListHist->Add(fHistNeventsJPsi);
313 fHistTPCsignalJPsi = new TH2D("fHistTPCsignalJPsi","fHistTPCsignalJPsi",240,0,120,240,0,120);
314 fListHist->Add(fHistTPCsignalJPsi);
316 fHistDiLeptonPtJPsi = new TH2D("fHistDiLeptonPtJPsi","fHistDiLeptonPtJPsi",350,0,3.5,350,0,3.5);
317 fListHist->Add(fHistDiLeptonPtJPsi);
319 fHistDiElectronMass = new TH1D("fHistDiElectronMass","Invariant mass of J/#psi candidates",100,2,5);
320 fHistDiElectronMass->GetXaxis()->SetTitle("Invariant mass(e^{+}e^{-}) (GeV/c)");
321 fListHist->Add(fHistDiElectronMass);
323 fHistDiMuonMass = new TH1D("fHistDiMuonMass","Invariant mass of J/#psi candidates",100,2,5);
324 fHistDiMuonMass->GetXaxis()->SetTitle("Invariant mass(#mu^{+}#mu^{-}) (GeV/c)");
325 fListHist->Add(fHistDiMuonMass);
327 TString CutNamePsi2s[14] = {"Analyzed","Triggered","Vertex cut","V0 decision","Neutron ZDC cut","Four good tracks",
328 "DiLepton - DiPion","Like sign leptons","Like sign pions","Like sign both","Oposite sign","PID","Dimuom","Dielectron"};
330 fHistNeventsPsi2s = new TH1D("fHistNeventsPsi2s","fHistNeventsPsi2s",14,0.5,14.5);
331 for (Int_t i = 0; i<14; i++) fHistNeventsPsi2s->GetXaxis()->SetBinLabel(i+1,CutNamePsi2s[i].Data());
332 fListHist->Add(fHistNeventsPsi2s);
334 fHistPsi2sMassVsPt = new TH2D("fHistPsi2sMassVsPt","Mass vs p_{T} of #psi(2s) candidates",100,3,6,50,0,5);
335 fHistPsi2sMassVsPt->GetXaxis()->SetTitle("Invariant mass(l^{+}l^{-}#pi^{+}#pi^{-}) (GeV/c)");
336 fHistPsi2sMassVsPt->GetYaxis()->SetTitle("#it{p}_{T} (GeV/#it{c})");
337 fListHist->Add(fHistPsi2sMassVsPt);
339 fHistPsi2sMassCoherent = new TH1D("fHistPsi2sMassAllCoherent","Invariant mass of coherent #psi(2s) candidates",100,3,6);
340 fHistPsi2sMassCoherent->GetXaxis()->SetTitle("Invariant mass(l^{+}l^{-}#pi^{+}#pi^{-}) (GeV/c)");
341 fListHist->Add(fHistPsi2sMassCoherent);
343 PostData(1, fJPsiTree);
344 PostData(2, fPsi2sTree);
345 PostData(3, fListTrig);
346 PostData(4, fListHist);
348 }//UserCreateOutputObjects
350 //_____________________________________________________________________________
351 void AliAnalysisTaskUpcPsi2s::UserExec(Option_t *)
354 //cout<<"#################### Next event ##################"<<endl;
358 if(fRunHist) RunESDhist();
359 if(fRunTree) RunESDtree();
364 if(fRunHist) RunAODhist();
365 if(fRunTree) RunAODtree();
369 //_____________________________________________________________________________
370 void AliAnalysisTaskUpcPsi2s::RunAODtrig()
374 AliAODEvent *aod = (AliAODEvent*) InputEvent();
377 fRunNum = aod ->GetRunNumber();
379 TString trigger = aod->GetFiredTriggerClasses();
381 if(trigger.Contains("CCUP4-B")) fHistCcup4TriggersPerRun->Fill(fRunNum); //CCUP4 triggers
382 if(trigger.Contains("CCUP7-B")) fHistCcup7TriggersPerRun->Fill(fRunNum); //CCUP7 triggers
383 if(trigger.Contains("CCUP2-B")) fHistCcup2TriggersPerRun->Fill(fRunNum); //CCUP2 triggers
385 if(trigger.Contains("CVLN_B2-B")) fHistCvlnTriggersPerRun->Fill(fRunNum); //CVLN triggers - synchronously downscaled
386 if(trigger.Contains("CVLN_R1-B")) fHistCvlnTriggersPerRun->Fill(fRunNum); //CVLN triggers - randomly downscaled
388 fL1inputs = aod->GetHeader()->GetL1TriggerInputs();
389 if(fL1inputs & (1 << 18)) fHistZedTriggersPerRun->Fill(fRunNum); //1ZED trigger inputs
391 //MB, Central and SemiCentral triggers
392 AliCentrality *centrality = aod->GetCentrality();
393 UInt_t selectionMask = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
395 Double_t percentile = centrality->GetCentralityPercentileUnchecked("V0M");
396 //Double_t percentile = centrality->GetCentralityPercentile("V0M");
398 if(((selectionMask & AliVEvent::kMB) == AliVEvent::kMB) && percentile<=80 && percentile>=0) fHistMBTriggersPerRun->Fill(fRunNum);
400 if(((selectionMask & AliVEvent::kCentral) == AliVEvent::kCentral) && percentile<=6 && percentile>=0 && (trigger.Contains("CVHN_R2-B"))) fHistCentralTriggersPerRun->Fill(fRunNum);
402 if(((selectionMask & AliVEvent::kSemiCentral) == AliVEvent::kSemiCentral) && percentile<=50 && percentile>=15) fHistSemiCentralTriggersPerRun->Fill(fRunNum);
404 PostData(3, fListTrig);
407 //_____________________________________________________________________________
408 void AliAnalysisTaskUpcPsi2s::RunAODhist()
411 TDatabasePDG *pdgdat = TDatabasePDG::Instance();
413 TParticlePDG *partMuon = pdgdat->GetParticle( 13 );
414 Double_t muonMass = partMuon->Mass();
416 TParticlePDG *partElectron = pdgdat->GetParticle( 11 );
417 Double_t electronMass = partElectron->Mass();
419 TParticlePDG *partPion = pdgdat->GetParticle( 211 );
420 Double_t pionMass = partPion->Mass();
423 AliAODEvent *aod = (AliAODEvent*) InputEvent();
426 fHistNeventsJPsi->Fill(1);
427 fHistNeventsPsi2s->Fill(1);
430 TString trigger = aod->GetFiredTriggerClasses();
432 if(!isMC && !trigger.Contains("CCUP") ) return;
434 fHistNeventsJPsi->Fill(2);
435 fHistNeventsPsi2s->Fill(2);
438 AliAODVertex *fAODVertex = aod->GetPrimaryVertex();
439 fVtxContrib = fAODVertex->GetNContributors();
440 if(fVtxContrib < 2) return;
442 fHistNeventsJPsi->Fill(3);
443 fHistNeventsPsi2s->Fill(3);
446 AliAODVZERO *fV0data = aod ->GetVZEROData();
447 AliAODZDC *fZDCdata = aod->GetZDCData();
449 fV0Adecision = fV0data->GetV0ADecision();
450 fV0Cdecision = fV0data->GetV0CDecision();
451 if(fV0Adecision != AliAODVZERO::kV0Empty || fV0Cdecision != AliAODVZERO::kV0Empty) return;
453 fHistNeventsJPsi->Fill(4);
454 fHistNeventsPsi2s->Fill(4);
456 fZDCAenergy = fZDCdata->GetZNATowerEnergy()[0];
457 fZDCCenergy = fZDCdata->GetZNCTowerEnergy()[0];
459 if( fZDCAenergy > 8200 || fZDCCenergy > 8200) return;
461 fHistNeventsJPsi->Fill(5);
462 fHistNeventsPsi2s->Fill(5);
465 Int_t nGoodTracks = 0;
466 Int_t TrackIndex[5] = {-1,-1,-1,-1,-1};
468 TLorentzVector vLepton[4], vPion[4], vCandidate, vDilepton;
469 Short_t qLepton[4], qPion[4];
470 UInt_t nLepton=0, nPion=0, nHighPt=0;
471 Double_t fRecTPCsignal[5];
472 Int_t mass[3]={-1,-1,-1};
476 for(Int_t itr=0; itr<aod ->GetNumberOfTracks(); itr++) {
477 AliAODTrack *trk = aod->GetTrack(itr);
479 if(!(trk->TestFilterBit(1<<0))) continue;
481 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
482 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
483 if(trk->GetTPCNcls() < 50)continue;
484 if(trk->Chi2perNDF() > 4)continue;
485 Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};
486 AliAODTrack* trk_clone=(AliAODTrack*)trk->Clone("trk_clone");
487 if(!trk_clone->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
490 if(TMath::Abs(dca[1]) > 2) continue;
492 TrackIndex[nGoodTracks] = itr;
495 if(nGoodTracks > 4) break;
498 nLepton=0; nPion=0; nHighPt=0;
499 mass[0]= -1; mass[1]= -1, mass[2]= -1;
501 if(nGoodTracks == 4){
502 fHistNeventsPsi2s->Fill(6);
503 for(Int_t i=0; i<4; i++){
504 AliAODTrack *trk = aod->GetTrack(TrackIndex[i]);
507 fRecTPCsignal[nLepton] = trk->GetTPCsignal();
508 qLepton[nLepton] = trk->Charge();
509 if(fRecTPCsignal[nLepton] > 40 && fRecTPCsignal[nLepton] < 70){
510 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), muonMass);
513 if(fRecTPCsignal[nLepton] > 70 && fRecTPCsignal[nLepton] < 100){
514 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), electronMass);
520 qPion[nPion] = trk->Charge();
521 vPion[nPion].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), pionMass);
525 if(nLepton > 2 || nPion > 2) break;
527 if((nLepton == 2) && (nPion == 2)){
528 fHistNeventsPsi2s->Fill(7);
529 if(qLepton[0]*qLepton[1] > 0) fHistNeventsPsi2s->Fill(8);
530 if(qPion[0]*qPion[1] > 0) fHistNeventsPsi2s->Fill(9);
531 if((qLepton[0]*qLepton[1] > 0) && (qPion[0]*qPion[1] > 0)) fHistNeventsPsi2s->Fill(10);
532 if((qLepton[0]*qLepton[1] < 0) && (qPion[0]*qPion[1] < 0)){
533 fHistNeventsPsi2s->Fill(11);
534 if(mass[0] == mass[1]) {
535 fHistNeventsPsi2s->Fill(12);
536 vCandidate = vLepton[0]+vLepton[1]+vPion[0]+vPion[1];
537 vDilepton = vLepton[0]+vLepton[1];
538 fHistPsi2sMassVsPt->Fill(vCandidate.M(),vCandidate.Pt());
539 if(vCandidate.Pt() < 0.15) fHistPsi2sMassCoherent->Fill(vCandidate.M());
540 if(mass[0] == 0) fHistNeventsPsi2s->Fill(13);
541 if(mass[0] == 1) fHistNeventsPsi2s->Fill(14);
549 for(Int_t itr=0; itr<aod ->GetNumberOfTracks(); itr++) {
550 AliAODTrack *trk = aod->GetTrack(itr);
552 if(!(trk->TestFilterBit(1<<0))) continue;
554 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
555 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
556 if(trk->GetTPCNcls() < 70)continue;
557 if(trk->Chi2perNDF() > 4)continue;
558 if((!trk->HasPointOnITSLayer(0))&&(!trk->HasPointOnITSLayer(1))) continue;
559 Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};
560 AliAODTrack* trk_clone=(AliAODTrack*)trk->Clone("trk_clone");
561 if(!trk_clone->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
563 if(TMath::Abs(dca[1]) > 2) continue;
564 if(TMath::Abs(dca[0]) > 0.2) continue;
566 TrackIndex[nGoodTracks] = itr;
569 if(nGoodTracks > 2) break;
572 nLepton=0; nPion=0; nHighPt=0;
573 mass[0]= -1; mass[1]= -1, mass[2]= -1;
575 if(nGoodTracks == 2){
576 fHistNeventsJPsi->Fill(6);
577 for(Int_t i=0; i<2; i++){
578 AliAODTrack *trk = aod->GetTrack(TrackIndex[i]);
579 if(trk->Pt() > 1) nHighPt++;
580 fRecTPCsignal[nLepton] = trk->GetTPCsignal();
581 qLepton[nLepton] = trk->Charge();
582 if(fRecTPCsignal[nLepton] > 40 && fRecTPCsignal[nLepton] < 70){
583 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), muonMass);
586 if(fRecTPCsignal[nLepton] > 70 && fRecTPCsignal[nLepton] < 100){
587 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), electronMass);
593 if(qLepton[0]*qLepton[1] > 0) fHistNeventsJPsi->Fill(7);
594 if(qLepton[0]*qLepton[1] < 0){
595 fHistNeventsJPsi->Fill(8);
597 fHistNeventsJPsi->Fill(9);
598 fHistTPCsignalJPsi->Fill(fRecTPCsignal[0],fRecTPCsignal[1]);
599 if(nHighPt == 2) fHistNeventsJPsi->Fill(10);
600 if(mass[0] == mass[1] && mass[0] != -1) {
601 fHistNeventsJPsi->Fill(11);
602 vCandidate = vLepton[0]+vLepton[1];
603 if( vCandidate.M() > 2.8 && vCandidate.M() < 3.2) fHistDiLeptonPtJPsi->Fill(vLepton[0].Pt(),vLepton[1].Pt());
605 fHistDiMuonMass->Fill(vCandidate.M());
606 fHistNeventsJPsi->Fill(12);
609 fHistDiElectronMass->Fill(vCandidate.M());
610 fHistNeventsJPsi->Fill(13);
619 PostData(4, fListHist);
623 //_____________________________________________________________________________
624 void AliAnalysisTaskUpcPsi2s::RunAODtree()
627 AliAODEvent *aod = (AliAODEvent*) InputEvent();
630 if(isMC) RunAODMC(aod);
633 const char *filnam = ((TTree*) GetInputData(0))->GetCurrentFile()->GetName();
634 fDataFilnam->Clear();
635 fDataFilnam->SetString(filnam);
636 fEvtNum = ((TTree*) GetInputData(0))->GetTree()->GetReadEntry();
637 fRunNum = aod ->GetRunNumber();
640 TString trigger = aod->GetFiredTriggerClasses();
642 fTrigger[0] = trigger.Contains("CCUP4-B"); // Central UPC Pb-Pb 2011
643 fTrigger[1] = trigger.Contains("CCUP2-B"); // Double gap
644 fTrigger[2] = trigger.Contains("CCUP7-B"); // Central UPC p-Pb 2013
646 Bool_t isTriggered = kFALSE;
647 for(Int_t i=0; i<ntrg; i++) {
648 if( fTrigger[i] ) isTriggered = kTRUE;
650 if(!isMC && !isTriggered ) return;
653 fL0inputs = aod->GetHeader()->GetL0TriggerInputs();
654 fL1inputs = aod->GetHeader()->GetL1TriggerInputs();
656 //Event identification
657 fPerNum = aod ->GetPeriodNumber();
658 fOrbNum = aod ->GetOrbitNumber();
659 fBCrossNum = aod ->GetBunchCrossNumber();
662 AliAODVertex *fAODVertex = aod->GetPrimaryVertex();
663 fVtxContrib = fAODVertex->GetNContributors();
664 fVtxPos[0] = fAODVertex->GetX();
665 fVtxPos[1] = fAODVertex->GetY();
666 fVtxPos[2] = fAODVertex->GetZ();
668 fAODVertex->GetCovarianceMatrix(CovMatx);
669 fVtxErr[0] = CovMatx[0];
670 fVtxErr[1] = CovMatx[1];
671 fVtxErr[2] = CovMatx[2];
672 fVtxChi2 = fAODVertex->GetChi2();
673 fVtxNDF = fAODVertex->GetNDF();
676 fNtracklets = aod->GetTracklets()->GetNumberOfTracklets();
679 AliAODVZERO *fV0data = aod ->GetVZEROData();
680 AliAODZDC *fZDCdata = aod->GetZDCData();
682 fV0Adecision = fV0data->GetV0ADecision();
683 fV0Cdecision = fV0data->GetV0CDecision();
684 fZDCAenergy = fZDCdata->GetZNATowerEnergy()[0];
685 fZDCCenergy = fZDCdata->GetZNCTowerEnergy()[0];
689 //Track loop - loose cuts
690 for(Int_t itr=0; itr<aod ->GetNumberOfTracks(); itr++) {
691 AliAODTrack *trk = aod->GetTrack(itr);
693 if(!(trk->TestFilterBit(1<<0))) continue;
695 if(!(trk->GetStatus() & AliAODTrack::kTPCrefit) ) continue;
696 if(!(trk->GetStatus() & AliAODTrack::kITSrefit) ) continue;
697 if(trk->GetTPCNcls() < 20)continue;
699 }//Track loop -loose cuts
702 Int_t TrackIndex[5] = {-1,-1,-1,-1,-1};
705 for(Int_t itr=0; itr<aod ->GetNumberOfTracks(); itr++) {
706 AliAODTrack *trk = aod->GetTrack(itr);
708 if(!(trk->TestFilterBit(1<<0))) continue;
710 if(!(trk->GetStatus() & AliAODTrack::kTPCrefit) ) continue;
711 if(!(trk->GetStatus() & AliAODTrack::kITSrefit) ) continue;
712 if(trk->GetTPCNcls() < 70)continue;
713 if(trk->Chi2perNDF() > 4)continue;
714 if((!trk->HasPointOnITSLayer(0))&&(!trk->HasPointOnITSLayer(1))) continue;
715 Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};
716 AliAODTrack* trk_clone=(AliAODTrack*)trk->Clone("trk_clone");
717 if(!trk_clone->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
719 if(TMath::Abs(dca[1]) > 2) continue;
720 if(TMath::Abs(dca[0]) > 0.2) continue;
722 TrackIndex[nGoodTracks] = itr;
725 if(nGoodTracks > 2) break;
728 fJPsiAODTracks->Clear("C");
729 if(nGoodTracks == 2){
731 TDatabasePDG *pdgdat = TDatabasePDG::Instance();
732 TParticlePDG *partMuon = pdgdat->GetParticle( 13 );
733 Double_t muonMass = partMuon->Mass();
734 TParticlePDG *partElectron = pdgdat->GetParticle( 11 );
735 Double_t electronMass = partElectron->Mass();
736 TParticlePDG *partPion = pdgdat->GetParticle( 211 );
737 Double_t pionMass = partPion->Mass();
741 Double_t KFmass = pionMass;
742 Double_t fRecTPCsignal;
743 AliKFParticle *KFpart[2];
744 AliKFVertex *KFvtx = new AliKFVertex();
745 KFvtx->SetField(aod->GetMagneticField());
747 for(Int_t i=0; i<2; i++){
748 AliAODTrack *trk = aod->GetTrack(TrackIndex[i]);
750 Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};
751 AliAODTrack* trk_clone=(AliAODTrack*)trk->Clone("trk_clone");
752 if(!trk_clone->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
755 new((*fJPsiAODTracks)[i]) AliAODTrack(*trk);
756 ((AliAODTrack*)((*fJPsiAODTracks)[i]))->SetDCA(dca[0],dca[1]);//to get DCAxy trk->DCA(); to get DCAz trk->ZAtDCA();
758 fPIDMuon[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kMuon);
759 fPIDElectron[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kElectron);
760 fPIDPion[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kPion);
762 trk->GetPosition(KFpar);
763 trk->PxPyPz(KFpar+3);
764 trk->GetCovMatrix(KFcov);
767 fRecTPCsignal = trk->GetTPCsignal();
768 if(fRecTPCsignal > 40 && fRecTPCsignal < 70) KFmass = muonMass;
769 if(fRecTPCsignal > 70 && fRecTPCsignal < 100)KFmass = electronMass;
771 else KFmass = pionMass;
773 KFpart[i] = new AliKFParticle();
774 KFpart[i]->SetField(aod->GetMagneticField());
775 KFpart[i]->AliKFParticleBase::Initialize(KFpar,KFcov,(Int_t) trk->Charge(), KFmass);
776 KFvtx->AddDaughter(*KFpart[i]);
779 Double_t pos[3]={0,0,0};
780 AliExternalTrackParam *parTrk = new AliExternalTrackParam();
781 parTrk->CopyFromVTrack((AliVTrack*) trk);
782 if(!parTrk->GetXYZAt(378,aod->GetMagneticField(),pos)) fTOFphi[i] = -666;
784 fTOFphi[i] = TMath::ATan2(pos[1],pos[0])*TMath::RadToDeg();
785 if(fTOFphi[i] < 0) fTOFphi[i]+=(2*TMath::Pi()*TMath::RadToDeg());
789 fKfVtxPos[0]= KFvtx->GetX();
790 fKfVtxPos[1]= KFvtx->GetY();
791 fKfVtxPos[2]= KFvtx->GetZ();
792 for(UInt_t i=0; i<2; i++)delete KFpart[i];
795 if(!isMC) fJPsiTree ->Fill();
800 for(Int_t itr=0; itr<aod ->GetNumberOfTracks(); itr++) {
801 AliAODTrack *trk = aod->GetTrack(itr);
803 if(!(trk->TestFilterBit(1<<0))) continue;
805 if(!(trk->GetStatus() & AliAODTrack::kTPCrefit) ) continue;
806 if(!(trk->GetStatus() & AliAODTrack::kITSrefit) ) continue;
807 if(trk->GetTPCNcls() < 50)continue;
808 if(trk->Chi2perNDF() > 4)continue;
809 Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};
810 AliAODTrack* trk_clone=(AliAODTrack*)trk->Clone("trk_clone");
811 if(!trk_clone->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
813 if(!trk->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
814 if(TMath::Abs(dca[1]) > 2) continue;
816 TrackIndex[nGoodTracks] = itr;
819 if(nGoodTracks > 4) break;
822 fPsi2sAODTracks->Clear("C");
823 if(nGoodTracks == 4){
825 TDatabasePDG *pdgdat = TDatabasePDG::Instance();
826 TParticlePDG *partMuon = pdgdat->GetParticle( 13 );
827 Double_t muonMass = partMuon->Mass();
828 TParticlePDG *partElectron = pdgdat->GetParticle( 11 );
829 Double_t electronMass = partElectron->Mass();
830 TParticlePDG *partPion = pdgdat->GetParticle( 211 );
831 Double_t pionMass = partPion->Mass();
835 Double_t KFmass = pionMass;
836 Double_t fRecTPCsignal;
837 AliKFParticle *KFpart[4];
838 AliKFVertex *KFvtx = new AliKFVertex();
839 KFvtx->SetField(aod->GetMagneticField());
841 for(Int_t i=0; i<4; i++){
842 AliAODTrack *trk = aod->GetTrack(TrackIndex[i]);
844 Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};
845 AliAODTrack* trk_clone=(AliAODTrack*)trk->Clone("trk_clone");
846 if(!trk_clone->PropagateToDCA(fAODVertex,aod->GetMagneticField(),300.,dca,cov)) continue;
849 new((*fPsi2sAODTracks)[i]) AliAODTrack(*trk);
850 ((AliAODTrack*)((*fPsi2sAODTracks)[i]))->SetDCA(dca[0],dca[1]);//to get DCAxy trk->DCA(); to get DCAz trk->ZAtDCA();
853 fPIDMuon[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kMuon);
854 fPIDElectron[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kElectron);
855 fPIDPion[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kPion);
857 trk->GetPosition(KFpar);
858 trk->PxPyPz(KFpar+3);
859 trk->GetCovMatrix(KFcov);
862 fRecTPCsignal = trk->GetTPCsignal();
863 if(fRecTPCsignal > 40 && fRecTPCsignal < 70) KFmass = muonMass;
864 if(fRecTPCsignal > 70 && fRecTPCsignal < 100)KFmass = electronMass;
866 else KFmass = pionMass;
868 KFpart[i] = new AliKFParticle();
869 KFpart[i]->SetField(aod->GetMagneticField());
870 KFpart[i]->AliKFParticleBase::Initialize(KFpar,KFcov,(Int_t) trk->Charge(), KFmass);
871 KFvtx->AddDaughter(*KFpart[i]);
873 Double_t pos[3]={0,0,0};
874 AliExternalTrackParam *parTrk = new AliExternalTrackParam();
875 parTrk->CopyFromVTrack((AliVTrack*) trk);
876 if(!parTrk->GetXYZAt(378,aod->GetMagneticField(),pos)) fTOFphi[i] = -666;
878 fTOFphi[i] = TMath::ATan2(pos[1],pos[0])*TMath::RadToDeg();
879 if(fTOFphi[i] < 0) fTOFphi[i]+=(2*TMath::Pi()*TMath::RadToDeg());
883 fKfVtxPos[0]= KFvtx->GetX();
884 fKfVtxPos[1]= KFvtx->GetY();
885 fKfVtxPos[2]= KFvtx->GetZ();
886 for(UInt_t i=0; i<4; i++)delete KFpart[i];
888 if(!isMC) fPsi2sTree ->Fill();
896 PostData(1, fJPsiTree);
897 PostData(2, fPsi2sTree);
902 //_____________________________________________________________________________
903 void AliAnalysisTaskUpcPsi2s::RunAODMC(AliAODEvent *aod)
906 fGenPart->Clear("C");
908 TClonesArray *arrayMC = (TClonesArray*) aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
912 //loop over mc particles
913 for(Int_t imc=0; imc<arrayMC->GetEntriesFast(); imc++) {
914 AliAODMCParticle *mcPart = (AliAODMCParticle*) arrayMC->At(imc);
915 if(!mcPart) continue;
917 if(mcPart->GetMother() >= 0) continue;
919 TParticle *part = (TParticle*) fGenPart->ConstructedAt(nmc++);
920 part->SetMomentum(mcPart->Px(), mcPart->Py(), mcPart->Pz(), mcPart->E());
921 part->SetPdgCode(mcPart->GetPdgCode());
922 part->SetUniqueID(imc);
923 }//loop over mc particles
928 //_____________________________________________________________________________
929 void AliAnalysisTaskUpcPsi2s::RunESDtrig()
933 AliESDEvent *esd = (AliESDEvent*) InputEvent();
936 fRunNum = esd ->GetRunNumber();
938 TString trigger = esd->GetFiredTriggerClasses();
940 if(trigger.Contains("CCUP4-B")) fHistCcup4TriggersPerRun->Fill(fRunNum); //CCUP4 triggers
941 if(trigger.Contains("CCUP7-B")) fHistCcup7TriggersPerRun->Fill(fRunNum); //CCUP7 triggers
942 if(trigger.Contains("CCUP2-B")) fHistCcup2TriggersPerRun->Fill(fRunNum); //CCUP2 triggers
944 if(trigger.Contains("CVLN_B2-B")) fHistCvlnTriggersPerRun->Fill(fRunNum); //CVLN triggers - synchronously downscaled
945 if(trigger.Contains("CVLN_R1-B")) fHistCvlnTriggersPerRun->Fill(fRunNum); //CVLN triggers - randomly downscaled
947 if(esd->GetHeader()->IsTriggerInputFired("1ZED")) fHistZedTriggersPerRun->Fill(fRunNum); //1ZED trigger inputs
949 //MB, Central and SemiCentral triggers
950 AliCentrality *centrality = esd->GetCentrality();
951 UInt_t selectionMask = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
953 //Double_t percentile = centrality->GetCentralityPercentile("V0M");
954 Double_t percentile = centrality->GetCentralityPercentileUnchecked("V0M");
956 if(((selectionMask & AliVEvent::kMB) == AliVEvent::kMB) && percentile<=80 && percentile>=0) fHistMBTriggersPerRun->Fill(fRunNum);
958 if(((selectionMask & AliVEvent::kCentral) == AliVEvent::kCentral) && percentile<=6 && percentile>=0 && (trigger.Contains("CVHN_R2-B"))) fHistCentralTriggersPerRun->Fill(fRunNum);
960 if(((selectionMask & AliVEvent::kSemiCentral) == AliVEvent::kSemiCentral) && percentile<=50 && percentile>=15) fHistSemiCentralTriggersPerRun->Fill(fRunNum);
963 PostData(3, fListTrig);
966 //_____________________________________________________________________________
967 void AliAnalysisTaskUpcPsi2s::RunESDhist()
970 TDatabasePDG *pdgdat = TDatabasePDG::Instance();
972 TParticlePDG *partMuon = pdgdat->GetParticle( 13 );
973 Double_t muonMass = partMuon->Mass();
975 TParticlePDG *partElectron = pdgdat->GetParticle( 11 );
976 Double_t electronMass = partElectron->Mass();
978 TParticlePDG *partPion = pdgdat->GetParticle( 211 );
979 Double_t pionMass = partPion->Mass();
982 AliESDEvent *esd = (AliESDEvent*) InputEvent();
985 fHistNeventsJPsi->Fill(1);
986 fHistNeventsPsi2s->Fill(1);
989 TString trigger = esd->GetFiredTriggerClasses();
991 if(!isMC && !trigger.Contains("CCUP") ) return;
993 fHistNeventsJPsi->Fill(2);
994 fHistNeventsPsi2s->Fill(2);
997 AliESDVertex *fESDVertex = (AliESDVertex*) esd->GetPrimaryVertex();
998 fVtxContrib = fESDVertex->GetNContributors();
999 if(fVtxContrib < 2) return;
1001 fHistNeventsJPsi->Fill(3);
1002 fHistNeventsPsi2s->Fill(3);
1005 AliESDVZERO *fV0data = esd->GetVZEROData();
1006 AliESDZDC *fZDCdata = esd->GetESDZDC();
1008 fV0Adecision = fV0data->GetV0ADecision();
1009 fV0Cdecision = fV0data->GetV0CDecision();
1010 if(fV0Adecision != AliESDVZERO::kV0Empty || fV0Cdecision != AliESDVZERO::kV0Empty) return;
1012 fZDCAenergy = fZDCdata->GetZN2TowerEnergy()[0];
1013 fZDCCenergy = fZDCdata->GetZN1TowerEnergy()[0];
1014 if( fZDCAenergy > 8200 || fZDCCenergy > 8200) return;
1016 fHistNeventsJPsi->Fill(4);
1017 fHistNeventsPsi2s->Fill(4);
1019 Int_t nGoodTracks=0;
1021 Int_t TrackIndex[5] = {-1,-1,-1,-1,-1};
1023 TLorentzVector vLepton[4], vPion[4], vCandidate, vDilepton;
1024 Short_t qLepton[4], qPion[4];
1025 UInt_t nLepton=0, nPion=0, nHighPt=0;
1026 Double_t fRecTPCsignal[5];
1027 Int_t mass[3]={-1,-1,-1};
1030 for(Int_t itr=0; itr<esd ->GetNumberOfTracks(); itr++) {
1031 AliESDtrack *trk = esd->GetTrack(itr);
1032 if( !trk ) continue;
1034 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
1035 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
1036 if(trk->GetTPCNcls() < 70)continue;
1037 if(trk->GetTPCchi2()/trk->GetTPCNcls() > 4)continue;
1038 if((!trk->HasPointOnITSLayer(0))&&(!trk->HasPointOnITSLayer(1))) continue;
1039 Float_t dca[2] = {0.0,0.0}; AliExternalTrackParam cParam;
1040 if(!trk->RelateToVertex(fESDVertex, esd->GetMagneticField(),300.,&cParam)) continue;
1041 trk->GetImpactParameters(dca[0],dca[1]);
1042 if(TMath::Abs(dca[1]) > 2) continue;
1043 if(TMath::Abs(dca[1]) > 0.2) continue;
1045 TrackIndex[nGoodTracks] = itr;
1047 if(nGoodTracks > 2) break;
1051 if(nGoodTracks == 2){
1052 fHistNeventsJPsi->Fill(5);
1053 for(Int_t i=0; i<2; i++){
1054 AliESDtrack *trk = esd->GetTrack(TrackIndex[i]);
1055 if(trk->Pt() > 1) nHighPt++;
1056 fRecTPCsignal[nLepton] = trk->GetTPCsignal();
1057 qLepton[nLepton] = trk->Charge();
1058 if(fRecTPCsignal[nLepton] > 40 && fRecTPCsignal[nLepton] < 70){
1059 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), muonMass);
1062 if(fRecTPCsignal[nLepton] > 70 && fRecTPCsignal[nLepton] < 100){
1063 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), electronMass);
1069 if(qLepton[0]*qLepton[1] > 0) fHistNeventsJPsi->Fill(6);
1070 if(qLepton[0]*qLepton[1] < 0){
1071 fHistNeventsJPsi->Fill(7);
1073 fHistNeventsJPsi->Fill(8);
1074 fHistTPCsignalJPsi->Fill(fRecTPCsignal[0],fRecTPCsignal[1]);
1075 if(nHighPt == 2) fHistNeventsJPsi->Fill(9);
1076 if(mass[0] == mass[1] && mass[0] != -1) {
1077 fHistNeventsJPsi->Fill(10);
1078 vCandidate = vLepton[0]+vLepton[1];
1079 if( vCandidate.M() > 2.8 && vCandidate.M() < 3.2) fHistDiLeptonPtJPsi->Fill(vLepton[0].Pt(),vLepton[1].Pt());
1081 fHistDiMuonMass->Fill(vCandidate.M());
1082 fHistNeventsJPsi->Fill(11);
1085 fHistDiElectronMass->Fill(vCandidate.M());
1086 fHistNeventsJPsi->Fill(12);
1093 nGoodTracks = 0; nLepton=0; nPion=0; nHighPt=0;
1094 mass[0]= -1; mass[1]= -1, mass[2]= -1;
1097 for(Int_t itr=0; itr<esd ->GetNumberOfTracks(); itr++) {
1098 AliESDtrack *trk = esd->GetTrack(itr);
1099 if( !trk ) continue;
1101 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
1102 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
1103 if(trk->GetTPCNcls() < 50)continue;
1104 if(trk->GetTPCchi2()/trk->GetTPCNcls() > 4)continue;
1105 Float_t dca[2] = {0.0,0.0}; AliExternalTrackParam cParam;
1106 if(!trk->RelateToVertex(fESDVertex, esd->GetMagneticField(),300.,&cParam)) continue;
1107 trk->GetImpactParameters(dca[0],dca[1]);
1108 if(TMath::Abs(dca[1]) > 2) continue;
1110 TrackIndex[nGoodTracks] = itr;
1112 if(nGoodTracks > 4) break;
1115 if(nGoodTracks == 4){
1116 fHistNeventsPsi2s->Fill(5);
1117 for(Int_t i=0; i<4; i++){
1118 AliESDtrack *trk = esd->GetTrack(TrackIndex[i]);
1121 fRecTPCsignal[nLepton] = trk->GetTPCsignal();
1122 qLepton[nLepton] = trk->Charge();
1123 if(fRecTPCsignal[nLepton] > 40 && fRecTPCsignal[nLepton] < 70){
1124 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), muonMass);
1127 if(fRecTPCsignal[nLepton] > 70 && fRecTPCsignal[nLepton] < 100){
1128 vLepton[nLepton].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), electronMass);
1134 qPion[nPion] = trk->Charge();
1135 vPion[nPion].SetPtEtaPhiM(trk->Pt(), trk->Eta(), trk->Phi(), pionMass);
1139 if(nLepton > 2 || nPion > 2) break;
1141 if((nLepton == 2) && (nPion == 2)){
1142 fHistNeventsPsi2s->Fill(6);
1143 if(qLepton[0]*qLepton[1] > 0) fHistNeventsPsi2s->Fill(7);
1144 if(qPion[0]*qPion[1] > 0) fHistNeventsPsi2s->Fill(8);
1145 if((qLepton[0]*qLepton[1] > 0) && (qPion[0]*qPion[1] > 0)) fHistNeventsPsi2s->Fill(9);
1146 if((qLepton[0]*qLepton[1] < 0) && (qPion[0]*qPion[1] < 0)){
1147 fHistNeventsPsi2s->Fill(10);
1148 if(mass[0] == mass[1]) {
1149 fHistNeventsPsi2s->Fill(11);
1150 vCandidate = vLepton[0]+vLepton[1]+vPion[0]+vPion[1];
1151 vDilepton = vLepton[0]+vLepton[1];
1152 fHistPsi2sMassVsPt->Fill(vCandidate.M(),vCandidate.Pt());
1153 if(vCandidate.Pt() < 0.15) fHistPsi2sMassCoherent->Fill(vCandidate.M());
1154 if(mass[0] == 0) fHistNeventsPsi2s->Fill(12);
1155 if(mass[0] == 1) fHistNeventsPsi2s->Fill(13);
1161 PostData(4, fListHist);
1165 //_____________________________________________________________________________
1166 void AliAnalysisTaskUpcPsi2s::RunESDtree()
1170 AliESDEvent *esd = (AliESDEvent*) InputEvent();
1173 if(isMC) RunESDMC(esd);
1176 const char *filnam = ((TTree*) GetInputData(0))->GetCurrentFile()->GetName();
1177 fDataFilnam->Clear();
1178 fDataFilnam->SetString(filnam);
1179 fEvtNum = ((TTree*) GetInputData(0))->GetTree()->GetReadEntry();
1180 fRunNum = esd->GetRunNumber();
1183 TString trigger = esd->GetFiredTriggerClasses();
1185 fTrigger[0] = trigger.Contains("CCUP4-B"); // Central UPC Pb-Pb 2011
1186 fTrigger[1] = trigger.Contains("CCUP2-B"); // Double gap
1187 fTrigger[2] = trigger.Contains("CCUP7-B"); // Central UPC p-Pb 2013
1189 Bool_t isTriggered = kFALSE;
1190 for(Int_t i=0; i<ntrg; i++) {
1191 if( fTrigger[i] ) isTriggered = kTRUE;
1193 if(!isMC && !isTriggered ) return;
1196 fL0inputs = esd->GetHeader()->GetL0TriggerInputs();
1197 fL1inputs = esd->GetHeader()->GetL1TriggerInputs();
1199 //TOF trigger info (0OMU)
1200 fTOFtrig1 = esd->GetHeader()->IsTriggerInputFired("0OMU");
1201 fTOFtrig2 = esd->GetHeader()->GetActiveTriggerInputs().Contains("0OMU") ? ((esd->GetHeader()) ? esd->GetHeader()->IsTriggerInputFired("0OMU") : kFALSE) : TESTBIT(esd->GetHeader()->GetL0TriggerInputs(), (kFALSE) ? 21 : 9);
1203 //Event identification
1204 fPerNum = esd->GetPeriodNumber();
1205 fOrbNum = esd->GetOrbitNumber();
1206 fBCrossNum = esd->GetBunchCrossNumber();
1209 AliESDVertex *fESDVertex = (AliESDVertex*) esd->GetPrimaryVertex();
1210 fVtxContrib = fESDVertex->GetNContributors();
1211 fVtxPos[0] = fESDVertex->GetX();
1212 fVtxPos[1] = fESDVertex->GetY();
1213 fVtxPos[2] = fESDVertex->GetZ();
1214 Double_t CovMatx[6];
1215 fESDVertex->GetCovarianceMatrix(CovMatx);
1216 fVtxErr[0] = CovMatx[0];
1217 fVtxErr[1] = CovMatx[1];
1218 fVtxErr[2] = CovMatx[2];
1219 fVtxChi2 = fESDVertex->GetChi2();
1220 fVtxNDF = fESDVertex->GetNDF();
1223 fNtracklets = esd->GetMultiplicity()->GetNumberOfTracklets();
1226 AliESDVZERO *fV0data = esd->GetVZEROData();
1227 AliESDZDC *fZDCdata = esd->GetESDZDC();
1229 fV0Adecision = fV0data->GetV0ADecision();
1230 fV0Cdecision = fV0data->GetV0CDecision();
1231 fZDCAenergy = fZDCdata->GetZN2TowerEnergy()[0];
1232 fZDCCenergy = fZDCdata->GetZN1TowerEnergy()[0];
1236 //Track loop - loose cuts
1237 for(Int_t itr=0; itr<esd ->GetNumberOfTracks(); itr++) {
1238 AliESDtrack *trk = esd->GetTrack(itr);
1239 if( !trk ) continue;
1241 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
1242 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
1243 if(trk->GetTPCNcls() < 20)continue;
1245 }//Track loop -loose cuts
1247 Int_t nGoodTracks=0;
1248 Int_t TrackIndex[5] = {-1,-1,-1,-1,-1};
1251 for(Int_t itr=0; itr<esd ->GetNumberOfTracks(); itr++) {
1252 AliESDtrack *trk = esd->GetTrack(itr);
1253 if( !trk ) continue;
1255 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
1256 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
1257 if(trk->GetTPCNcls() < 70)continue;
1258 if(trk->GetTPCchi2()/trk->GetTPCNcls() > 4)continue;
1259 if((!trk->HasPointOnITSLayer(0))&&(!trk->HasPointOnITSLayer(1))) continue;
1260 Float_t dca[2] = {0.0,0.0}; AliExternalTrackParam cParam;
1261 if(!trk->RelateToVertex(fESDVertex, esd->GetMagneticField(),300.,&cParam)) continue;
1262 trk->GetImpactParameters(dca[0],dca[1]);
1263 if(TMath::Abs(dca[1]) > 2) continue;
1264 if(TMath::Abs(dca[1]) > 0.2) continue;
1266 TrackIndex[nGoodTracks] = itr;
1268 if(nGoodTracks > 2) break;
1271 fJPsiESDTracks->Clear("C");
1272 if(nGoodTracks == 2){
1273 for(Int_t i=0; i<2; i++){
1274 AliESDtrack *trk = esd->GetTrack(TrackIndex[i]);
1276 AliExternalTrackParam cParam;
1277 trk->RelateToVertex(fESDVertex, esd->GetMagneticField(),300.,&cParam);// to get trk->GetImpactParameters(DCAxy,DCAz);
1279 new((*fJPsiESDTracks)[i]) AliESDtrack(*trk);
1281 fPIDMuon[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kMuon);
1282 fPIDElectron[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kElectron);
1283 fPIDPion[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kPion);
1285 Double_t pos[3]={0,0,0};
1286 if(!trk->GetXYZAt(378,esd->GetMagneticField(),pos)) fTOFphi[i] = -666;
1288 fTOFphi[i] = TMath::ATan2(pos[1],pos[0])*TMath::RadToDeg();
1289 if(fTOFphi[i] < 0) fTOFphi[i]+=(2*TMath::Pi()*TMath::RadToDeg());
1292 if(!isMC) fJPsiTree ->Fill();
1297 for(Int_t itr=0; itr<esd ->GetNumberOfTracks(); itr++) {
1298 AliESDtrack *trk = esd->GetTrack(itr);
1299 if( !trk ) continue;
1301 if(!(trk->GetStatus() & AliESDtrack::kTPCrefit) ) continue;
1302 if(!(trk->GetStatus() & AliESDtrack::kITSrefit) ) continue;
1303 if(trk->GetTPCNcls() < 50)continue;
1304 if(trk->GetTPCchi2()/trk->GetTPCNcls() > 4)continue;
1305 Float_t dca[2] = {0.0,0.0}; AliExternalTrackParam cParam;
1306 if(!trk->RelateToVertex(fESDVertex, esd->GetMagneticField(),300.,&cParam)) continue;
1307 trk->GetImpactParameters(dca[0],dca[1]);
1308 if(TMath::Abs(dca[1]) > 2) continue;
1310 TrackIndex[nGoodTracks] = itr;
1312 if(nGoodTracks > 4) break;
1315 fPsi2sESDTracks->Clear("C");
1316 if(nGoodTracks == 4){
1317 for(Int_t i=0; i<4; i++){
1318 AliESDtrack *trk = esd->GetTrack(TrackIndex[i]);
1320 AliExternalTrackParam cParam;
1321 trk->RelateToVertex(fESDVertex, esd->GetMagneticField(),300.,&cParam);// to get trk->GetImpactParameters(DCAxy,DCAz);
1323 new((*fPsi2sESDTracks)[i]) AliESDtrack(*trk);
1325 fPIDMuon[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kMuon);
1326 fPIDElectron[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kElectron);
1327 fPIDPion[i] = fPIDResponse->NumberOfSigmasTPC(trk,AliPID::kPion);
1329 Double_t pos[3]={0,0,0};
1330 if(!trk->GetXYZAt(378,esd->GetMagneticField(),pos)) fTOFphi[i] = -666;
1332 fTOFphi[i] = TMath::ATan2(pos[1],pos[0])*TMath::RadToDeg();
1333 if(fTOFphi[i] < 0) fTOFphi[i]+=(2*TMath::Pi()*TMath::RadToDeg());
1336 if(!isMC) fPsi2sTree ->Fill();
1341 fPsi2sTree ->Fill();
1344 PostData(1, fJPsiTree);
1345 PostData(2, fPsi2sTree);
1350 //_____________________________________________________________________________
1351 void AliAnalysisTaskUpcPsi2s::RunESDMC(AliESDEvent* esd)
1354 AliTriggerAnalysis *fTrigAna = new AliTriggerAnalysis();
1355 fTrigAna->SetAnalyzeMC(isMC);
1357 if(fTrigAna->SPDFiredChips(esd,1,kFALSE,2) > 1) fTriggerInputsMC[0] = kTRUE;
1358 if(fTrigAna->SPDFiredChips(esd,1,kFALSE,2) < 2) fTriggerInputsMC[0] = kFALSE;
1360 fTriggerInputsMC[1] = esd->GetHeader()->IsTriggerInputFired("0OMU");
1361 fTriggerInputsMC[2] = esd->GetHeader()->IsTriggerInputFired("0VBA");
1362 fTriggerInputsMC[3] = esd->GetHeader()->IsTriggerInputFired("0VBC");
1364 fGenPart->Clear("C");
1366 AliMCEvent *mc = MCEvent();
1370 //loop over mc particles
1371 for(Int_t imc=0; imc<mc->GetNumberOfTracks(); imc++) {
1372 AliMCParticle *mcPart = (AliMCParticle*) mc->GetTrack(imc);
1373 if(!mcPart) continue;
1375 if(mcPart->GetMother() >= 0) continue;
1377 TParticle *part = (TParticle*) fGenPart->ConstructedAt(nmc++);
1378 part->SetMomentum(mcPart->Px(), mcPart->Py(), mcPart->Pz(), mcPart->E());
1379 part->SetPdgCode(mcPart->PdgCode());
1380 part->SetUniqueID(imc);
1381 }//loop over mc particles
1387 //_____________________________________________________________________________
1388 void AliAnalysisTaskUpcPsi2s::Terminate(Option_t *)
1391 cout<<"Analysis complete."<<endl;