]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliD0toKpiAnalysis.cxx
Primary vertex reconstruction and standalone ITS tracking in the reconstruction chain
[u/mrichter/AliRoot.git] / ANALYSIS / AliD0toKpiAnalysis.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //    Implementation of the D0toKpi reconstruction and analysis class
18 //
19 // Note: the two decay tracks are labelled: 0 (positive track)
20 //                                          1 (negative track)
21 //
22 //            Origin: A. Dainese    andrea.dainese@pd.infn.it            
23 //----------------------------------------------------------------------------
24 #include <Riostream.h>
25 #include <TH1.h>
26 #include <TH2.h>
27 #include <TKey.h>
28 #include <TFile.h>
29 #include <TTree.h>
30 #include <TString.h>
31 #include <TSystem.h>
32 #include <TParticle.h>
33 #include "AliESD.h"
34 #include "AliStack.h"
35 #include "AliRunLoader.h"
36 #include "AliITStrackV2.h"
37 #include "AliITSVertexerTracks.h"
38 #include "AliESDVertex.h"
39 #include "AliV0vertexer.h"
40 #include "AliV0vertex.h"
41 #include "AliD0toKpi.h"
42 #include "AliD0toKpiAnalysis.h"
43
44 typedef struct {
45   Int_t lab;
46   Int_t pdg;
47   Int_t mumlab;
48   Int_t mumpdg;
49   Float_t Vx,Vy,Vz;
50   Float_t Px,Py,Pz;
51 } REFTRACK;
52
53 ClassImp(AliD0toKpiAnalysis)
54
55 //----------------------------------------------------------------------------
56 AliD0toKpiAnalysis::AliD0toKpiAnalysis() {
57   // Default constructor
58
59   SetBz();
60   SetPtCut();
61   Setd0Cut();
62   SetMassCut();
63   SetD0Cuts();
64   SetVertex1();
65   SetPID();
66   fVertexOnTheFly = kFALSE;
67   fSim = kFALSE;
68   fOnlySignal = kFALSE;
69   fDebug = kFALSE;
70 }
71 //----------------------------------------------------------------------------
72 AliD0toKpiAnalysis::~AliD0toKpiAnalysis() {}
73 //----------------------------------------------------------------------------
74 void AliD0toKpiAnalysis::ApplySelection(const Char_t *inName,const Char_t *outName) const {
75   // select candidates that pass fD0Cuts and write them to a new file
76
77   TFile *inFile = TFile::Open(inName);
78
79   TTree *treeD0in=(TTree*)inFile->Get("TreeD0");
80   AliD0toKpiAnalysis *inAnalysis = (AliD0toKpiAnalysis*)inFile->Get("D0toKpiAnalysis");
81   printf("+++\n+++  I N P U T   S T A T U S:\n+++\n");
82   inAnalysis->PrintStatus();
83
84
85   AliD0toKpi *d = 0; 
86   treeD0in->SetBranchAddress("D0toKpi",&d);
87   Int_t entries = (Int_t)treeD0in->GetEntries();
88
89   printf("+++\n+++ Number of D0 in input tree:  %d\n+++\n",entries);
90
91   TTree *treeD0out = new TTree("TreeD0","Tree with selected D0 candidates");
92   treeD0out->Branch("D0toKpi","AliD0toKpi",&d,200000,0);
93
94
95   Int_t okD0=0,okD0bar=0;
96   Int_t nSel = 0;
97
98   for(Int_t i=0; i<entries; i++) {
99     // get event from tree
100     treeD0in->GetEvent(i);
101
102     if(fSim && fOnlySignal && !d->IsSignal()) continue; 
103
104     // check if candidate passes selection (as D0 or D0bar)
105     if(d->Select(fD0Cuts,okD0,okD0bar)) {
106       nSel++;
107       treeD0out->Fill();
108     }
109
110   }
111
112   AliD0toKpiAnalysis *outAnalysis = (AliD0toKpiAnalysis*)inAnalysis->Clone("D0toKpiAnalysis");
113   outAnalysis->SetD0Cuts(fD0Cuts);
114   printf("------------------------------------------\n");
115   printf("+++\n+++  O U T P U T   S T A T U S:\n+++\n");
116   outAnalysis->PrintStatus();
117
118   printf("+++\n+++ Number of D0 in output tree:  %d\n+++\n",nSel);
119
120   TFile* outFile = new TFile(outName,"recreate");
121   treeD0out->Write();
122   outAnalysis->Write();
123   outFile->Close();
124
125   return;
126 }
127 //----------------------------------------------------------------------------
128 Double_t AliD0toKpiAnalysis::CalculateTOFmass(Double_t mom,Double_t length,
129                                               Double_t time) const {
130   // calculated the mass from momentum, track length from vertex to TOF
131   // and time measured by the TOF
132   if(length==0.) return -1000.;
133   Double_t a = time*time/length/length;
134   if(a > 1.) {
135     a = TMath::Sqrt(a-1.);
136   } else {
137     a = -TMath::Sqrt(1.-a);
138   }
139
140   return mom*a;
141 }
142 //----------------------------------------------------------------------------
143 void AliD0toKpiAnalysis::FindCandidates(Int_t evFirst,Int_t evLast,
144                                         const Char_t *outName) {
145   // Find D0 candidates and calculate parameters
146
147   if(fBz<-9000.) {
148     printf("AliD0toKpiAnalysis::FindCandidates():  Set B!\n");
149     return;
150   }
151   AliKalmanTrack::SetConvConst(100/0.299792458/fBz);
152
153   TString trkName("AliITStracksV2.root");
154   if(gSystem->AccessPathName(trkName.Data(),kFileExists)) {
155     printf("AliD0toKpiAnalysis::FindCandidates(): No tracks file found!\n"); 
156     return;
157   }
158
159   TString outName1=outName;
160   TString outName2("nTotEvents.dat");
161
162   Int_t    ev;
163   Int_t    nTotEv=0,nD0rec=0,nD0rec1ev=0;
164   Double_t dca;
165   Double_t v2[3],mom[6],d0[2];
166   Double_t alphaP,alphaN,ptP,ptN,phiP,phiN;
167   Int_t    iTrkP,iTrkN,trkEntries;
168   Int_t    nTrksP=0,nTrksN=0;
169   Int_t    trkNum[2];
170   Int_t    okD0=0,okD0bar=0;
171   Char_t   trkTreeName[100],vtxName[100];
172   AliITStrackV2 *postrack = 0;
173   AliITStrackV2 *negtrack = 0;
174
175   // create the AliITSVertexerTracks object
176   // (it will be used only if fVertexOnTheFly=kTrue)
177   AliITSVertexerTracks *vertexer1 = new AliITSVertexerTracks;
178   vertexer1->SetMinTracks(2);
179   vertexer1->SetDebug(0);
180   Int_t  skipped[2];
181   Bool_t goodVtx1;
182   
183
184   // define the cuts for vertexing
185   Double_t vtxcuts[]={50., // max. allowed chi2
186                       0.0, // min. allowed negative daughter's impact param 
187                       0.0, // min. allowed positive daughter's impact param 
188                       1.0, // max. allowed DCA between the daughter tracks
189                      -1.0, // min. allowed cosine of pointing angle
190                       0.0, // min. radius of the fiducial volume
191                       2.9};// max. radius of the fiducial volume
192   
193   // create the AliV0vertexer object
194   AliV0vertexer *vertexer2 = new AliV0vertexer(vtxcuts);
195
196   // create tree for reconstructed D0s
197   AliD0toKpi *ioD0toKpi=0;
198   TTree *treeD0 = new TTree("TreeD0","Tree with D0 candidates");
199   treeD0->Branch("D0toKpi","AliD0toKpi",&ioD0toKpi,200000,0);
200
201   // open file with tracks
202   TFile *trkFile = TFile::Open(trkName.Data());
203
204   // loop on events in file
205   for(ev=evFirst; ev<=evLast; ev++) {
206     printf(" --- Looking for D0->Kpi in event  %d ---\n",ev);
207
208     // retrieve primary vertex from file
209     sprintf(vtxName,"VertexTracks_%d",ev);
210     AliESDVertex *vertex1stored = (AliESDVertex*)trkFile->Get(vtxName);
211     vertex1stored->GetXYZ(fV1);
212     delete vertex1stored;
213
214     // retrieve tracks from file
215     sprintf(trkTreeName,"TreeT_ITS_%d",ev);
216     TTree *trkTree=(TTree*)trkFile->Get(trkTreeName);
217     if(!trkTree) { 
218       printf("AliD0toKpiAnalysis::FindCandidates():\n  tracks tree not found for evet %d\n",ev); 
219       continue; 
220     }
221     trkEntries = (Int_t)trkTree->GetEntries();
222     printf(" Number of tracks: %d\n",trkEntries);
223
224     // count the total number of events
225     nTotEv++;
226
227     // call function which applies sigle-track selection and
228     // separetes positives and negatives
229     TObjArray trksP(trkEntries/2);
230     Int_t *trkEntryP = new Int_t[trkEntries];
231     TObjArray trksN(trkEntries/2);
232     Int_t *trkEntryN = new Int_t[trkEntries];
233     SelectTracks(*trkTree,trksP,trkEntryP,nTrksP,trksN,trkEntryN,nTrksN);
234       
235     nD0rec1ev = 0;
236
237     // loop on positive tracks
238     for(iTrkP=0; iTrkP<nTrksP; iTrkP++) {
239       if(iTrkP%10==0) printf("  Processing positive track number %d of %d\n",iTrkP,nTrksP);
240           
241       // get track from track array
242       postrack = (AliITStrackV2*)trksP.At(iTrkP);
243       trkNum[0] = trkEntryP[iTrkP];      
244
245       // loop on negative tracks 
246       for(iTrkN=0; iTrkN<nTrksN; iTrkN++) {
247         // get track from tracks array
248         negtrack = (AliITStrackV2*)trksN.At(iTrkN);
249         trkNum[1] = trkEntryN[iTrkN];      
250
251         AliITStrackV2 nt(*negtrack), pt(*postrack), *pnt=&nt, *ppt=&pt;
252
253         //
254         // ----------- DCA MINIMIZATION ------------------
255         //
256         // find the DCA and propagate the tracks to the DCA 
257         dca = vertexer2->PropagateToDCA(pnt,ppt);
258
259         // define the AliV0vertex object
260         AliV0vertex *vertex2 = new AliV0vertex(*pnt,*ppt);
261           
262         // get position of the secondary vertex
263         vertex2->GetXYZ(v2[0],v2[1],v2[2]);
264         
265         delete vertex2;
266   
267         // momenta of the tracks at the vertex
268         ptP = 1./TMath::Abs(ppt->Get1Pt());
269         alphaP = ppt->GetAlpha();
270         phiP = alphaP+TMath::ASin(ppt->GetSnp());
271         mom[0] = ptP*TMath::Cos(phiP); 
272         mom[1] = ptP*TMath::Sin(phiP);
273         mom[2] = ptP*ppt->GetTgl();
274           
275         ptN = 1./TMath::Abs(pnt->Get1Pt());
276         alphaN = pnt->GetAlpha();
277         phiN = alphaN+TMath::ASin(pnt->GetSnp());
278         mom[3] = ptN*TMath::Cos(phiN); 
279         mom[4] = ptN*TMath::Sin(phiN);
280         mom[5] = ptN*pnt->GetTgl();
281           
282         goodVtx1 = kTRUE;
283         // no vertexing if DeltaMass > fMassCut 
284         if(fVertexOnTheFly) {
285           goodVtx1 = kFALSE;
286           if(SelectInvMass(mom)) {
287             // primary vertex from *other* tracks in event
288             vertexer1->SetVtxStart(fV1[0],fV1[1]);
289             skipped[0] = trkEntryP[iTrkP];
290             skipped[1] = trkEntryN[iTrkN];
291             vertexer1->SetSkipTracks(2,skipped);
292             AliESDVertex *vertex1onfly = 
293               (AliESDVertex*)vertexer1->VertexOnTheFly(*trkTree); 
294             if(vertex1onfly->GetNContributors()>0) goodVtx1 = kTRUE;
295             vertex1onfly->GetXYZ(fV1);
296             //vertex1onfly->PrintStatus();
297             delete vertex1onfly;        
298           }
299         }         
300
301         // impact parameters of the tracks w.r.t. the primary vertex
302         d0[0] =  10000.*ppt->GetD(fV1[0],fV1[1]);
303         d0[1] = -10000.*pnt->GetD(fV1[0],fV1[1]);
304
305         // create the object AliD0toKpi
306         AliD0toKpi theD0(ev,trkNum,fV1,v2,dca,mom,d0);
307
308         // select D0s
309         if(goodVtx1 && theD0.Select(fD0Cuts,okD0,okD0bar)) {
310               
311           // fill the tree
312           ioD0toKpi=&theD0;
313           treeD0->Fill();
314
315           nD0rec++; nD0rec1ev++;
316           ioD0toKpi=0;  
317         }
318         
319         negtrack = 0;
320       } // loop on negative tracks
321       postrack = 0;
322     }   // loop on positive tracks
323
324     trksP.Delete();
325     trksN.Delete();
326     delete [] trkEntryP;
327     delete [] trkEntryN;
328     delete trkTree;
329
330     printf(" Number of D0 candidates: %d\n",nD0rec1ev);
331   }    // loop on events in file
332
333
334   printf("\n+++\n+++ Total number of events: %d\n+++\n",nTotEv);
335   printf("\n+++\n+++ Total number of D0 candidates: %d\n+++\n",nD0rec);
336
337   delete vertexer1;
338   delete vertexer2;
339
340   trkFile->Close();
341
342   // create a copy of this class to be written to output file
343   AliD0toKpiAnalysis *copy = (AliD0toKpiAnalysis*)this->Clone("D0toKpiAnalysis");
344   copy->PrintStatus();
345
346   // add PDG codes to decay tracks in found candidates (in simulation mode)
347   // and store tree in the output file
348   if(!fSim) {
349     TFile *outroot = new TFile(outName1.Data(),"recreate");
350     treeD0->Write();
351     copy->Write();
352     outroot->Close();
353     delete outroot;
354   } else {
355     printf(" Now adding information from simulation (PDG codes) ...\n");
356     TTree *treeD0sim = new TTree("TreeD0","Tree with D0 candidates");
357     MakeTracksRefFile(evFirst,evLast);
358     SimulationInfo(treeD0,treeD0sim);
359     delete treeD0;
360     TFile *outroot = new TFile(outName1.Data(),"recreate");
361     treeD0sim->Write();
362     copy->Write();
363     outroot->Close();
364     delete outroot;
365   }
366
367   // write to a file the total number of events
368   FILE *outdat = fopen(outName2.Data(),"w");
369   fprintf(outdat,"%d\n",nTotEv);
370   fclose(outdat);
371
372   return;
373 }
374 //----------------------------------------------------------------------------
375 void AliD0toKpiAnalysis::FindCandidatesESD(Int_t evFirst,Int_t evLast,
376                                            const Char_t *outName) {
377   // Find D0 candidates and calculate parameters
378
379   if(fBz<-9000.) {
380     printf("AliD0toKpiAnalysis::FindCandidatesESD():  Set B!\n");
381     return;
382   }
383   AliKalmanTrack::SetConvConst(100/0.299792458/fBz);
384
385   TString esdName("AliESDs.root");
386   if(gSystem->AccessPathName(esdName.Data(),kFileExists)) {
387     printf("AliD0toKpiAnalysis::FindCandidatesESD(): No ESDs file found!\n"); 
388     return;
389   }
390
391   TString outName1=outName;
392   TString outName2("nTotEvents.dat");
393
394   Int_t    nTotEv=0,nD0rec=0,nD0rec1ev=0;
395   Double_t dca;
396   Double_t v2[3],mom[6],d0[2];
397   Double_t alphaP,alphaN,ptP,ptN,phiP,phiN;
398   Int_t    iTrkP,iTrkN,trkEntries;
399   Int_t    nTrksP=0,nTrksN=0;
400   Int_t    trkNum[2];
401   Double_t tofmass[2];
402   Int_t    okD0=0,okD0bar=0;
403   AliITStrackV2 *postrack = 0;
404   AliITStrackV2 *negtrack = 0;
405
406   // create the AliITSVertexerTracks object
407   // (it will be used only if fVertexOnTheFly=kTrue)
408   AliITSVertexerTracks *vertexer1 = new AliITSVertexerTracks;
409   vertexer1->SetMinTracks(2);
410   vertexer1->SetDebug(0);
411   Int_t  skipped[2];
412   Bool_t goodVtx1;
413   
414
415   // define the cuts for vertexing
416   Double_t vtxcuts[]={50., // max. allowed chi2
417                       0.0, // min. allowed negative daughter's impact param 
418                       0.0, // min. allowed positive daughter's impact param 
419                       1.0, // max. allowed DCA between the daughter tracks
420                      -1.0, // min. allowed cosine of pointing angle
421                       0.0, // min. radius of the fiducial volume
422                       2.9};// max. radius of the fiducial volume
423   
424   // create the AliV0vertexer object
425   AliV0vertexer *vertexer2 = new AliV0vertexer(vtxcuts);
426
427   // create tree for reconstructed D0s
428   AliD0toKpi *ioD0toKpi=0;
429   TTree *treeD0 = new TTree("TreeD0","Tree with D0 candidates");
430   treeD0->Branch("D0toKpi","AliD0toKpi",&ioD0toKpi,200000,0);
431
432   // open file with tracks
433   TFile *esdFile = TFile::Open(esdName.Data());
434
435   TKey *key=0;
436   TIter next(esdFile->GetListOfKeys());
437   // loop on events in file
438   while ((key=(TKey*)next())!=0) {
439     AliESD *event=(AliESD*)key->ReadObj();
440     Int_t ev = (Int_t)event->GetEventNumber();
441     if(ev<evFirst || ev>evLast) { delete event; continue; }
442     printf("--- Finding D0 -> Kpi in event %d\n",ev);
443
444     // retrieve primary vertex from file
445     //sprintf(vtxName,"Vertex_%d",ev);
446     //AliESDVertex *vertex1stored = (AliESDVertex*)trkFile->Get(vtxName);
447     //vertex1stored->GetXYZ(fV1);
448     //delete vertex1stored;
449     event->GetVertex()->GetXYZ(fV1);
450
451     trkEntries = event->GetNumberOfTracks();
452     printf(" Number of tracks: %d\n",trkEntries);
453
454     // count the total number of events
455     nTotEv++;
456
457     // call function which applies sigle-track selection and
458     // separetes positives and negatives
459     TObjArray trksP(trkEntries/2);
460     Int_t    *trkEntryP   = new Int_t[trkEntries];
461     TObjArray trksN(trkEntries/2);
462     Int_t    *trkEntryN   = new Int_t[trkEntries];
463     TTree *trkTree = new TTree();
464     if(fVertexOnTheFly) {
465       SelectTracksESDvtx(*event,trkTree,trksP,trkEntryP,nTrksP,
466                                         trksN,trkEntryN,nTrksN);
467     } else {
468       SelectTracksESD(*event,trksP,trkEntryP,nTrksP,
469                              trksN,trkEntryN,nTrksN);
470     }      
471
472     if(fDebug) printf(" pos. tracks: %d    neg .tracks: %d\n",nTrksP,nTrksN);
473
474     nD0rec1ev = 0;
475
476     // loop on positive tracks
477     for(iTrkP=0; iTrkP<nTrksP; iTrkP++) {
478       if((iTrkP%10==0) || fDebug) printf("  Processing positive track number %d of %d\n",iTrkP,nTrksP);
479           
480       // get track from track array
481       postrack = (AliITStrackV2*)trksP.UncheckedAt(iTrkP);
482       trkNum[0] = trkEntryP[iTrkP];      
483
484       // loop on negative tracks 
485       for(iTrkN=0; iTrkN<nTrksN; iTrkN++) {
486
487         // get track from tracks array
488         negtrack = (AliITStrackV2*)trksN.UncheckedAt(iTrkN);
489         trkNum[1] = trkEntryN[iTrkN];      
490
491         AliITStrackV2 nt(*negtrack), pt(*postrack), *pnt=&nt, *ppt=&pt;
492
493         //
494         // ----------- DCA MINIMIZATION ------------------
495         //
496         // find the DCA and propagate the tracks to the DCA 
497         dca = vertexer2->PropagateToDCA(pnt,ppt);
498
499         // define the AliV0vertex object
500         AliV0vertex *vertex2 = new AliV0vertex(*pnt,*ppt);
501           
502         // get position of the secondary vertex
503         vertex2->GetXYZ(v2[0],v2[1],v2[2]);
504         
505         delete vertex2;
506   
507         // momenta of the tracks at the vertex
508         ptP = 1./TMath::Abs(ppt->Get1Pt());
509         alphaP = ppt->GetAlpha();
510         phiP = alphaP+TMath::ASin(ppt->GetSnp());
511         mom[0] = ptP*TMath::Cos(phiP); 
512         mom[1] = ptP*TMath::Sin(phiP);
513         mom[2] = ptP*ppt->GetTgl();
514           
515         ptN = 1./TMath::Abs(pnt->Get1Pt());
516         alphaN = pnt->GetAlpha();
517         phiN = alphaN+TMath::ASin(pnt->GetSnp());
518         mom[3] = ptN*TMath::Cos(phiN); 
519         mom[4] = ptN*TMath::Sin(phiN);
520         mom[5] = ptN*pnt->GetTgl();
521           
522         goodVtx1 = kTRUE;
523         // no vertexing if DeltaMass > fMassCut 
524         if(fVertexOnTheFly) {
525           goodVtx1 = kFALSE;
526           if(SelectInvMass(mom)) {
527             // primary vertex from *other* tracks in event
528             vertexer1->SetVtxStart(fV1[0],fV1[1]);
529             skipped[0] = trkEntryP[iTrkP];
530             skipped[1] = trkEntryN[iTrkN];
531             vertexer1->SetSkipTracks(2,skipped);
532             AliESDVertex *vertex1onfly = 
533               (AliESDVertex*)vertexer1->VertexOnTheFly(*trkTree); 
534             if(vertex1onfly->GetNContributors()>0) goodVtx1 = kTRUE;
535             vertex1onfly->GetXYZ(fV1);
536             //vertex1onfly->PrintStatus();
537             delete vertex1onfly;        
538           }
539         }         
540
541         // impact parameters of the tracks w.r.t. the primary vertex
542         d0[0] =  10000.*ppt->GetD(fV1[0],fV1[1]);
543         d0[1] = -10000.*pnt->GetD(fV1[0],fV1[1]);
544
545         // create the object AliD0toKpi
546         AliD0toKpi theD0(ev,trkNum,fV1,v2,dca,mom,d0);
547
548         // select D0s
549         if(goodVtx1 && theD0.Select(fD0Cuts,okD0,okD0bar)) {
550               
551           // get PID info from ESD
552           AliESDtrack *t0 = (AliESDtrack*)event->GetTrack(trkNum[0]);
553           Double_t esdpid0[5];
554           t0->GetESDpid(esdpid0);
555           if(t0->GetStatus()&AliESDtrack::kTOFpid) {
556             tofmass[0] = CalculateTOFmass(t0->GetP(),
557                                           t0->GetIntegratedLength(),
558                                           t0->GetTOFsignal());
559           } else {
560             tofmass[0] = -1000.;
561           }
562           AliESDtrack *t1 = (AliESDtrack*)event->GetTrack(trkNum[1]);
563           Double_t esdpid1[5];
564           t1->GetESDpid(esdpid1);
565           if(t1->GetStatus()&AliESDtrack::kTOFpid) {
566             tofmass[1] = CalculateTOFmass(t1->GetP(),
567                                           t1->GetIntegratedLength(),
568                                           t1->GetTOFsignal());
569           } else {
570             tofmass[1] = -1000.;
571           }
572
573           theD0.SetPIDresponse(esdpid0,esdpid1);
574           theD0.SetTOFmasses(tofmass);
575
576           // fill the tree
577           ioD0toKpi=&theD0;
578           treeD0->Fill();
579
580           nD0rec++; nD0rec1ev++;
581           ioD0toKpi=0;  
582         }
583         
584         negtrack = 0;
585       } // loop on negative tracks
586       postrack = 0;
587     }   // loop on positive tracks
588
589     trksP.Delete();
590     trksN.Delete();
591     delete [] trkEntryP;
592     delete [] trkEntryN;
593     delete trkTree;
594     delete event;
595
596
597     printf(" Number of D0 candidates: %d\n",nD0rec1ev);
598   }    // loop on events in file
599
600
601   printf("\n+++\n+++ Total number of events: %d\n+++\n",nTotEv);
602   printf("\n+++\n+++ Total number of D0 candidates: %d\n+++\n",nD0rec);
603
604   delete vertexer1;
605   delete vertexer2;
606
607   esdFile->Close();
608
609   // create a copy of this class to be written to output file
610   AliD0toKpiAnalysis *copy = (AliD0toKpiAnalysis*)this->Clone("D0toKpiAnalysis");
611
612   // add PDG codes to decay tracks in found candidates (in simulation mode)
613   // and store tree in the output file
614   if(!fSim) {
615     TFile *outroot = new TFile(outName1.Data(),"recreate");
616     treeD0->Write();
617     copy->Write();
618     outroot->Close();
619     delete outroot;
620   } else {
621     printf(" Now adding information from simulation (PDG codes) ...\n");
622     TTree *treeD0sim = new TTree("TreeD0","Tree with D0 candidates");
623     MakeTracksRefFileESD();
624     SimulationInfo(treeD0,treeD0sim);
625     delete treeD0;
626     TFile *outroot = new TFile(outName1.Data(),"recreate");
627     treeD0sim->Write();
628     copy->Write();
629     outroot->Close();
630     delete outroot;
631   }
632
633   // write to a file the total number of events
634   FILE *outdat = fopen(outName2.Data(),"w");
635   fprintf(outdat,"%d\n",nTotEv);
636   fclose(outdat);
637
638   return;
639 }
640 //-----------------------------------------------------------------------------
641 void AliD0toKpiAnalysis::PrintStatus() const {
642   // Print parameters being used
643
644   printf(" fBz      = %f T\n",fBz);
645   printf("Preselections:\n");
646   printf("    fPtCut   = %f GeV\n",fPtCut);
647   printf("    fd0Cut   = %f micron\n",fd0Cut);
648   printf("    fMassCut = %f GeV\n",fMassCut);
649   if(fVertexOnTheFly) printf("Primary vertex on the fly\n");
650   if(fSim) { 
651     printf("Simulation mode\n");
652     if(fOnlySignal) printf("  Only signal goes to file\n");
653   }
654   printf("Cuts on candidates:\n");
655   printf("    |M-MD0| [GeV]    < %f\n",fD0Cuts[0]);
656   printf("    dca    [micron]  < %f\n",fD0Cuts[1]);
657   printf("    cosThetaStar     < %f\n",fD0Cuts[2]);
658   printf("    pTK     [GeV]    > %f\n",fD0Cuts[3]);
659   printf("    pTpi    [GeV]    > %f\n",fD0Cuts[4]);
660   printf("    |d0K|  [micron]  < %f\n",fD0Cuts[5]);
661   printf("    |d0pi| [micron]  < %f\n",fD0Cuts[6]);
662   printf("    d0d0  [micron^2] < %f\n",fD0Cuts[7]);
663   printf("    cosThetaPoint    > %f\n",fD0Cuts[8]);
664
665   return;
666 }
667 //-----------------------------------------------------------------------------
668 Bool_t AliD0toKpiAnalysis::SelectInvMass(const Double_t p[6]) const {
669   // Apply preselection in the invariant mass of the pair
670
671   Double_t mD0 = 1.8645;
672   Double_t mPi = 0.13957;
673   Double_t mKa = 0.49368;
674
675   Double_t energy[2];
676   Double_t mom2[2],momTot2;
677
678   mom2[0] = p[0]*p[0] + p[1]*p[1] + p[2]*p[2];
679   mom2[1] = p[3]*p[3] + p[4]*p[4] + p[5]*p[5];
680
681   momTot2 = (p[0]+p[3])*(p[0]+p[3])+
682             (p[1]+p[4])*(p[1]+p[4])+
683             (p[2]+p[5])*(p[2]+p[5]);
684
685   // D0 -> K- pi+
686   energy[1] = TMath::Sqrt(mKa*mKa+mom2[1]);
687   energy[0] = TMath::Sqrt(mPi*mPi+mom2[0]);
688
689   Double_t minvD0 = TMath::Sqrt((energy[0]+energy[1])*(energy[0]+energy[1])-momTot2);
690     
691   // D0bar -> K+ pi-
692   energy[0] = TMath::Sqrt(mKa*mKa+mom2[0]);
693   energy[1] = TMath::Sqrt(mPi*mPi+mom2[1]);
694
695   Double_t minvD0bar = TMath::Sqrt((energy[0]+energy[1])*(energy[0]+energy[1])-momTot2);
696
697   if(TMath::Abs(minvD0-mD0)    < fMassCut) return kTRUE;
698   if(TMath::Abs(minvD0bar-mD0) < fMassCut) return kTRUE;
699   return kFALSE;
700 }
701 //-----------------------------------------------------------------------------
702 void AliD0toKpiAnalysis::SelectTracks(TTree &trkTree,
703                   TObjArray &trksP,Int_t *trkEntryP,Int_t &nTrksP,
704                   TObjArray &trksN,Int_t *trkEntryN,Int_t &nTrksN) const {
705   // Create two TObjArrays with positive and negative tracks and 
706   // apply single-track preselection
707
708   nTrksP=0,nTrksN=0;
709  
710   Int_t entr = (Int_t)trkTree.GetEntries();
711
712   // trasfer tracks from tree to arrays
713   for(Int_t i=0; i<entr; i++) {
714
715     AliITStrackV2 *track = new AliITStrackV2; 
716     trkTree.SetBranchAddress("tracks",&track);
717
718     trkTree.GetEvent(i);
719
720     // single track selection
721     if(!SingleTrkCuts(*track)) { delete track; continue; }
722
723     if(track->Get1Pt()>0.) { // negative track
724       trksN.AddLast(track);
725       trkEntryN[nTrksN] = i;
726       nTrksN++;
727     } else {                 // positive track
728       trksP.AddLast(track);
729       trkEntryP[nTrksP] = i;
730       nTrksP++;
731     }
732
733   }
734
735   return;
736 }
737 //-----------------------------------------------------------------------------
738 void AliD0toKpiAnalysis::SelectTracksESD(AliESD &event,
739         TObjArray &trksP,Int_t *trkEntryP,Int_t &nTrksP,
740         TObjArray &trksN,Int_t *trkEntryN,Int_t &nTrksN) const {
741   // Create two TObjArrays with positive and negative tracks and 
742   // apply single-track preselection
743
744   nTrksP=0,nTrksN=0;
745
746   Int_t entr = event.GetNumberOfTracks();
747  
748   // transfer ITS tracks from ESD to arrays and to a tree
749   for(Int_t i=0; i<entr; i++) {
750     //if(fDebug) printf(" SelectTracksESD: %d/%d\n",i,entr);
751
752     AliESDtrack *esdtrack = (AliESDtrack*)event.GetTrack(i);
753     UInt_t status = esdtrack->GetStatus();
754
755     if(!(status&AliESDtrack::kITSrefit)) continue;
756
757     AliITStrackV2 *itstrack = new AliITStrackV2(*esdtrack); 
758
759     // single track selection
760     if(!SingleTrkCuts(*itstrack)) 
761       { delete itstrack; continue; }
762
763     if(itstrack->Get1Pt()>0.) { // negative track
764       trksN.AddLast(itstrack);
765       trkEntryN[nTrksN] = i;
766       nTrksN++;
767     } else {                 // positive track
768       trksP.AddLast(itstrack);
769       trkEntryP[nTrksP] = i;
770       nTrksP++;
771     }
772
773   } // loop on ESD tracks
774
775   return;
776 }
777 //-----------------------------------------------------------------------------
778 void AliD0toKpiAnalysis::SelectTracksESDvtx(AliESD &event,TTree *trkTree,
779         TObjArray &trksP,Int_t *trkEntryP,Int_t &nTrksP,
780         TObjArray &trksN,Int_t *trkEntryN,Int_t &nTrksN) const {
781   // Create two TObjArrays with positive and negative tracks and 
782   // apply single-track preselection
783
784   nTrksP=0,nTrksN=0;
785
786   Int_t entr = event.GetNumberOfTracks();
787  
788   AliITStrackV2 *itstrackfortree = 0;
789   trkTree->Branch("tracks","AliITStrackV2",&itstrackfortree,entr,0);
790
791
792   // transfer ITS tracks from ESD to arrays and to a tree
793   for(Int_t i=0; i<entr; i++) {
794
795     AliESDtrack *esdtrack = (AliESDtrack*)event.GetTrack(i);
796     UInt_t status = esdtrack->GetStatus();
797
798     if(!(status&AliESDtrack::kITSrefit)) continue;
799
800     AliITStrackV2 *itstrack = new AliITStrackV2(*esdtrack);
801
802     // store track in the tree to be used for primary vertex finding
803     itstrackfortree = new AliITStrackV2(*esdtrack);
804     trkTree->Fill();
805     itstrackfortree = 0;
806
807     // single track selection
808     if(!SingleTrkCuts(*itstrack)) 
809       { delete itstrack; continue; }
810
811     if(itstrack->Get1Pt()>0.) { // negative track
812       trksN.AddLast(itstrack);
813       trkEntryN[nTrksN] = i;
814       nTrksN++;
815     } else {                 // positive track
816       trksP.AddLast(itstrack);
817       trkEntryP[nTrksP] = i;
818       nTrksP++;
819     }
820
821   } // loop on esd tracks
822
823   delete itstrackfortree;
824
825   return;
826 }
827 //-----------------------------------------------------------------------------
828 void AliD0toKpiAnalysis::SetD0Cuts(Double_t cut0,Double_t cut1,
829                                    Double_t cut2,Double_t cut3,Double_t cut4,
830                                    Double_t cut5,Double_t cut6,
831                                    Double_t cut7,Double_t cut8) {
832   // Set the cuts for D0 selection
833   fD0Cuts[0] = cut0;
834   fD0Cuts[1] = cut1;
835   fD0Cuts[2] = cut2;
836   fD0Cuts[3] = cut3;
837   fD0Cuts[4] = cut4;
838   fD0Cuts[5] = cut5;
839   fD0Cuts[6] = cut6;
840   fD0Cuts[7] = cut7;
841   fD0Cuts[8] = cut8;
842
843   return;
844 }
845 //-----------------------------------------------------------------------------
846 void AliD0toKpiAnalysis::SetD0Cuts(const Double_t cuts[9]) {
847   // Set the cuts for D0 selection
848
849   for(Int_t i=0; i<9; i++) fD0Cuts[i] = cuts[i];
850
851   return;
852 }
853 //-----------------------------------------------------------------------------
854 Bool_t AliD0toKpiAnalysis::SingleTrkCuts(const AliITStrackV2& trk) const {
855   // Check if track passes some kinematical cuts  
856
857   if(TMath::Abs(1./trk.Get1Pt()) < fPtCut) 
858     return kFALSE;
859   if(TMath::Abs(10000.*trk.GetD(fV1[0],fV1[1])) < fd0Cut) 
860     return kFALSE;
861
862   return kTRUE;
863 }
864 //----------------------------------------------------------------------------
865 void AliD0toKpiAnalysis::MakeTracksRefFile(Int_t evFirst,Int_t evLast) 
866   const {
867   // Create a file with simulation info for the reconstructed tracks
868   
869   TFile *out = TFile::Open("ITStracksRefFile.root","recreate");
870   TFile *trk = TFile::Open("AliITStracksV2.root");
871   AliRunLoader *rl = AliRunLoader::Open("galice.root");
872
873   // load kinematics
874   rl->LoadKinematics();
875   
876   Int_t      label;
877   TParticle *part;  
878   TParticle *mumpart;
879   REFTRACK   reftrk;
880   
881   for(Int_t ev=evFirst; ev<=evLast; ev++){
882     rl->GetEvent(ev);  
883     AliStack *stack = rl->Stack();
884
885     trk->cd();
886
887     // Tree with tracks
888     char tname[100];
889     sprintf(tname,"TreeT_ITS_%d",ev);
890
891     TTree *tracktree=(TTree*)trk->Get(tname);
892     if(!tracktree) continue;
893     AliITStrackV2 *track = new AliITStrackV2; 
894     tracktree->SetBranchAddress("tracks",&track);
895     Int_t nentr=(Int_t)tracktree->GetEntries();
896
897     // Tree for true track parameters
898     char ttname[100];
899     sprintf(ttname,"Tree_Ref_%d",ev);
900     TTree *reftree = new TTree(ttname,"Tree with true track params");
901     reftree->Branch("rectracks",&reftrk,"lab/I:pdg:mumlab:mumpdg:Vx/F:Vy:Vz:Px:Py:Pz");
902
903     for(Int_t i=0; i<nentr; i++) {
904       tracktree->GetEvent(i);
905       label = TMath::Abs(track->GetLabel());
906
907       part = (TParticle*)stack->Particle(label);
908       reftrk.lab = label;
909       reftrk.pdg = part->GetPdgCode();
910       reftrk.mumlab = part->GetFirstMother();
911       if(part->GetFirstMother()>=0) {
912         mumpart = (TParticle*)stack->Particle(part->GetFirstMother());
913         reftrk.mumpdg = mumpart->GetPdgCode();
914       } else {
915         reftrk.mumpdg=-1;
916       }
917       reftrk.Vx = part->Vx();
918       reftrk.Vy = part->Vy();
919       reftrk.Vz = part->Vz();
920       reftrk.Px = part->Px();
921       reftrk.Py = part->Py();
922       reftrk.Pz = part->Pz();
923       
924       reftree->Fill();
925     } // loop on tracks   
926
927     out->cd();
928     reftree->Write();
929
930     delete track;
931     delete reftree;
932     delete tracktree;
933     delete stack;
934   } // loop on events
935
936   trk->Close();
937   out->Close();
938
939   return;
940 }
941 //----------------------------------------------------------------------------
942 void AliD0toKpiAnalysis::MakeTracksRefFileESD() const {
943   // Create a file with simulation info for the reconstructed tracks
944   
945   TFile *outFile = TFile::Open("ITStracksRefFile.root","recreate");
946   TFile *esdFile = TFile::Open("AliESDs.root");
947   AliRunLoader *rl = AliRunLoader::Open("galice.root");
948
949   // load kinematics
950   rl->LoadKinematics();
951   
952   Int_t      label;
953   TParticle *part;  
954   TParticle *mumpart;
955   REFTRACK   reftrk;
956   
957   TKey *key=0;
958   TIter next(esdFile->GetListOfKeys());
959   // loop on events in file
960   while ((key=(TKey*)next())!=0) {
961     AliESD *event=(AliESD*)key->ReadObj();
962     Int_t ev = (Int_t)event->GetEventNumber();
963
964     rl->GetEvent(ev);  
965     AliStack *stack = rl->Stack();
966
967     Int_t nentr=(Int_t)event->GetNumberOfTracks();
968
969     // Tree for true track parameters
970     char ttname[100];
971     sprintf(ttname,"Tree_Ref_%d",ev);
972     TTree *reftree = new TTree(ttname,"Tree with true track params");
973     reftree->Branch("rectracks",&reftrk,"lab/I:pdg:mumlab:mumpdg:Vx/F:Vy:Vz:Px:Py:Pz");
974
975     for(Int_t i=0; i<nentr; i++) {
976       AliESDtrack *esdtrack = (AliESDtrack*)event->GetTrack(i);
977       label = TMath::Abs(esdtrack->GetLabel());
978
979       part = (TParticle*)stack->Particle(label);
980       reftrk.lab = label;
981       reftrk.pdg = part->GetPdgCode();
982       reftrk.mumlab = part->GetFirstMother();
983       if(part->GetFirstMother()>=0) {
984         mumpart = (TParticle*)stack->Particle(part->GetFirstMother());
985         reftrk.mumpdg = mumpart->GetPdgCode();
986       } else {
987         reftrk.mumpdg=-1;
988       }
989       reftrk.Vx = part->Vx();
990       reftrk.Vy = part->Vy();
991       reftrk.Vz = part->Vz();
992       reftrk.Px = part->Px();
993       reftrk.Py = part->Py();
994       reftrk.Pz = part->Pz();
995       
996       reftree->Fill();
997
998     } // loop on tracks   
999
1000     outFile->cd();
1001     reftree->Write();
1002
1003     delete reftree;
1004     delete event;
1005     delete stack;
1006   } // loop on events
1007
1008   esdFile->Close();
1009   outFile->Close();
1010
1011   return;
1012 }
1013 //-----------------------------------------------------------------------------
1014 void AliD0toKpiAnalysis::SimulationInfo(TTree *treeD0in,TTree *treeD0out) const {
1015   // add pdg codes to candidate decay tracks (for sim)
1016
1017   TString refFileName("ITStracksRefFile.root");
1018   if(fSim && gSystem->AccessPathName(refFileName.Data(),kFileExists)) { 
1019     printf("AliD0toKpiAnalysis::SimulationInfo: no reference file found!\n"); 
1020     return;
1021   }
1022   TFile *refFile = TFile::Open(refFileName.Data());
1023
1024   Char_t refTreeName[100];
1025   Int_t  event;
1026   Int_t  pdg[2],mumpdg[2],mumlab[2];
1027   REFTRACK reftrk;
1028
1029   // read-in reference tree for event 0 (the only event for Pb-Pb)
1030   sprintf(refTreeName,"Tree_Ref_%d",0);
1031   TTree *refTree0 = (TTree*)refFile->Get(refTreeName);
1032   refTree0->SetBranchAddress("rectracks",&reftrk);
1033
1034   AliD0toKpi *theD0 = 0; 
1035   treeD0in->SetBranchAddress("D0toKpi",&theD0);
1036   treeD0out->Branch("D0toKpi","AliD0toKpi",&theD0,200000,0);
1037
1038   Int_t entries = (Int_t)treeD0in->GetEntries();
1039
1040   for(Int_t i=0; i<entries; i++) {
1041     if(i%100==0) printf("  done %d candidates of %d\n",i,entries);    
1042
1043     treeD0in->GetEvent(i);
1044     event = theD0->EventNo();
1045
1046     if(event==0) { // always true for Pb-Pb (avoid to read-in tree every time)
1047       refTree0->GetEvent(theD0->GetTrkNum(0));
1048       pdg[0]    = reftrk.pdg;
1049       mumpdg[0] = reftrk.mumpdg;
1050       mumlab[0] = reftrk.mumlab;
1051       refTree0->GetEvent(theD0->GetTrkNum(1));
1052       pdg[1]    = reftrk.pdg;
1053       mumpdg[1] = reftrk.mumpdg;
1054       mumlab[1] = reftrk.mumlab;
1055     } else {
1056       sprintf(refTreeName,"Tree_Ref_%d",event);
1057       TTree *refTree = (TTree*)refFile->Get(refTreeName);
1058       refTree->SetBranchAddress("rectracks",&reftrk);
1059       refTree->GetEvent(theD0->GetTrkNum(0));
1060       pdg[0]    = reftrk.pdg;
1061       mumpdg[0] = reftrk.mumpdg;
1062       mumlab[0] = reftrk.mumlab;
1063       refTree->GetEvent(theD0->GetTrkNum(1));
1064       pdg[1]    = reftrk.pdg;
1065       mumpdg[1] = reftrk.mumpdg;
1066       mumlab[1] = reftrk.mumlab;
1067       delete refTree;
1068     }
1069     
1070     theD0->SetPdgCodes(pdg);
1071     theD0->SetMumPdgCodes(mumpdg);
1072     
1073     if(TMath::Abs(mumpdg[0])==421 && TMath::Abs(mumpdg[1])==421 
1074        && mumlab[0]==mumlab[1]) theD0->SetSignal();
1075     
1076     if(!fOnlySignal || theD0->IsSignal()) treeD0out->Fill();
1077
1078   }
1079
1080   delete refTree0;
1081
1082   refFile->Close();
1083
1084   return;
1085 }
1086
1087
1088
1089
1090
1091