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 "AliAODEvent.h"
30 #include "AliAODVertex.h"
31 #include "AliESDtrack.h"
32 #include "AliAODTrack.h"
33 #include "AliAODMCHeader.h"
34 #include "AliAODMCParticle.h"
35 #include "AliAODRecoDecayHF2Prong.h"
36 #include "AliAODRecoCascadeHF.h"
37 #include "AliAnalysisVertexingHF.h"
38 #include "AliAnalysisTaskSE.h"
39 #include "AliAnalysisTaskSECompareHF.h"
41 ClassImp(AliAnalysisTaskSECompareHF)
44 //________________________________________________________________________
45 AliAnalysisTaskSECompareHF::AliAnalysisTaskSECompareHF():
53 // Default constructor
55 // NO DefineOutput() HERE (ONLY IN STANDARD CONSTRUCTOR)
58 //________________________________________________________________________
59 AliAnalysisTaskSECompareHF::AliAnalysisTaskSECompareHF(const char *name):
60 AliAnalysisTaskSE(name),
67 // Standard constructor
69 // Output slot #1 writes into a TList container
70 DefineOutput(1,TList::Class()); //My private output
71 // Output slot #2 writes into a TNtuple container
72 DefineOutput(2,TNtuple::Class()); //My private output
75 //________________________________________________________________________
76 AliAnalysisTaskSECompareHF::~AliAnalysisTaskSECompareHF()
89 //________________________________________________________________________
90 void AliAnalysisTaskSECompareHF::Init()
94 if(fDebug > 1) printf("AnalysisTaskSECompareHF::Init() \n");
96 gROOT->LoadMacro("ConfigVertexingHF.C");
98 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
104 //________________________________________________________________________
105 void AliAnalysisTaskSECompareHF::UserCreateOutputObjects()
107 // Create the output container
109 if(fDebug > 1) printf("AnalysisTaskSECompareHF::UserCreateOutputObjects() \n");
111 // Several histograms are more conveniently managed in a TList
112 fOutput = new TList();
115 fHistMass = new TH1F("fHistMass", "D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
117 fHistMass->SetMinimum(0);
118 fOutput->Add(fHistMass);
120 fHistNEvents = new TH1F("fHistNEvents", "Number of processed events; ; Events",3,-1.5,1.5);
121 fHistNEvents->Sumw2();
122 fHistNEvents->SetMinimum(0);
123 fOutput->Add(fHistNEvents);
125 OpenFile(2); // 2 is the slot number of the ntuple
126 fNtupleCmp = new TNtuple("fNtupleCmp","Charm comparison","pdg:nprongs:VxRec:VxTrue:ErrVx:VyRec:VyTrue:ErrVy:VzRec:VzTrue:ErrVz:Chi2toNDF:PtRec:Mrec:CPta:Prodd0");
131 //________________________________________________________________________
132 void AliAnalysisTaskSECompareHF::UserExec(Option_t */*option*/)
134 // Execute analysis for current event:
135 // heavy flavor candidates association to MC truth
138 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
140 fHistNEvents->Fill(0); // count event
141 // Post the data already here
145 TClonesArray *inputArrayVertices =
146 (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
147 if(!inputArrayVertices) {
148 printf("AliAnalysisTaskSECompareHF::UserExec: Vertices branch not found!\n");
152 // load D0->Kpi candidates
153 TClonesArray *inputArrayD0toKpi =
154 (TClonesArray*)aod->GetList()->FindObject("D0toKpi");
155 if(!inputArrayD0toKpi) {
156 printf("AliAnalysisTaskSECompareHF::UserExec: D0toKpi branch not found!\n");
161 // load D*+ candidates
162 TClonesArray *inputArrayDstar =
163 (TClonesArray*)aod->GetList()->FindObject("Dstar");
164 if(!inputArrayDstar) {
165 printf("AliAnalysisTaskSECompareHF::UserExec: Dstar branch not found!\n");
170 // AOD primary vertex
171 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
175 TClonesArray *mcArray =
176 (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
178 printf("AliAnalysisTaskSECompareHF::UserExec: MC particles branch not found!\n");
183 AliAODMCHeader *mcHeader =
184 (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
186 printf("AliAnalysisTaskSECompareHF::UserExec: MC header branch not found!\n");
190 Int_t nprongs,lab,okD0,okD0bar,pdg;
192 Double_t invmass,posRec[3],posTrue[3],covRec[6],errx,erry,errz;
193 AliAODRecoDecayHF2Prong *d2=0;
194 AliAODRecoDecayHF3Prong *d3=0;
195 AliAODRecoDecayHF4Prong *d4=0;
197 // loop over vertices
198 Int_t nVertices = inputArrayVertices->GetEntriesFast();
199 printf("Number of vertices: %d\n",nVertices);
201 for (Int_t iVtx = 0; iVtx < nVertices; iVtx++) {
202 AliAODVertex *vtx = (AliAODVertex*)inputArrayVertices->UncheckedAt(iVtx);
205 vtx->GetCovarianceMatrix(covRec);
206 errx=1.; erry=1.; errz=1.;
207 if(covRec[0]>0) errx = TMath::Sqrt(covRec[0]);
208 if(covRec[2]>0) erry = TMath::Sqrt(covRec[2]);
209 if(covRec[5]>0) errz = TMath::Sqrt(covRec[5]);
212 // get parent AliAODRecoDecayHF
213 nprongs= vtx->GetNDaughters();
214 // check that the daughters have kITSrefit and kTPCrefit
215 Bool_t allDgOK=kTRUE;
216 for(Int_t idg=0; idg<nprongs; idg++) {
217 AliAODTrack *track = (AliAODTrack*)vtx->GetDaughter(idg);
218 if(!(track->GetStatus()&AliESDtrack::kITSrefit)) allDgOK=kFALSE;
219 if(!(track->GetStatus()&AliESDtrack::kTPCrefit)) allDgOK=kFALSE;
221 if(!allDgOK) continue;
225 case 2: // look for D0->Kpi
226 d2 = (AliAODRecoDecayHF2Prong*)vtx->GetParent();
227 if(d2->IsLikeSign()) continue;
228 if(d2->Charge() != 0) continue; // these are D*
229 lab = d2->MatchToMC(421,mcArray);
232 if(!d2->GetOwnPrimaryVtx()) {
233 d2->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
237 if(d2->SelectD0(fVHF->GetD0toKpiCuts(),okD0,okD0bar)) { // beware: cuts may bias the resolution!
238 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
239 pdg = dMC->GetPdgCode();
240 invmass = (pdg==421 ? d2->InvMassD0() : d2->InvMassD0bar());
241 fHistMass->Fill(invmass);
242 // Post the data already here
244 // get a daughter for true pos of decay vertex
245 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
246 dg0MC->XvYvZv(posTrue);
247 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
248 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
249 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
250 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
251 (Float_t)d2->GetReducedChi2(),(Float_t)d2->Pt(),(Float_t)invmass,
252 (Float_t)d2->CosPointingAngle(),(Float_t)d2->Prodd0d0()};
253 fNtupleCmp->Fill(tmp);
254 PostData(2,fNtupleCmp);
256 if(unsetvtx) d2->UnsetOwnPrimaryVtx();
259 case 3: // look for D+
260 d3 = (AliAODRecoDecayHF3Prong*)vtx->GetParent();
261 if(d3->IsLikeSign()) continue;
262 lab = d3->MatchToMC(411,mcArray);
265 if(!d3->GetOwnPrimaryVtx()) {
266 d3->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
269 //if(d3->SelectDplus(fVHF->GetDplusCuts())) {
270 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
271 pdg = dMC->GetPdgCode();
272 invmass = d3->InvMassDplus();
273 // get a daughter for true pos of decay vertex
274 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
275 dg0MC->XvYvZv(posTrue);
276 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
277 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
278 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
279 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
280 (Float_t)d3->GetReducedChi2(),(Float_t)d3->Pt(),(Float_t)invmass,
281 (Float_t)d3->CosPointingAngle(),(Float_t)(d3->Getd0Prong(0)*d3->Getd0Prong(1)*d3->Getd0Prong(2))};
282 fNtupleCmp->Fill(tmp);
283 PostData(2,fNtupleCmp);
285 if(unsetvtx) d3->UnsetOwnPrimaryVtx();
288 case 4: // look for D0->Kpipipi
289 d4 = (AliAODRecoDecayHF4Prong*)vtx->GetParent();
290 if(d4->IsLikeSign()) continue;
291 lab = d4->MatchToMC(421,mcArray);
294 if(!d4->GetOwnPrimaryVtx()) {
295 d4->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
299 //if(d4->SelectD0(fVHF->GetD0to4ProngsCuts(),okD0,okD0bar)) {
300 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
301 pdg = dMC->GetPdgCode();
302 //invmass = (pdg==421 ? d->InvMassD0() : d->InvMassD0bar());
303 invmass = 10.; // dummy
304 // get a daughter for true pos of decay vertex
305 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
306 dg0MC->XvYvZv(posTrue);
307 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
308 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
309 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
310 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
311 (Float_t)d4->GetReducedChi2(),(Float_t)d4->Pt(),(Float_t)invmass,
312 (Float_t)d4->CosPointingAngle(),(Float_t)(d4->Getd0Prong(0)*d4->Getd0Prong(1)*d4->Getd0Prong(2)*d4->Getd0Prong(3))};
313 fNtupleCmp->Fill(tmp);
314 PostData(2,fNtupleCmp);
316 if(unsetvtx) d4->UnsetOwnPrimaryVtx();
321 } // end loop on vertices
324 // loop over D*+ candidates
326 for (Int_t iDstar = 0; iDstar < inputArrayDstar->GetEntries(); iDstar++) {
327 AliAODRecoCascadeHF *c = (AliAODRecoCascadeHF*)inputArrayDstar->UncheckedAt(iDstar);
328 Int_t labDstar = c->MatchToMC(413,421,mcArray);
329 if(labDstar>=0) printf("GOOD MATCH FOR D*+\n");
333 // Post the data already here
335 PostData(2,fNtupleCmp);
339 //________________________________________________________________________
340 void AliAnalysisTaskSECompareHF::Terminate(Option_t */*option*/)
342 // Terminate analysis
344 if(fDebug > 1) printf("AnalysisTaskSECompareHF: Terminate() \n");
346 fOutput = dynamic_cast<TList*> (GetOutputData(1));
348 printf("ERROR: fOutput not available\n");
352 fHistMass = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMass"));
353 fHistNEvents = dynamic_cast<TH1F*>(fOutput->FindObject("fHistNEvents"));
355 //fNtupleCmp = dynamic_cast<TNtuple*> (GetOutputData(2));