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