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>
28 #include "AliAODEvent.h"
29 #include "AliAODVertex.h"
30 #include "AliAODTrack.h"
31 #include "AliAODRecoDecayHF2Prong.h"
32 #include "AliAnalysisTaskSE.h"
33 #include "AliAnalysisTaskSEBkgLikeSignJPSI.h"
35 ClassImp(AliAnalysisTaskSEBkgLikeSignJPSI)
37 //________________________________________________________________________
38 AliAnalysisTaskSEBkgLikeSignJPSI::AliAnalysisTaskSEBkgLikeSignJPSI():
59 // Default constructor
61 // Output slot #1 writes into a TList container
62 DefineOutput(1,TList::Class()); //My private output
65 //________________________________________________________________________
66 AliAnalysisTaskSEBkgLikeSignJPSI::AliAnalysisTaskSEBkgLikeSignJPSI(const char *name):
67 AliAnalysisTaskSE(name),
87 // Default constructor
89 // Output slot #1 writes into a TList container
90 DefineOutput(1,TList::Class()); //My private output
93 //________________________________________________________________________
94 AliAnalysisTaskSEBkgLikeSignJPSI::~AliAnalysisTaskSEBkgLikeSignJPSI()
108 //________________________________________________________________________
109 void AliAnalysisTaskSEBkgLikeSignJPSI::Init()
113 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI::Init() \n");
115 gROOT->LoadMacro("ConfigVertexingHF.C");
117 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
125 //________________________________________________________________________
126 void AliAnalysisTaskSEBkgLikeSignJPSI::UserCreateOutputObjects()
128 // Create the output container
130 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI::UserCreateOutputObjects() \n");
132 // Several histograms are more conveniently managed in a TList
133 fOutput = new TList();
136 fHistMassJPSI = new TH1F("fHistMassJPSI", "J/#Psi invariant mass; M [GeV]; Entries",200,2.8,3.25);
137 fHistMassJPSI->Sumw2();
138 fHistMassJPSI->SetMinimum(0);
139 fOutput->Add(fHistMassJPSI);
141 fHistMassLS = new TH1F("fHistMassLS", "Like sign pairs invariant mass; M [GeV]; Entries",200,2.8,3.25);
142 fHistMassLS->Sumw2();
143 fHistMassLS->SetMinimum(0);
144 fOutput->Add(fHistMassLS);
146 fHistCtsJPSI = new TH1F("fHistCtsJPSI", "J/#Psi cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
147 fHistCtsJPSI->Sumw2();
148 fHistCtsJPSI->SetMinimum(0);
149 fOutput->Add(fHistCtsJPSI);
151 fHistCtsLS = new TH1F("fHistCtsLS", "Like sign pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
153 fHistCtsLS->SetMinimum(0);
154 fOutput->Add(fHistCtsLS);
156 fHistCtsLSpos = new TH1F("fHistCtsLSpos", "Like sign ++ pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
157 fHistCtsLSpos->Sumw2();
158 fHistCtsLSpos->SetMinimum(0);
159 fOutput->Add(fHistCtsLSpos);
161 fHistCtsLSneg = new TH1F("fHistCtsLSneg", "Like sign -- pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
162 fHistCtsLSneg->Sumw2();
163 fHistCtsLSneg->SetMinimum(0);
164 fOutput->Add(fHistCtsLSneg);
166 fHistCPtaJPSI = new TH1F("fHistCPtaJPSI", "J/#Psi cosine of pointing angle; Cos#Theta_{point}; Entries",200,-1.,1.);
167 fHistCPtaJPSI->Sumw2();
168 fHistCPtaJPSI->SetMinimum(0);
169 fOutput->Add(fHistCPtaJPSI);
171 fHistCPtaLS = new TH1F("fHistCPtaLS", "Like sign pairs cosine of pointing angle; Cos#Theta_{point}; Entries",200,-1.,1.);
172 fHistCPtaLS->Sumw2();
173 fHistCPtaLS->SetMinimum(0);
174 fOutput->Add(fHistCPtaLS);
176 fHistd0d0JPSI = new TH1F("fHistd0d0JPSI", "J/#Psi product of impact parameters; d0xd0 [#mu m^{2}]; Entries",200,-100000.,100000.);
177 fHistd0d0JPSI->Sumw2();
178 fHistd0d0JPSI->SetMinimum(0);
179 fOutput->Add(fHistd0d0JPSI);
181 fHistd0d0LS = new TH1F("fHistd0d0LS", "Like sign pairs product of impact parameters; d0xd0 [#mu m^{2}]; Entries",200,-100000.,100000.);
182 fHistd0d0LS->Sumw2();
183 fHistd0d0LS->SetMinimum(0);
184 fOutput->Add(fHistd0d0LS);
186 fHistDCAJPSI = new TH1F("fHistDCAJPSI", "J/#Psi distance of closest approach; dca [10^{2}#mu m]; Entries",100,0.,5.);
187 fHistDCAJPSI->Sumw2();
188 fHistDCAJPSI->SetMinimum(0);
189 fOutput->Add(fHistDCAJPSI);
191 fHistDCALS = new TH1F("fHistDCALS", "Like sign pairs distance of closest approach; dca [10^{2}#mu m]; Entries",100,0.,5.);
193 fHistDCALS->SetMinimum(0);
194 fOutput->Add(fHistDCALS);
199 //________________________________________________________________________
200 void AliAnalysisTaskSEBkgLikeSignJPSI::UserExec(Option_t */*option*/)
202 // Execute analysis for current event:
203 // heavy flavor candidates association to MC truth
205 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
207 // load heavy flavour vertices
208 TClonesArray *arrayVerticesHF =
209 (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
210 if(!arrayVerticesHF) {
211 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: VerticesHF branch not found!\n");
215 // load JPSI->ee candidates
216 TClonesArray *arrayJPSItoEle =
217 (TClonesArray*)aod->GetList()->FindObject("JPSItoEle");
218 if(!arrayJPSItoEle) {
219 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: JPSItoEle branch not found!\n");
223 // load like sign candidates
224 TClonesArray *arrayLikeSign =
225 (TClonesArray*)aod->GetList()->FindObject("LikeSign2Prong");
227 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: LikeSign2Prong branch not found!\n");
231 // AOD primary vertex
232 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
234 // make trkIDtoEntry register (temporary)
235 Int_t trkIDtoEntry[100000];
236 for(Int_t it=0;it<aod->GetNumberOfTracks();it++) {
237 AliAODTrack *track = aod->GetTrack(it);
238 trkIDtoEntry[track->GetID()]=it;
241 // loop over Like sign candidates
242 Int_t nPosPairs=0,nNegPairs=0;
243 Int_t nLikeSign = arrayLikeSign->GetEntriesFast();
244 printf("+++\n+++Number of like sign pairs ---> %d \n+++\n", nLikeSign);
246 for(Int_t iLikeSign = 0; iLikeSign < nLikeSign; iLikeSign++) {
247 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayLikeSign->UncheckedAt(iLikeSign);
248 Bool_t unsetvtx=kFALSE;
249 if(!d->GetOwnPrimaryVtx()) {
250 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
254 if(d->SelectBtoJPSI(fVHF->GetBtoJPSICuts(),okBtoJPSIls)) {
255 fHistMassLS->Fill(d->InvMassJPSIee());
256 fHistCPtaLS->Fill(d->CosPointingAngle());
257 fHistd0d0LS->Fill(1e8*d->Prodd0d0());
258 fHistCtsLS->Fill(d->CosThetaStarJPSI());
259 fHistDCALS->Fill(100*d->GetDCA());
261 AliAODTrack *trk0 = (AliAODTrack*)d->GetDaughter(0);
263 trk0=aod->GetTrack(trkIDtoEntry[d->GetProngID(0)]);
264 printf("references to standard AOD not available \n");
266 if((trk0->Charge())==1) {
268 fHistCtsLSpos->Fill(d->CosThetaStarJPSI());
272 fHistCtsLSneg->Fill(d->CosThetaStarJPSI());
277 if(unsetvtx) d->UnsetOwnPrimaryVtx();
280 printf("------------ N. of positive pairs in Event ----- %d \n", nPosPairs);
281 printf("------------ N. of negative pairs in Event ----- %d \n", nNegPairs);
283 fTotPosPairs += nPosPairs;
284 fTotNegPairs += nNegPairs;
286 // loop over JPSI candidates
287 Int_t nBtoJpsiToEle = arrayJPSItoEle->GetEntriesFast();
288 printf("Number of like JPSI -> ee candidates ---> %d \n", nBtoJpsiToEle);
290 for (Int_t iBtoJpsiToEle = 0; iBtoJpsiToEle < nBtoJpsiToEle; iBtoJpsiToEle++) {
291 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayJPSItoEle->UncheckedAt(iBtoJpsiToEle);
292 Bool_t unsetvtx=kFALSE;
293 if(!d->GetOwnPrimaryVtx()) {
294 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
298 if(d->SelectBtoJPSI(fVHF->GetBtoJPSICuts(),okBtoJPSI)) {
299 fHistMassJPSI->Fill(d->InvMassJPSIee());
300 fHistCtsJPSI->Fill(d->CosThetaStarJPSI());
301 fHistd0d0JPSI->Fill(1e8*d->Prodd0d0());
302 fHistCPtaJPSI->Fill(d->CosPointingAngle());
303 fHistDCAJPSI->Fill(100*d->GetDCA());
306 if(unsetvtx) d->UnsetOwnPrimaryVtx();
312 //________________________________________________________________________
313 void AliAnalysisTaskSEBkgLikeSignJPSI::Terminate(Option_t */*option*/)
315 // Terminate analysis
317 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI: Terminate() \n");
319 fOutput = dynamic_cast<TList*> (GetOutputData(1));
321 printf("ERROR: fOutput not available\n");
325 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs*fTotNegPairs);
327 fHistMassJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMassJPSI"));
328 fHistMassLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMassLS"));
329 fHistCtsJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsJPSI"));
330 fHistCtsLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLS"));
331 fHistCtsLSpos = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLSpos"));
332 fHistCtsLSneg = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLSneg"));
333 fHistCPtaJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCPtaJPSI"));
334 fHistCPtaLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCPtaLS"));
335 fHistd0d0JPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistd0d0JPSI"));
336 fHistd0d0LS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistd0d0LS"));
337 fHistDCAJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistDCAJPSI"));
338 fHistDCALS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistDCALS"));
340 fHistMassLS->Scale((1/fLsNormalization)*fHistMassLS->GetEntries());
341 fHistCtsLS->Scale((1/fLsNormalization)*fHistCtsLS->GetEntries());
342 fHistCtsLSpos->Scale((1/fLsNormalization)*fHistCtsLSpos->GetEntries());
343 fHistCtsLSneg->Scale((1/fLsNormalization)*fHistCtsLSneg->GetEntries());
344 fHistCPtaLS->Scale((1/fLsNormalization)*fHistCPtaLS->GetEntries());
345 fHistd0d0LS->Scale((1/fLsNormalization)*fHistd0d0LS->GetEntries());
346 fHistDCALS->Scale((1/fLsNormalization)*fHistDCALS->GetEntries());