]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/vertexingHF/AliAnalysisTaskSECompareHF.cxx
Warning fixed
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliAnalysisTaskSECompareHF.cxx
CommitLineData
3a219f60 1/**************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/////////////////////////////////////////////////////////////
17//
18// AliAnalysisTaskSE for the comparison of heavy flavor
19// decay candidates with the MC truth.
20//
21// Author: A.Dainese, andrea.dainese@lnl.infn.it
22/////////////////////////////////////////////////////////////
23
24#include <TClonesArray.h>
25#include <TNtuple.h>
26#include <TList.h>
27#include <TH1F.h>
28
29#include "AliAODEvent.h"
30#include "AliAODVertex.h"
607719af 31#include "AliESDtrack.h"
3a219f60 32#include "AliAODTrack.h"
33#include "AliAODMCHeader.h"
34#include "AliAODMCParticle.h"
35#include "AliAODRecoDecayHF2Prong.h"
b3999999 36#include "AliAODRecoCascadeHF.h"
2ae9c7cb 37#include "AliAnalysisVertexingHF.h"
3a219f60 38#include "AliAnalysisTaskSE.h"
39#include "AliAnalysisTaskSECompareHF.h"
40
41ClassImp(AliAnalysisTaskSECompareHF)
42
43
44//________________________________________________________________________
45AliAnalysisTaskSECompareHF::AliAnalysisTaskSECompareHF():
46AliAnalysisTaskSE(),
47fOutput(0),
2a041947 48fNtupleCmp(0),
12bfc069 49fHistMass(0),
a273622f 50fHistNEvents(0),
12bfc069 51fVHF(0)
3a219f60 52{
53 // Default constructor
3a219f60 54
9d6d35b0 55 // NO DefineOutput() HERE (ONLY IN STANDARD CONSTRUCTOR)
3a219f60 56}
57
58//________________________________________________________________________
59AliAnalysisTaskSECompareHF::AliAnalysisTaskSECompareHF(const char *name):
60AliAnalysisTaskSE(name),
61fOutput(0),
2a041947 62fNtupleCmp(0),
12bfc069 63fHistMass(0),
a273622f 64fHistNEvents(0),
12bfc069 65fVHF(0)
3a219f60 66{
9d6d35b0 67 // Standard constructor
3a219f60 68
69 // Output slot #1 writes into a TList container
70 DefineOutput(1,TList::Class()); //My private output
10bdd1ae 71 // Output slot #2 writes into a TNtuple container
72 DefineOutput(2,TNtuple::Class()); //My private output
3a219f60 73}
74
75//________________________________________________________________________
76AliAnalysisTaskSECompareHF::~AliAnalysisTaskSECompareHF()
77{
78 // Destructor
79 if (fOutput) {
80 delete fOutput;
81 fOutput = 0;
82 }
12bfc069 83 if (fVHF) {
84 delete fVHF;
85 fVHF = 0;
86 }
3a219f60 87}
88
89//________________________________________________________________________
90void AliAnalysisTaskSECompareHF::Init()
91{
92 // Initialization
93
94 if(fDebug > 1) printf("AnalysisTaskSECompareHF::Init() \n");
9d6d35b0 95
12bfc069 96 gROOT->LoadMacro("ConfigVertexingHF.C");
97
98 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
99 fVHF->PrintStatus();
9d6d35b0 100
3a219f60 101 return;
102}
103
104//________________________________________________________________________
105void AliAnalysisTaskSECompareHF::UserCreateOutputObjects()
106{
107 // Create the output container
108 //
109 if(fDebug > 1) printf("AnalysisTaskSECompareHF::UserCreateOutputObjects() \n");
110
111 // Several histograms are more conveniently managed in a TList
112 fOutput = new TList();
113 fOutput->SetOwner();
114
115 fHistMass = new TH1F("fHistMass", "D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
116 fHistMass->Sumw2();
117 fHistMass->SetMinimum(0);
118 fOutput->Add(fHistMass);
119
a273622f 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);
124
2a5c77c6 125 OpenFile(2); // 2 is the slot number of the ntuple
607719af 126 fNtupleCmp = new TNtuple("fNtupleCmp","Charm comparison","pdg:nprongs:VxRec:VxTrue:ErrVx:VyRec:VyTrue:ErrVy:VzRec:VzTrue:ErrVz:Chi2toNDF:PtRec:Mrec:CPta:Prodd0");
3a219f60 127
128 return;
129}
130
131//________________________________________________________________________
132void AliAnalysisTaskSECompareHF::UserExec(Option_t */*option*/)
133{
134 // Execute analysis for current event:
135 // heavy flavor candidates association to MC truth
9d6d35b0 136
3a219f60 137
138 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
9d6d35b0 139
a273622f 140 fHistNEvents->Fill(0); // count event
141 // Post the data already here
142 PostData(1,fOutput);
2a041947 143
144 // load HF vertices
145 TClonesArray *inputArrayVertices =
146 (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
147 if(!inputArrayVertices) {
148 printf("AliAnalysisTaskSECompareHF::UserExec: Vertices branch not found!\n");
149 return;
150 }
151
3a219f60 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");
157 return;
158 }
159
9d6d35b0 160
b3999999 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");
166 return;
167 }
9d6d35b0 168
b3999999 169
3a219f60 170 // AOD primary vertex
171 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
172 //vtx1->Print();
173
174 // load MC particles
175 TClonesArray *mcArray =
176 (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
177 if(!mcArray) {
178 printf("AliAnalysisTaskSECompareHF::UserExec: MC particles branch not found!\n");
179 return;
180 }
181
182 // load MC header
183 AliAODMCHeader *mcHeader =
184 (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
185 if(!mcHeader) {
186 printf("AliAnalysisTaskSECompareHF::UserExec: MC header branch not found!\n");
187 return;
188 }
189
2a041947 190 Int_t nprongs,lab,okD0,okD0bar,pdg;
191 Bool_t unsetvtx;
90426506 192 Double_t invmass,posRec[3],posTrue[3],covRec[6],errx,erry,errz;
952f1e33 193 AliAODRecoDecayHF2Prong *d2=0;
194 AliAODRecoDecayHF3Prong *d3=0;
195 AliAODRecoDecayHF4Prong *d4=0;
2a041947 196
c8112d2f 197 Int_t pdgDgD0toKpi[2]={321,211};
198 Int_t pdgDgDplustoKpipi[3]={321,211,211};
199 Int_t pdgDgD0toKpipipi[4]={321,211,211,211};
200
2a041947 201 // loop over vertices
202 Int_t nVertices = inputArrayVertices->GetEntriesFast();
3b1598da 203 if(fDebug>1) printf("Number of vertices: %d\n",nVertices);
2a041947 204
205 for (Int_t iVtx = 0; iVtx < nVertices; iVtx++) {
206 AliAODVertex *vtx = (AliAODVertex*)inputArrayVertices->UncheckedAt(iVtx);
207
90426506 208 vtx->GetXYZ(posRec);
209 vtx->GetCovarianceMatrix(covRec);
210 errx=1.; erry=1.; errz=1.;
211 if(covRec[0]>0) errx = TMath::Sqrt(covRec[0]);
212 if(covRec[2]>0) erry = TMath::Sqrt(covRec[2]);
213 if(covRec[5]>0) errz = TMath::Sqrt(covRec[5]);
214
215
2a041947 216 // get parent AliAODRecoDecayHF
217 nprongs= vtx->GetNDaughters();
607719af 218 // check that the daughters have kITSrefit and kTPCrefit
219 Bool_t allDgOK=kTRUE;
220 for(Int_t idg=0; idg<nprongs; idg++) {
221 AliAODTrack *track = (AliAODTrack*)vtx->GetDaughter(idg);
222 if(!(track->GetStatus()&AliESDtrack::kITSrefit)) allDgOK=kFALSE;
223 if(!(track->GetStatus()&AliESDtrack::kTPCrefit)) allDgOK=kFALSE;
224 }
225 if(!allDgOK) continue;
226
2a041947 227
228 switch(nprongs) {
229 case 2: // look for D0->Kpi
952f1e33 230 d2 = (AliAODRecoDecayHF2Prong*)vtx->GetParent();
27917274 231 if(d2->IsLikeSign()) continue;
a273622f 232 if(d2->Charge() != 0) continue; // these are D*
c8112d2f 233 lab = d2->MatchToMC(421,mcArray,2,pdgDgD0toKpi);
2a041947 234 if(lab>=0) {
235 unsetvtx=kFALSE;
236 if(!d2->GetOwnPrimaryVtx()) {
237 d2->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
238 unsetvtx=kTRUE;
239 }
240 okD0=0; okD0bar=0;
607719af 241 if(d2->SelectD0(fVHF->GetD0toKpiCuts(),okD0,okD0bar)) { // beware: cuts may bias the resolution!
242 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
243 pdg = dMC->GetPdgCode();
244 invmass = (pdg==421 ? d2->InvMassD0() : d2->InvMassD0bar());
607719af 245 // get a daughter for true pos of decay vertex
246 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
247 dg0MC->XvYvZv(posTrue);
caf4ca8b 248 fHistMass->Fill(invmass);
249 // Post the data already here
250 PostData(1,fOutput);
607719af 251 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
252 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
253 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
254 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
255 (Float_t)d2->GetReducedChi2(),(Float_t)d2->Pt(),(Float_t)invmass,
256 (Float_t)d2->CosPointingAngle(),(Float_t)d2->Prodd0d0()};
257 fNtupleCmp->Fill(tmp);
258 PostData(2,fNtupleCmp);
259 }
2a041947 260 if(unsetvtx) d2->UnsetOwnPrimaryVtx();
3a219f60 261 }
2a041947 262 break;
263 case 3: // look for D+
952f1e33 264 d3 = (AliAODRecoDecayHF3Prong*)vtx->GetParent();
27917274 265 if(d3->IsLikeSign()) continue;
c8112d2f 266 lab = d3->MatchToMC(411,mcArray,3,pdgDgDplustoKpipi);
2a041947 267 if(lab>=0) {
268 unsetvtx=kFALSE;
269 if(!d3->GetOwnPrimaryVtx()) {
270 d3->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
271 unsetvtx=kTRUE;
272 }
90426506 273 //if(d3->SelectDplus(fVHF->GetDplusCuts())) {
274 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
275 pdg = dMC->GetPdgCode();
276 invmass = d3->InvMassDplus();
277 // get a daughter for true pos of decay vertex
278 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
279 dg0MC->XvYvZv(posTrue);
607719af 280 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
281 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
282 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
283 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
284 (Float_t)d3->GetReducedChi2(),(Float_t)d3->Pt(),(Float_t)invmass,
285 (Float_t)d3->CosPointingAngle(),(Float_t)(d3->Getd0Prong(0)*d3->Getd0Prong(1)*d3->Getd0Prong(2))};
286 fNtupleCmp->Fill(tmp);
a7c1934a 287 PostData(2,fNtupleCmp);
90426506 288 //}
2a041947 289 if(unsetvtx) d3->UnsetOwnPrimaryVtx();
290 }
291 break;
292 case 4: // look for D0->Kpipipi
952f1e33 293 d4 = (AliAODRecoDecayHF4Prong*)vtx->GetParent();
27917274 294 if(d4->IsLikeSign()) continue;
c8112d2f 295 lab = d4->MatchToMC(421,mcArray,4,pdgDgD0toKpipipi);
2a041947 296 if(lab>=0) {
297 unsetvtx=kFALSE;
298 if(!d4->GetOwnPrimaryVtx()) {
299 d4->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
300 unsetvtx=kTRUE;
301 }
302 okD0=0; okD0bar=0;
90426506 303 //if(d4->SelectD0(fVHF->GetD0to4ProngsCuts(),okD0,okD0bar)) {
304 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
305 pdg = dMC->GetPdgCode();
306 //invmass = (pdg==421 ? d->InvMassD0() : d->InvMassD0bar());
307 invmass = 10.; // dummy
308 // get a daughter for true pos of decay vertex
309 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
310 dg0MC->XvYvZv(posTrue);
607719af 311 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
312 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
313 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
314 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
315 (Float_t)d4->GetReducedChi2(),(Float_t)d4->Pt(),(Float_t)invmass,
316 (Float_t)d4->CosPointingAngle(),(Float_t)(d4->Getd0Prong(0)*d4->Getd0Prong(1)*d4->Getd0Prong(2)*d4->Getd0Prong(3))};
317 fNtupleCmp->Fill(tmp);
a7c1934a 318 PostData(2,fNtupleCmp);
90426506 319 //}
2a041947 320 if(unsetvtx) d4->UnsetOwnPrimaryVtx();
321 }
322 break;
3a219f60 323 }
3a219f60 324
2a041947 325 } // end loop on vertices
3a219f60 326
9d6d35b0 327
b3999999 328 // loop over D*+ candidates
4dd17b11 329 /*
b3999999 330 for (Int_t iDstar = 0; iDstar < inputArrayDstar->GetEntries(); iDstar++) {
331 AliAODRecoCascadeHF *c = (AliAODRecoCascadeHF*)inputArrayDstar->UncheckedAt(iDstar);
332 Int_t labDstar = c->MatchToMC(413,421,mcArray);
333 if(labDstar>=0) printf("GOOD MATCH FOR D*+\n");
334 }
4dd17b11 335 */
2a041947 336
4dd17b11 337 // Post the data already here
338 PostData(1,fOutput);
2a041947 339 PostData(2,fNtupleCmp);
b3999999 340
3a219f60 341 return;
342}
3a219f60 343//________________________________________________________________________
344void AliAnalysisTaskSECompareHF::Terminate(Option_t */*option*/)
345{
346 // Terminate analysis
347 //
348 if(fDebug > 1) printf("AnalysisTaskSECompareHF: Terminate() \n");
349
350 fOutput = dynamic_cast<TList*> (GetOutputData(1));
351 if (!fOutput) {
352 printf("ERROR: fOutput not available\n");
353 return;
354 }
355
3a219f60 356 fHistMass = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMass"));
a273622f 357 fHistNEvents = dynamic_cast<TH1F*>(fOutput->FindObject("fHistNEvents"));
3a219f60 358
27917274 359 //fNtupleCmp = dynamic_cast<TNtuple*> (GetOutputData(2));
10bdd1ae 360
3a219f60 361 return;
362}
363