1 /**************************************************************************
2 * Copyright(c) 1998-2009, 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 **************************************************************************/
16 ///////////////////////////////////////////////////////////////////////////
18 // AliAnalysisTaskSE for reading both reconstructed JPSI -> ee candidates
19 // and like sign pairs and for drawing corresponding distributions
21 // Author: C.Di Giglio, carmelo.digiglio@ba.infn.it
22 ///////////////////////////////////////////////////////////////////////////
26 #include <TClonesArray.h>
31 #include "AliAnalysisManager.h"
32 #include "AliAODHandler.h"
33 #include "AliAODEvent.h"
34 #include "AliAODVertex.h"
35 #include "AliAODTrack.h"
36 #include "AliAODRecoDecayHF2Prong.h"
37 #include "AliAnalysisVertexingHF.h"
38 #include "AliAnalysisTaskSE.h"
39 #include "AliAnalysisTaskSEBkgLikeSignJPSI.h"
41 ClassImp(AliAnalysisTaskSEBkgLikeSignJPSI)
43 //________________________________________________________________________
44 AliAnalysisTaskSEBkgLikeSignJPSI::AliAnalysisTaskSEBkgLikeSignJPSI():
65 // Default constructor
69 //________________________________________________________________________
70 AliAnalysisTaskSEBkgLikeSignJPSI::AliAnalysisTaskSEBkgLikeSignJPSI(const char *name):
71 AliAnalysisTaskSE(name),
91 // Standard constructor
93 // Output slot #1 writes into a TList container
94 DefineOutput(1,TList::Class()); //My private output
97 //________________________________________________________________________
98 AliAnalysisTaskSEBkgLikeSignJPSI::~AliAnalysisTaskSEBkgLikeSignJPSI()
112 //________________________________________________________________________
113 void AliAnalysisTaskSEBkgLikeSignJPSI::Init()
117 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI::Init() \n");
119 gROOT->LoadMacro("ConfigVertexingHF.C");
121 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
127 //________________________________________________________________________
128 void AliAnalysisTaskSEBkgLikeSignJPSI::UserCreateOutputObjects()
130 // Create the output container
132 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI::UserCreateOutputObjects() \n");
134 // Several histograms are more conveniently managed in a TList
135 fOutput = new TList();
138 fHistMassJPSI = new TH1F("fHistMassJPSI", "J/#Psi invariant mass; M [GeV]; Entries",200,2.8,3.25);
139 fHistMassJPSI->Sumw2();
140 fHistMassJPSI->SetMinimum(0);
141 fOutput->Add(fHistMassJPSI);
143 fHistMassLS = new TH1F("fHistMassLS", "Like sign pairs invariant mass; M [GeV]; Entries",200,2.8,3.25);
144 fHistMassLS->Sumw2();
145 fHistMassLS->SetMinimum(0);
146 fOutput->Add(fHistMassLS);
148 fHistCtsJPSI = new TH1F("fHistCtsJPSI", "J/#Psi cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
149 fHistCtsJPSI->Sumw2();
150 fHistCtsJPSI->SetMinimum(0);
151 fOutput->Add(fHistCtsJPSI);
153 fHistCtsLS = new TH1F("fHistCtsLS", "Like sign pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
155 fHistCtsLS->SetMinimum(0);
156 fOutput->Add(fHistCtsLS);
158 fHistCtsLSpos = new TH1F("fHistCtsLSpos", "Like sign ++ pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
159 fHistCtsLSpos->Sumw2();
160 fHistCtsLSpos->SetMinimum(0);
161 fOutput->Add(fHistCtsLSpos);
163 fHistCtsLSneg = new TH1F("fHistCtsLSneg", "Like sign -- pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
164 fHistCtsLSneg->Sumw2();
165 fHistCtsLSneg->SetMinimum(0);
166 fOutput->Add(fHistCtsLSneg);
168 fHistCPtaJPSI = new TH1F("fHistCPtaJPSI", "J/#Psi cosine of pointing angle; Cos#Theta_{point}; Entries",200,-1.,1.);
169 fHistCPtaJPSI->Sumw2();
170 fHistCPtaJPSI->SetMinimum(0);
171 fOutput->Add(fHistCPtaJPSI);
173 fHistCPtaLS = new TH1F("fHistCPtaLS", "Like sign pairs cosine of pointing angle; Cos#Theta_{point}; Entries",200,-1.,1.);
174 fHistCPtaLS->Sumw2();
175 fHistCPtaLS->SetMinimum(0);
176 fOutput->Add(fHistCPtaLS);
178 fHistd0d0JPSI = new TH1F("fHistd0d0JPSI", "J/#Psi product of impact parameters; d0xd0 [#mu m^{2}]; Entries",200,-100000.,100000.);
179 fHistd0d0JPSI->Sumw2();
180 fHistd0d0JPSI->SetMinimum(0);
181 fOutput->Add(fHistd0d0JPSI);
183 fHistd0d0LS = new TH1F("fHistd0d0LS", "Like sign pairs product of impact parameters; d0xd0 [#mu m^{2}]; Entries",200,-100000.,100000.);
184 fHistd0d0LS->Sumw2();
185 fHistd0d0LS->SetMinimum(0);
186 fOutput->Add(fHistd0d0LS);
188 fHistDCAJPSI = new TH1F("fHistDCAJPSI", "J/#Psi distance of closest approach; dca [10^{2}#mu m]; Entries",100,0.,5.);
189 fHistDCAJPSI->Sumw2();
190 fHistDCAJPSI->SetMinimum(0);
191 fOutput->Add(fHistDCAJPSI);
193 fHistDCALS = new TH1F("fHistDCALS", "Like sign pairs distance of closest approach; dca [10^{2}#mu m]; Entries",100,0.,5.);
195 fHistDCALS->SetMinimum(0);
196 fOutput->Add(fHistDCALS);
201 //________________________________________________________________________
202 void AliAnalysisTaskSEBkgLikeSignJPSI::UserExec(Option_t */*option*/)
204 // Execute analysis for current event:
205 // heavy flavor candidates association to MC truth
207 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
209 TClonesArray *arrayJPSItoEle = 0;
210 TClonesArray *arrayLikeSign = 0;
212 if(!aod && AODEvent() && IsStandardAOD()) {
213 // In case there is an AOD handler writing a standard AOD, use the AOD
214 // event in memory rather than the input (ESD) event.
215 aod = dynamic_cast<AliAODEvent*> (AODEvent());
216 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
217 // have to taken from the AOD event hold by the AliAODExtension
218 AliAODHandler* aodHandler = (AliAODHandler*)
219 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
220 if(aodHandler->GetExtensions()) {
221 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
222 AliAODEvent *aodFromExt = ext->GetAOD();
223 // load Jpsi candidates
224 arrayJPSItoEle=(TClonesArray*)aodFromExt->GetList()->FindObject("JPSItoEle");
225 // load like sign candidates
226 arrayLikeSign=(TClonesArray*)aodFromExt->GetList()->FindObject("LikeSign2Prong");
229 // load Jpsi candidates
230 arrayJPSItoEle=(TClonesArray*)aod->GetList()->FindObject("JPSItoEle");
231 // load like sign candidates
232 arrayLikeSign=(TClonesArray*)aod->GetList()->FindObject("LikeSign2Prong");
236 if(!arrayJPSItoEle) {
237 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: JPSItoEle branch not found!\n");
241 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: LikeSign2Prong branch not found!\n");
245 // fix for temporary bug in ESDfilter
246 // the AODs with null vertex pointer didn't pass the PhysSel
247 if(!aod->GetPrimaryVertex()) return;
249 // AOD primary vertex
250 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
252 // make trkIDtoEntry register (temporary)
253 Int_t trkIDtoEntry[100000];
254 for(Int_t it=0;it<aod->GetNumberOfTracks();it++) {
255 AliAODTrack *track = aod->GetTrack(it);
256 trkIDtoEntry[track->GetID()]=it;
259 // loop over Like sign candidates
260 Int_t nPosPairs=0,nNegPairs=0;
261 Int_t nLikeSign = arrayLikeSign->GetEntriesFast();
262 if(fDebug>1) printf("+++\n+++Number of like sign pairs ---> %d \n+++\n", nLikeSign);
264 for(Int_t iLikeSign = 0; iLikeSign < nLikeSign; iLikeSign++) {
265 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayLikeSign->UncheckedAt(iLikeSign);
266 Bool_t unsetvtx=kFALSE;
267 if(!d->GetOwnPrimaryVtx()) {
268 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
272 if(d->SelectBtoJPSI(fVHF->GetBtoJPSICuts(),okBtoJPSIls)) {
273 AliAODTrack *trk0 = (AliAODTrack*)d->GetDaughter(0);
274 fHistMassLS->Fill(d->InvMassJPSIee());
275 fHistCPtaLS->Fill(d->CosPointingAngle());
276 fHistd0d0LS->Fill(1e8*d->Prodd0d0());
277 fHistDCALS->Fill(100*d->GetDCA());
278 //PostData(1,fOutput);
280 trk0=aod->GetTrack(trkIDtoEntry[d->GetProngID(0)]);
281 printf("references to standard AOD not available \n");
283 if((trk0->Charge())==1) {
285 fHistCtsLS->Fill(d->CosThetaStar(0,443,11,11));
286 fHistCtsLSpos->Fill(d->CosThetaStar(0,443,11,11));
287 //PostData(1,fOutput);
290 fHistCtsLS->Fill(d->CosThetaStarJPSI());
291 fHistCtsLSneg->Fill(d->CosThetaStarJPSI());
292 //PostData(1,fOutput);
296 if(unsetvtx) d->UnsetOwnPrimaryVtx();
299 if(fDebug>1) printf("------------ N. of positive pairs in Event ----- %d \n", nPosPairs);
300 if(fDebug>1) printf("------------ N. of negative pairs in Event ----- %d \n", nNegPairs);
302 fTotPosPairs += nPosPairs;
303 fTotNegPairs += nNegPairs;
305 // loop over JPSI candidates
306 Int_t nBtoJpsiToEle = arrayJPSItoEle->GetEntriesFast();
307 if(fDebug>1) printf("Number of like JPSI -> ee candidates ---> %d \n", nBtoJpsiToEle);
309 for (Int_t iBtoJpsiToEle = 0; iBtoJpsiToEle < nBtoJpsiToEle; iBtoJpsiToEle++) {
310 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayJPSItoEle->UncheckedAt(iBtoJpsiToEle);
311 Bool_t unsetvtx=kFALSE;
312 if(!d->GetOwnPrimaryVtx()) {
313 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
317 if(d->SelectBtoJPSI(fVHF->GetBtoJPSICuts(),okBtoJPSI)) {
318 fHistMassJPSI->Fill(d->InvMassJPSIee());
319 fHistCtsJPSI->Fill(d->CosThetaStarJPSI());
320 fHistd0d0JPSI->Fill(1e8*d->Prodd0d0());
321 fHistCPtaJPSI->Fill(d->CosPointingAngle());
322 fHistDCAJPSI->Fill(100*d->GetDCA());
325 if(unsetvtx) d->UnsetOwnPrimaryVtx();
331 //________________________________________________________________________
332 void AliAnalysisTaskSEBkgLikeSignJPSI::Terminate(Option_t */*option*/)
334 // Terminate analysis
336 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI: Terminate() \n");
338 fOutput = dynamic_cast<TList*> (GetOutputData(1));
340 printf("ERROR: fOutput not available\n");
344 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs*fTotNegPairs);
346 fHistMassJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMassJPSI"));
347 fHistMassLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMassLS"));
348 fHistCtsJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsJPSI"));
349 fHistCtsLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLS"));
350 fHistCtsLSpos = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLSpos"));
351 fHistCtsLSneg = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLSneg"));
352 fHistCPtaJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCPtaJPSI"));
353 fHistCPtaLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCPtaLS"));
354 fHistd0d0JPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistd0d0JPSI"));
355 fHistd0d0LS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistd0d0LS"));
356 fHistDCAJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistDCAJPSI"));
357 fHistDCALS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistDCALS"));
359 if(fLsNormalization>0.) {
360 fHistMassLS->Scale((1/fLsNormalization)*fHistMassLS->GetEntries());
361 fHistCtsLS->Scale((1/fLsNormalization)*fHistCtsLS->GetEntries());
362 fHistCtsLSpos->Scale((1/fLsNormalization)*fHistCtsLSpos->GetEntries());
363 fHistCtsLSneg->Scale((1/fLsNormalization)*fHistCtsLSneg->GetEntries());
364 fHistCPtaLS->Scale((1/fLsNormalization)*fHistCPtaLS->GetEntries());
365 fHistd0d0LS->Scale((1/fLsNormalization)*fHistd0d0LS->GetEntries());
366 fHistDCALS->Scale((1/fLsNormalization)*fHistDCALS->GetEntries());