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 **************************************************************************/
16 /////////////////////////////////////////////////////////////
18 // AliAnalysisTaskSE for the comparison of heavy flavor
19 // decay candidates with the MC truth.
21 // Author: A.Dainese, andrea.dainese@lnl.infn.it
22 /////////////////////////////////////////////////////////////
24 #include <TClonesArray.h>
29 #include "AliAnalysisManager.h"
30 #include "AliAODHandler.h"
31 #include "AliAODEvent.h"
32 #include "AliAODVertex.h"
33 #include "AliESDtrack.h"
34 #include "AliAODTrack.h"
35 #include "AliAODMCHeader.h"
36 #include "AliAODMCParticle.h"
37 #include "AliAODRecoDecayHF2Prong.h"
38 #include "AliAODRecoCascadeHF.h"
39 #include "AliAnalysisVertexingHF.h"
40 #include "AliAnalysisTaskSE.h"
41 #include "AliAnalysisTaskSECompareHF.h"
43 ClassImp(AliAnalysisTaskSECompareHF)
46 //________________________________________________________________________
47 AliAnalysisTaskSECompareHF::AliAnalysisTaskSECompareHF():
55 // Default constructor
57 // NO DefineOutput() HERE (ONLY IN STANDARD CONSTRUCTOR)
60 //________________________________________________________________________
61 AliAnalysisTaskSECompareHF::AliAnalysisTaskSECompareHF(const char *name):
62 AliAnalysisTaskSE(name),
69 // Standard constructor
71 // Output slot #1 writes into a TList container
72 DefineOutput(1,TList::Class()); //My private output
73 // Output slot #2 writes into a TNtuple container
74 DefineOutput(2,TNtuple::Class()); //My private output
77 //________________________________________________________________________
78 AliAnalysisTaskSECompareHF::~AliAnalysisTaskSECompareHF()
91 //________________________________________________________________________
92 void AliAnalysisTaskSECompareHF::Init()
96 if(fDebug > 1) printf("AnalysisTaskSECompareHF::Init() \n");
98 gROOT->LoadMacro("ConfigVertexingHF.C");
100 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
106 //________________________________________________________________________
107 void AliAnalysisTaskSECompareHF::UserCreateOutputObjects()
109 // Create the output container
111 if(fDebug > 1) printf("AnalysisTaskSECompareHF::UserCreateOutputObjects() \n");
113 // Several histograms are more conveniently managed in a TList
114 fOutput = new TList();
117 fHistMass = new TH1F("fHistMass", "D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
119 fHistMass->SetMinimum(0);
120 fOutput->Add(fHistMass);
122 fHistNEvents = new TH1F("fHistNEvents", "Number of processed events; ; Events",3,-1.5,1.5);
123 fHistNEvents->Sumw2();
124 fHistNEvents->SetMinimum(0);
125 fOutput->Add(fHistNEvents);
127 OpenFile(2); // 2 is the slot number of the ntuple
128 fNtupleCmp = new TNtuple("fNtupleCmp","Charm comparison","pdg:nprongs:VxRec:VxTrue:ErrVx:VyRec:VyTrue:ErrVy:VzRec:VzTrue:ErrVz:Chi2toNDF:PtRec:Mrec:CPta:Prodd0");
133 //________________________________________________________________________
134 void AliAnalysisTaskSECompareHF::UserExec(Option_t */*option*/)
136 // Execute analysis for current event:
137 // heavy flavor candidates association to MC truth
140 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
142 TClonesArray *inputArrayVertices = 0;
143 TClonesArray *inputArrayD0toKpi = 0;
144 TClonesArray *inputArrayDstar = 0;
146 if(!aod && AODEvent() && IsStandardAOD()) {
147 // In case there is an AOD handler writing a standard AOD, use the AOD
148 // event in memory rather than the input (ESD) event.
149 aod = dynamic_cast<AliAODEvent*> (AODEvent());
150 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
151 // have to taken from the AOD event hold by the AliAODExtension
152 AliAODHandler* aodHandler = (AliAODHandler*)
153 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
154 if(aodHandler->GetExtensions()) {
155 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
156 AliAODEvent *aodFromExt = ext->GetAOD();
158 inputArrayVertices = (TClonesArray*)aodFromExt->GetList()->FindObject("VerticesHF");
159 // load D0->Kpi candidates
160 inputArrayD0toKpi = (TClonesArray*)aodFromExt->GetList()->FindObject("D0toKpi");
161 // load D*+ candidates
162 inputArrayDstar = (TClonesArray*)aodFromExt->GetList()->FindObject("Dstar");
166 inputArrayVertices = (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
167 // load D0->Kpi candidates
168 inputArrayD0toKpi = (TClonesArray*)aod->GetList()->FindObject("D0toKpi");
169 // load D*+ candidates
170 inputArrayDstar = (TClonesArray*)aod->GetList()->FindObject("Dstar");
174 if(!inputArrayVertices) {
175 printf("AliAnalysisTaskSECompareHF::UserExec: Vertices branch not found!\n");
178 if(!inputArrayD0toKpi) {
179 printf("AliAnalysisTaskSECompareHF::UserExec: D0toKpi branch not found!\n");
182 if(!inputArrayDstar) {
183 printf("AliAnalysisTaskSECompareHF::UserExec: Dstar branch not found!\n");
188 fHistNEvents->Fill(0); // count event
189 // Post the data already here
192 // AOD primary vertex
193 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
197 TClonesArray *mcArray =
198 (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
200 printf("AliAnalysisTaskSECompareHF::UserExec: MC particles branch not found!\n");
205 AliAODMCHeader *mcHeader =
206 (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
208 printf("AliAnalysisTaskSECompareHF::UserExec: MC header branch not found!\n");
212 Int_t nprongs,lab,okD0,okD0bar,pdg;
214 Double_t invmass,posRec[3],posTrue[3],covRec[6],errx,erry,errz;
215 AliAODRecoDecayHF2Prong *d2=0;
216 AliAODRecoDecayHF3Prong *d3=0;
217 AliAODRecoDecayHF4Prong *d4=0;
219 Int_t pdgDgD0toKpi[2]={321,211};
220 Int_t pdgDgDplustoKpipi[3]={321,211,211};
221 Int_t pdgDgD0toKpipipi[4]={321,211,211,211};
223 // loop over vertices
224 Int_t nVertices = inputArrayVertices->GetEntriesFast();
225 if(fDebug>1) printf("Number of vertices: %d\n",nVertices);
227 for (Int_t iVtx = 0; iVtx < nVertices; iVtx++) {
228 AliAODVertex *vtx = (AliAODVertex*)inputArrayVertices->UncheckedAt(iVtx);
231 vtx->GetCovarianceMatrix(covRec);
232 errx=1.; erry=1.; errz=1.;
233 if(covRec[0]>0) errx = TMath::Sqrt(covRec[0]);
234 if(covRec[2]>0) erry = TMath::Sqrt(covRec[2]);
235 if(covRec[5]>0) errz = TMath::Sqrt(covRec[5]);
238 // get parent AliAODRecoDecayHF
239 nprongs= vtx->GetNDaughters();
240 // check that the daughters have kITSrefit and kTPCrefit
241 Bool_t allDgOK=kTRUE;
242 for(Int_t idg=0; idg<nprongs; idg++) {
243 AliAODTrack *track = (AliAODTrack*)vtx->GetDaughter(idg);
244 if(!(track->GetStatus()&AliESDtrack::kITSrefit)) allDgOK=kFALSE;
245 if(!(track->GetStatus()&AliESDtrack::kTPCrefit)) allDgOK=kFALSE;
247 if(!allDgOK) continue;
251 case 2: // look for D0->Kpi
252 d2 = (AliAODRecoDecayHF2Prong*)vtx->GetParent();
253 if(d2->IsLikeSign()) continue;
254 if(d2->Charge() != 0) continue; // these are D*
255 lab = d2->MatchToMC(421,mcArray,2,pdgDgD0toKpi);
258 if(!d2->GetOwnPrimaryVtx()) {
259 d2->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
263 if(d2->SelectD0(fVHF->GetD0toKpiCuts(),okD0,okD0bar)) { // beware: cuts may bias the resolution!
264 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
265 pdg = dMC->GetPdgCode();
266 invmass = (pdg==421 ? d2->InvMassD0() : d2->InvMassD0bar());
267 // get a daughter for true pos of decay vertex
268 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
269 dg0MC->XvYvZv(posTrue);
270 fHistMass->Fill(invmass);
271 // Post the data already here
273 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
274 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
275 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
276 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
277 (Float_t)d2->GetReducedChi2(),(Float_t)d2->Pt(),(Float_t)invmass,
278 (Float_t)d2->CosPointingAngle(),(Float_t)d2->Prodd0d0()};
279 fNtupleCmp->Fill(tmp);
280 PostData(2,fNtupleCmp);
282 if(unsetvtx) d2->UnsetOwnPrimaryVtx();
285 case 3: // look for D+
286 d3 = (AliAODRecoDecayHF3Prong*)vtx->GetParent();
287 if(d3->IsLikeSign()) continue;
288 lab = d3->MatchToMC(411,mcArray,3,pdgDgDplustoKpipi);
291 if(!d3->GetOwnPrimaryVtx()) {
292 d3->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
295 //if(d3->SelectDplus(fVHF->GetDplusCuts())) {
296 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
297 pdg = dMC->GetPdgCode();
298 invmass = d3->InvMassDplus();
299 // get a daughter for true pos of decay vertex
300 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
301 dg0MC->XvYvZv(posTrue);
302 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
303 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
304 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
305 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
306 (Float_t)d3->GetReducedChi2(),(Float_t)d3->Pt(),(Float_t)invmass,
307 (Float_t)d3->CosPointingAngle(),(Float_t)(d3->Getd0Prong(0)*d3->Getd0Prong(1)*d3->Getd0Prong(2))};
308 fNtupleCmp->Fill(tmp);
309 PostData(2,fNtupleCmp);
311 if(unsetvtx) d3->UnsetOwnPrimaryVtx();
314 case 4: // look for D0->Kpipipi
315 d4 = (AliAODRecoDecayHF4Prong*)vtx->GetParent();
316 if(d4->IsLikeSign()) continue;
317 lab = d4->MatchToMC(421,mcArray,4,pdgDgD0toKpipipi);
320 if(!d4->GetOwnPrimaryVtx()) {
321 d4->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
325 //if(d4->SelectD0(fVHF->GetD0to4ProngsCuts(),okD0,okD0bar)) {
326 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
327 pdg = dMC->GetPdgCode();
328 //invmass = (pdg==421 ? d->InvMassD0() : d->InvMassD0bar());
329 invmass = 10.; // dummy
330 // get a daughter for true pos of decay vertex
331 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
332 dg0MC->XvYvZv(posTrue);
333 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
334 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
335 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
336 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
337 (Float_t)d4->GetReducedChi2(),(Float_t)d4->Pt(),(Float_t)invmass,
338 (Float_t)d4->CosPointingAngle(),(Float_t)(d4->Getd0Prong(0)*d4->Getd0Prong(1)*d4->Getd0Prong(2)*d4->Getd0Prong(3))};
339 fNtupleCmp->Fill(tmp);
340 PostData(2,fNtupleCmp);
342 if(unsetvtx) d4->UnsetOwnPrimaryVtx();
347 } // end loop on vertices
350 // loop over D*+ candidates
352 for (Int_t iDstar = 0; iDstar < inputArrayDstar->GetEntries(); iDstar++) {
353 AliAODRecoCascadeHF *c = (AliAODRecoCascadeHF*)inputArrayDstar->UncheckedAt(iDstar);
354 Int_t labDstar = c->MatchToMC(413,421,mcArray);
355 if(labDstar>=0) printf("GOOD MATCH FOR D*+\n");
359 // Post the data already here
361 PostData(2,fNtupleCmp);
365 //________________________________________________________________________
366 void AliAnalysisTaskSECompareHF::Terminate(Option_t */*option*/)
368 // Terminate analysis
370 if(fDebug > 1) printf("AnalysisTaskSECompareHF: Terminate() \n");
372 fOutput = dynamic_cast<TList*> (GetOutputData(1));
374 printf("ERROR: fOutput not available\n");
378 fHistMass = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMass"));
379 fHistNEvents = dynamic_cast<TH1F*>(fOutput->FindObject("fHistNEvents"));
381 //fNtupleCmp = dynamic_cast<TNtuple*> (GetOutputData(2));