]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/vertexingHF/AliAnalysisTaskSECompareHF.cxx
Added method to add a default set of jet fidners for a deltaAOD
[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
8931c313 140 // In case there is an AOD handler writing a standard AOD, use the AOD
141 // event in memory rather than the input (ESD) event.
142 if (!aod && AODEvent() && IsStandardAOD()) aod = dynamic_cast<AliAODEvent*> (AODEvent());
143
a273622f 144 fHistNEvents->Fill(0); // count event
145 // Post the data already here
146 PostData(1,fOutput);
2a041947 147
148 // load HF vertices
149 TClonesArray *inputArrayVertices =
150 (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
151 if(!inputArrayVertices) {
152 printf("AliAnalysisTaskSECompareHF::UserExec: Vertices branch not found!\n");
153 return;
154 }
155
3a219f60 156 // load D0->Kpi candidates
157 TClonesArray *inputArrayD0toKpi =
158 (TClonesArray*)aod->GetList()->FindObject("D0toKpi");
159 if(!inputArrayD0toKpi) {
160 printf("AliAnalysisTaskSECompareHF::UserExec: D0toKpi branch not found!\n");
161 return;
162 }
163
9d6d35b0 164
b3999999 165 // load D*+ candidates
166 TClonesArray *inputArrayDstar =
167 (TClonesArray*)aod->GetList()->FindObject("Dstar");
168 if(!inputArrayDstar) {
169 printf("AliAnalysisTaskSECompareHF::UserExec: Dstar branch not found!\n");
170 return;
171 }
9d6d35b0 172
b3999999 173
3a219f60 174 // AOD primary vertex
175 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
176 //vtx1->Print();
177
178 // load MC particles
179 TClonesArray *mcArray =
180 (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
181 if(!mcArray) {
182 printf("AliAnalysisTaskSECompareHF::UserExec: MC particles branch not found!\n");
183 return;
184 }
185
186 // load MC header
187 AliAODMCHeader *mcHeader =
188 (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
189 if(!mcHeader) {
190 printf("AliAnalysisTaskSECompareHF::UserExec: MC header branch not found!\n");
191 return;
192 }
193
2a041947 194 Int_t nprongs,lab,okD0,okD0bar,pdg;
195 Bool_t unsetvtx;
90426506 196 Double_t invmass,posRec[3],posTrue[3],covRec[6],errx,erry,errz;
952f1e33 197 AliAODRecoDecayHF2Prong *d2=0;
198 AliAODRecoDecayHF3Prong *d3=0;
199 AliAODRecoDecayHF4Prong *d4=0;
2a041947 200
c8112d2f 201 Int_t pdgDgD0toKpi[2]={321,211};
202 Int_t pdgDgDplustoKpipi[3]={321,211,211};
203 Int_t pdgDgD0toKpipipi[4]={321,211,211,211};
204
2a041947 205 // loop over vertices
206 Int_t nVertices = inputArrayVertices->GetEntriesFast();
3b1598da 207 if(fDebug>1) printf("Number of vertices: %d\n",nVertices);
2a041947 208
209 for (Int_t iVtx = 0; iVtx < nVertices; iVtx++) {
210 AliAODVertex *vtx = (AliAODVertex*)inputArrayVertices->UncheckedAt(iVtx);
211
90426506 212 vtx->GetXYZ(posRec);
213 vtx->GetCovarianceMatrix(covRec);
214 errx=1.; erry=1.; errz=1.;
215 if(covRec[0]>0) errx = TMath::Sqrt(covRec[0]);
216 if(covRec[2]>0) erry = TMath::Sqrt(covRec[2]);
217 if(covRec[5]>0) errz = TMath::Sqrt(covRec[5]);
218
219
2a041947 220 // get parent AliAODRecoDecayHF
221 nprongs= vtx->GetNDaughters();
607719af 222 // check that the daughters have kITSrefit and kTPCrefit
223 Bool_t allDgOK=kTRUE;
224 for(Int_t idg=0; idg<nprongs; idg++) {
225 AliAODTrack *track = (AliAODTrack*)vtx->GetDaughter(idg);
226 if(!(track->GetStatus()&AliESDtrack::kITSrefit)) allDgOK=kFALSE;
227 if(!(track->GetStatus()&AliESDtrack::kTPCrefit)) allDgOK=kFALSE;
228 }
229 if(!allDgOK) continue;
230
2a041947 231
232 switch(nprongs) {
233 case 2: // look for D0->Kpi
952f1e33 234 d2 = (AliAODRecoDecayHF2Prong*)vtx->GetParent();
27917274 235 if(d2->IsLikeSign()) continue;
a273622f 236 if(d2->Charge() != 0) continue; // these are D*
c8112d2f 237 lab = d2->MatchToMC(421,mcArray,2,pdgDgD0toKpi);
2a041947 238 if(lab>=0) {
239 unsetvtx=kFALSE;
240 if(!d2->GetOwnPrimaryVtx()) {
241 d2->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
242 unsetvtx=kTRUE;
243 }
244 okD0=0; okD0bar=0;
607719af 245 if(d2->SelectD0(fVHF->GetD0toKpiCuts(),okD0,okD0bar)) { // beware: cuts may bias the resolution!
246 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
247 pdg = dMC->GetPdgCode();
248 invmass = (pdg==421 ? d2->InvMassD0() : d2->InvMassD0bar());
607719af 249 // get a daughter for true pos of decay vertex
250 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
251 dg0MC->XvYvZv(posTrue);
caf4ca8b 252 fHistMass->Fill(invmass);
253 // Post the data already here
254 PostData(1,fOutput);
607719af 255 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
256 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
257 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
258 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
259 (Float_t)d2->GetReducedChi2(),(Float_t)d2->Pt(),(Float_t)invmass,
260 (Float_t)d2->CosPointingAngle(),(Float_t)d2->Prodd0d0()};
261 fNtupleCmp->Fill(tmp);
262 PostData(2,fNtupleCmp);
263 }
2a041947 264 if(unsetvtx) d2->UnsetOwnPrimaryVtx();
3a219f60 265 }
2a041947 266 break;
267 case 3: // look for D+
952f1e33 268 d3 = (AliAODRecoDecayHF3Prong*)vtx->GetParent();
27917274 269 if(d3->IsLikeSign()) continue;
c8112d2f 270 lab = d3->MatchToMC(411,mcArray,3,pdgDgDplustoKpipi);
2a041947 271 if(lab>=0) {
272 unsetvtx=kFALSE;
273 if(!d3->GetOwnPrimaryVtx()) {
274 d3->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
275 unsetvtx=kTRUE;
276 }
90426506 277 //if(d3->SelectDplus(fVHF->GetDplusCuts())) {
278 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
279 pdg = dMC->GetPdgCode();
280 invmass = d3->InvMassDplus();
281 // get a daughter for true pos of decay vertex
282 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
283 dg0MC->XvYvZv(posTrue);
607719af 284 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
285 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
286 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
287 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
288 (Float_t)d3->GetReducedChi2(),(Float_t)d3->Pt(),(Float_t)invmass,
289 (Float_t)d3->CosPointingAngle(),(Float_t)(d3->Getd0Prong(0)*d3->Getd0Prong(1)*d3->Getd0Prong(2))};
290 fNtupleCmp->Fill(tmp);
a7c1934a 291 PostData(2,fNtupleCmp);
90426506 292 //}
2a041947 293 if(unsetvtx) d3->UnsetOwnPrimaryVtx();
294 }
295 break;
296 case 4: // look for D0->Kpipipi
952f1e33 297 d4 = (AliAODRecoDecayHF4Prong*)vtx->GetParent();
27917274 298 if(d4->IsLikeSign()) continue;
c8112d2f 299 lab = d4->MatchToMC(421,mcArray,4,pdgDgD0toKpipipi);
2a041947 300 if(lab>=0) {
301 unsetvtx=kFALSE;
302 if(!d4->GetOwnPrimaryVtx()) {
303 d4->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
304 unsetvtx=kTRUE;
305 }
306 okD0=0; okD0bar=0;
90426506 307 //if(d4->SelectD0(fVHF->GetD0to4ProngsCuts(),okD0,okD0bar)) {
308 AliAODMCParticle *dMC = (AliAODMCParticle*)mcArray->At(lab);
309 pdg = dMC->GetPdgCode();
310 //invmass = (pdg==421 ? d->InvMassD0() : d->InvMassD0bar());
311 invmass = 10.; // dummy
312 // get a daughter for true pos of decay vertex
313 AliAODMCParticle *dg0MC = (AliAODMCParticle*)mcArray->At(dMC->GetDaughter(0));
314 dg0MC->XvYvZv(posTrue);
607719af 315 Float_t tmp[16]={(Float_t)pdg,(Float_t)nprongs,
316 (Float_t)posRec[0],(Float_t)posTrue[0],(Float_t)errx,
317 (Float_t)posRec[1],(Float_t)posTrue[1],(Float_t)erry,
318 (Float_t)posRec[2],(Float_t)posTrue[2],(Float_t)errz,
319 (Float_t)d4->GetReducedChi2(),(Float_t)d4->Pt(),(Float_t)invmass,
320 (Float_t)d4->CosPointingAngle(),(Float_t)(d4->Getd0Prong(0)*d4->Getd0Prong(1)*d4->Getd0Prong(2)*d4->Getd0Prong(3))};
321 fNtupleCmp->Fill(tmp);
a7c1934a 322 PostData(2,fNtupleCmp);
90426506 323 //}
2a041947 324 if(unsetvtx) d4->UnsetOwnPrimaryVtx();
325 }
326 break;
3a219f60 327 }
3a219f60 328
2a041947 329 } // end loop on vertices
3a219f60 330
9d6d35b0 331
b3999999 332 // loop over D*+ candidates
4dd17b11 333 /*
b3999999 334 for (Int_t iDstar = 0; iDstar < inputArrayDstar->GetEntries(); iDstar++) {
335 AliAODRecoCascadeHF *c = (AliAODRecoCascadeHF*)inputArrayDstar->UncheckedAt(iDstar);
336 Int_t labDstar = c->MatchToMC(413,421,mcArray);
337 if(labDstar>=0) printf("GOOD MATCH FOR D*+\n");
338 }
4dd17b11 339 */
2a041947 340
4dd17b11 341 // Post the data already here
342 PostData(1,fOutput);
2a041947 343 PostData(2,fNtupleCmp);
b3999999 344
3a219f60 345 return;
346}
3a219f60 347//________________________________________________________________________
348void AliAnalysisTaskSECompareHF::Terminate(Option_t */*option*/)
349{
350 // Terminate analysis
351 //
352 if(fDebug > 1) printf("AnalysisTaskSECompareHF: Terminate() \n");
353
354 fOutput = dynamic_cast<TList*> (GetOutputData(1));
355 if (!fOutput) {
356 printf("ERROR: fOutput not available\n");
357 return;
358 }
359
3a219f60 360 fHistMass = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMass"));
a273622f 361 fHistNEvents = dynamic_cast<TH1F*>(fOutput->FindObject("fHistNEvents"));
3a219f60 362
27917274 363 //fNtupleCmp = dynamic_cast<TNtuple*> (GetOutputData(2));
10bdd1ae 364
3a219f60 365 return;
366}
367