]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliAnalysisVertexingHF.cxx
Undo change of rev. 30509, that was committed by mistake
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliAnalysisVertexingHF.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 heavy-flavour vertexing analysis class
18 // Candidates are stored in the AOD as objects deriving from AliAODRecoDecay.
19 // To be used as a task of AliAnalysisManager by means of the interface
20 // class AliAnalysisTaskSEVertexingHF. 
21 // An example of usage in the macro AliAnalysisTaskSEVertexingHFTest.C.
22 //
23 //  Origin: E.Bruna, G.E.Bruno, A.Dainese, F.Prino, R.Romita, X.M.Zhang
24 //  Contact: andrea.dainese@lnl.infn.it
25 //----------------------------------------------------------------------------
26 #include <TFile.h>
27 #include <TDatabasePDG.h>
28 #include <TString.h>
29 #include "AliLog.h"
30 #include "AliVEvent.h"
31 #include "AliVVertex.h"
32 #include "AliVTrack.h"
33 #include "AliVertexerTracks.h"
34 #include "AliKFVertex.h"
35 #include "AliESDEvent.h"
36 #include "AliESDVertex.h"
37 #include "AliESDtrack.h"
38 #include "AliESDtrackCuts.h"
39 #include "AliAODEvent.h"
40 #include "AliAODRecoDecay.h"
41 #include "AliAODRecoDecayHF.h"
42 #include "AliAODRecoDecayHF2Prong.h"
43 #include "AliAODRecoDecayHF3Prong.h"
44 #include "AliAODRecoDecayHF4Prong.h"
45 #include "AliAODRecoCascadeHF.h"
46 #include "AliAnalysisFilter.h"
47 #include "AliAnalysisVertexingHF.h"
48
49 ClassImp(AliAnalysisVertexingHF)
50
51 //----------------------------------------------------------------------------
52 AliAnalysisVertexingHF::AliAnalysisVertexingHF():
53 fInputAOD(kFALSE),
54 fAODMap(0),
55 fBzkG(0.),
56 fSecVtxWithKF(kFALSE),
57 fRecoPrimVtxSkippingTrks(kFALSE),
58 fRmTrksFromPrimVtx(kFALSE),
59 fV1(0x0),
60 fD0toKpi(kTRUE),
61 fJPSItoEle(kTRUE),
62 f3Prong(kTRUE),
63 f4Prong(kTRUE),
64 fDstar(kTRUE),
65 fTrackFilter(0x0),
66 fTrackFilterSoftPi(0x0)
67 {
68   // Default constructor
69
70   SetD0toKpiCuts();
71   SetBtoJPSICuts();
72   SetDplusCuts();
73   SetDsCuts();
74   SetLcCuts();
75   SetDstarCuts();
76
77 }
78 //--------------------------------------------------------------------------
79 AliAnalysisVertexingHF::AliAnalysisVertexingHF(const AliAnalysisVertexingHF &source) : 
80 TNamed(source),
81 fInputAOD(source.fInputAOD),
82 fAODMap(source.fAODMap),
83 fBzkG(source.fBzkG),
84 fSecVtxWithKF(source.fSecVtxWithKF),
85 fRecoPrimVtxSkippingTrks(source.fRecoPrimVtxSkippingTrks),
86 fRmTrksFromPrimVtx(source.fRmTrksFromPrimVtx),
87 fV1(source.fV1),
88 fD0toKpi(source.fD0toKpi),
89 fJPSItoEle(source.fJPSItoEle),
90 f3Prong(source.f3Prong),
91 f4Prong(source.f4Prong),
92 fDstar(source.fDstar),
93 fTrackFilter(source.fTrackFilter),
94 fTrackFilterSoftPi(source.fTrackFilterSoftPi)
95 {
96   //
97   // Copy constructor
98   //
99   for(Int_t i=0; i<9; i++)  fD0toKpiCuts[i]=source.fD0toKpiCuts[i];
100   for(Int_t i=0; i<9; i++)  fBtoJPSICuts[i]=source.fBtoJPSICuts[i];
101   for(Int_t i=0; i<12; i++) fDplusCuts[i]=source.fDplusCuts[i];
102   for(Int_t i=0; i<13; i++) fDsCuts[i]=source.fDsCuts[i];
103   for(Int_t i=0; i<12; i++) fLcCuts[i]=source.fLcCuts[i];
104   for(Int_t i=0; i<5; i++)  fDstarCuts[i]=source.fDstarCuts[i];
105 }
106 //--------------------------------------------------------------------------
107 AliAnalysisVertexingHF &AliAnalysisVertexingHF::operator=(const AliAnalysisVertexingHF &source)
108 {
109   //
110   // assignment operator
111   //
112   if(&source == this) return *this;
113   fInputAOD = source.fInputAOD;
114   fBzkG = source.fBzkG;
115   fSecVtxWithKF = source.fSecVtxWithKF;
116   fRecoPrimVtxSkippingTrks = source.fRecoPrimVtxSkippingTrks;
117   fRmTrksFromPrimVtx = source.fRmTrksFromPrimVtx;
118   fV1 = source.fV1;
119   fD0toKpi = source.fD0toKpi;
120   fJPSItoEle = source.fJPSItoEle;
121   f3Prong = source.f3Prong;
122   f4Prong = source.f4Prong;
123   fDstar = source.fDstar;
124   fTrackFilter = source.fTrackFilter;
125   fTrackFilterSoftPi = source.fTrackFilterSoftPi;
126
127   for(Int_t i=0; i<9; i++)  fD0toKpiCuts[i]=source.fD0toKpiCuts[i];
128   for(Int_t i=0; i<9; i++)  fBtoJPSICuts[i]=source.fBtoJPSICuts[i];
129   for(Int_t i=0; i<12; i++) fDplusCuts[i]=source.fDplusCuts[i];
130   for(Int_t i=0; i<13; i++) fDsCuts[i]=source.fDsCuts[i];
131   for(Int_t i=0; i<12; i++) fLcCuts[i]=source.fLcCuts[i];
132   for(Int_t i=0; i<5; i++)  fDstarCuts[i]=source.fDstarCuts[i];
133
134   return *this;
135 }
136 //----------------------------------------------------------------------------
137 AliAnalysisVertexingHF::~AliAnalysisVertexingHF() {
138   // Destructor
139   if(fV1) { delete fV1; fV1=0; }
140   if(fTrackFilter) { delete fTrackFilter; fTrackFilter=0; }
141   if(fTrackFilterSoftPi) { delete fTrackFilterSoftPi; fTrackFilterSoftPi=0; }
142 }
143 //----------------------------------------------------------------------------
144 void AliAnalysisVertexingHF::FindCandidates(AliVEvent *event,
145                                             TClonesArray *aodVerticesHFTClArr,
146                                             TClonesArray *aodD0toKpiTClArr,
147                                             TClonesArray *aodJPSItoEleTClArr,
148                                             TClonesArray *aodCharm3ProngTClArr,
149                                             TClonesArray *aodCharm4ProngTClArr,
150                                             TClonesArray *aodDstarTClArr)
151 {
152   // Find heavy-flavour vertex candidates
153   // Input:  ESD or AOD
154   // Output: AOD (additional branches added)
155
156
157   TString evtype = event->IsA()->GetName();
158   fInputAOD = ((evtype=="AliAODEvent") ? kTRUE : kFALSE);
159
160   if(!aodVerticesHFTClArr) {
161     printf("ERROR: no aodVerticesHFTClArr");
162     return;
163   }
164   if(fD0toKpi && !aodD0toKpiTClArr) {
165     printf("ERROR: no aodD0toKpiTClArr");
166     return;
167   }
168   if(fJPSItoEle && !aodJPSItoEleTClArr) {
169     printf("ERROR: no aodJPSItoEleTClArr");
170     return;
171   }
172   if(f3Prong && !aodCharm3ProngTClArr) {
173     printf("ERROR: no aodCharm3ProngTClArr");
174     return;
175   }
176   if(f4Prong && !aodCharm4ProngTClArr) {
177     printf("ERROR: no aodCharm4ProngTClArr");
178     return;
179   }
180   if(fDstar && !aodDstarTClArr) {
181     printf("ERROR: no aodDstarTClArr");
182     return;
183   }
184
185   // delete candidates from previous event and create references
186   Int_t iVerticesHF=0,iD0toKpi=0,iJPSItoEle=0,i3Prong=0,i4Prong=0,iDstar=0;
187   aodVerticesHFTClArr->Delete();
188   iVerticesHF = aodVerticesHFTClArr->GetEntriesFast();
189   TClonesArray &verticesHFRef = *aodVerticesHFTClArr;
190   if(fD0toKpi)   {
191     aodD0toKpiTClArr->Delete();
192     iD0toKpi = aodD0toKpiTClArr->GetEntriesFast();
193   }
194   if(fJPSItoEle) {
195     aodJPSItoEleTClArr->Delete();
196     iJPSItoEle = aodJPSItoEleTClArr->GetEntriesFast();
197   }
198   if(f3Prong) {   
199     aodCharm3ProngTClArr->Delete();
200     i3Prong = aodCharm3ProngTClArr->GetEntriesFast();
201   }
202   if(f4Prong) {
203     aodCharm4ProngTClArr->Delete();
204     i4Prong = aodCharm4ProngTClArr->GetEntriesFast();
205   }
206   if(fDstar) {
207     aodDstarTClArr->Delete();
208     iDstar = aodDstarTClArr->GetEntriesFast();
209   }
210   TClonesArray &aodD0toKpiRef      = *aodD0toKpiTClArr;
211   TClonesArray &aodJPSItoEleRef    = *aodJPSItoEleTClArr;
212   TClonesArray &aodCharm3ProngRef  = *aodCharm3ProngTClArr;
213   TClonesArray &aodCharm4ProngRef  = *aodCharm4ProngTClArr;
214   TClonesArray &aodDstarRef        = *aodDstarTClArr;
215   
216
217   AliAODRecoDecayHF2Prong *io2Prong  = 0;
218   AliAODRecoDecayHF3Prong *io3Prong  = 0;
219   AliAODRecoDecayHF4Prong *io4Prong  = 0;
220   AliAODRecoCascadeHF     *ioCascade = 0;
221
222   Int_t    iTrkP1,iTrkP2,iTrkN1,iTrkN2,iTrkSoftPi,trkEntries;
223   Double_t xdummy,ydummy,dcap1n1,dcap1n2,dcap2n1,dcap1p2,dcan1n2,dcaCasc;
224   Bool_t   okD0=kFALSE,okJPSI=kFALSE,ok3Prong=kFALSE,ok4Prong=kFALSE;
225   Bool_t   okDstar=kFALSE,okD0fromDstar=kFALSE;
226   AliESDtrack *postrack1 = 0;
227   AliESDtrack *postrack2 = 0;
228   AliESDtrack *negtrack1 = 0;
229   AliESDtrack *negtrack2 = 0;
230   AliESDtrack *trackPi   = 0;
231   Double_t dcaMax = fD0toKpiCuts[1];
232   if(dcaMax < fBtoJPSICuts[1]) dcaMax=fBtoJPSICuts[1];
233   if(dcaMax < fDplusCuts[11])  dcaMax=fDplusCuts[11];
234   AliDebug(2,Form(" dca cut set to %f cm",dcaMax));
235
236
237   // get Bz
238   fBzkG = (Double_t)event->GetMagneticField(); 
239
240   trkEntries = (Int_t)event->GetNumberOfTracks();
241   AliDebug(1,Form(" Number of tracks: %d",trkEntries));
242
243   if(trkEntries<2) {
244     AliDebug(1,Form(" Not enough tracks: %d",trkEntries));
245     return;
246   }
247     
248   const AliVVertex *primary = event->GetPrimaryVertex();
249   if(!primary) {
250     AliDebug(1," No primary vertex from tracks");
251     return;
252   }
253   TString primTitle = primary->GetTitle();
254   if(!primTitle.Contains("VertexerTracks")) {
255     AliDebug(1," No primary vertex from tracks");
256     return;
257   }
258
259   // call function that applies sigle-track selection,
260   // for displaced tracks and soft pions (both charges) for D*,
261   // and retrieves primary vertex
262   TObjArray seleTrksArray(trkEntries);
263   UChar_t  *seleFlags = new UChar_t[trkEntries]; // bit 0: displaced, bit 1: softpi
264   Int_t     nSeleTrks=0;
265   SelectTracksAndCopyVertex(event,seleTrksArray,nSeleTrks,seleFlags);
266     
267   AliDebug(1,Form(" Selected tracks: %d",nSeleTrks));
268     
269   TObjArray *twoTrackArray1    = new TObjArray(2);
270   TObjArray *twoTrackArray2    = new TObjArray(2);
271   TObjArray *twoTrackArrayCasc = new TObjArray(2);
272   TObjArray *threeTrackArray   = new TObjArray(3);
273   TObjArray *fourTrackArray    = new TObjArray(4);
274   
275   Double_t dispersion;
276
277   AliAODRecoDecayHF   *rd = 0;
278   AliAODRecoCascadeHF *rc = 0;
279
280   // LOOP ON  POSITIVE  TRACKS
281   for(iTrkP1=0; iTrkP1<nSeleTrks; iTrkP1++) {
282     if(iTrkP1%1==0) AliDebug(1,Form("  1st loop on pos: track number %d of %d",iTrkP1,nSeleTrks));  
283     // get track from tracks array
284     postrack1 = (AliESDtrack*)seleTrksArray.UncheckedAt(iTrkP1);
285     if(postrack1->Charge()<0 || !TESTBIT(seleFlags[iTrkP1],kBitDispl)) continue;
286     // LOOP ON  NEGATIVE  TRACKS
287     for(iTrkN1=0; iTrkN1<nSeleTrks; iTrkN1++) {
288       if(iTrkN1==iTrkP1) continue;
289       if(iTrkN1%1==0) AliDebug(1,Form("    1st loop on neg: track number %d of %d",iTrkN1,nSeleTrks));  
290       // get track from tracks array
291       negtrack1 = (AliESDtrack*)seleTrksArray.UncheckedAt(iTrkN1);
292       if(negtrack1->Charge()>0 || !TESTBIT(seleFlags[iTrkN1],kBitDispl)) continue;
293       // back to primary vertex
294       postrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
295       negtrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
296       // DCA between the two tracks
297       dcap1n1 = postrack1->GetDCA(negtrack1,fBzkG,xdummy,ydummy);
298       if(dcap1n1>dcaMax) { negtrack1=0; continue; }
299       // Vertexing
300       twoTrackArray1->AddAt(postrack1,0);
301       twoTrackArray1->AddAt(negtrack1,1);
302       AliAODVertex *vertexp1n1 = ReconstructSecondaryVertex(twoTrackArray1,dispersion);
303       if(!vertexp1n1) { 
304         twoTrackArray1->Clear();
305         negtrack1=0; 
306         continue; 
307       }
308       if(fD0toKpi || fJPSItoEle || fDstar) { 
309         io2Prong = Make2Prong(twoTrackArray1,event,vertexp1n1,dcap1n1,okD0,okJPSI,okD0fromDstar);
310         if((fD0toKpi && okD0) || (fJPSItoEle && okJPSI)) {
311           // add the vertex and the decay to the AOD
312           AliAODVertex *v2Prong = new(verticesHFRef[iVerticesHF++])AliAODVertex(*vertexp1n1);
313           if(fInputAOD) AddDaughterRefs(v2Prong,event,twoTrackArray1);
314           if(okD0) {  
315             rd = new(aodD0toKpiRef[iD0toKpi++])AliAODRecoDecayHF2Prong(*io2Prong);
316             rd->SetSecondaryVtx(v2Prong);
317             v2Prong->SetParent(rd);
318           }
319           if(okJPSI) {
320             rd = new(aodJPSItoEleRef[iJPSItoEle++])AliAODRecoDecayHF2Prong(*io2Prong);
321             rd->SetSecondaryVtx(v2Prong);
322             if(!okD0) v2Prong->SetParent(rd); // it cannot have two mothers ...
323           }
324         }
325         if(fDstar && okD0fromDstar) {
326           if(fInputAOD) { // write references in io2Prong
327             AddDaughterRefs(vertexp1n1,event,twoTrackArray1);
328             io2Prong->SetSecondaryVtx(vertexp1n1);
329           }
330           // create a track from the D0
331           //AliESDtrack *trackD0 = new AliESDtrack(io2Prong); // to be uncommented after commit of AliAODRecoDecay : public AliVTrack
332           AliESDtrack *trackD0 = new AliESDtrack(); // temporary, just to allow compilation
333           // LOOP ON TRACKS THAT PASSED THE SOFT PION CUTS
334           for(iTrkSoftPi=0; iTrkSoftPi<nSeleTrks; iTrkSoftPi++) {
335             if(iTrkSoftPi==iTrkP1 || iTrkSoftPi==iTrkN1) continue;
336             if(!TESTBIT(seleFlags[iTrkSoftPi],kBitSoftPi)) continue;
337             if(iTrkSoftPi%1==0) AliDebug(1,Form("    1st loop on pi_s: track number %d of %d",iTrkSoftPi,nSeleTrks));  
338             // get track from tracks array
339             trackPi = (AliESDtrack*)seleTrksArray.UncheckedAt(iTrkSoftPi);
340
341             trackPi->PropagateToDCA(fV1,fBzkG,kVeryBig);
342             trackD0->PropagateToDCA(fV1,fBzkG,kVeryBig);
343             // DCA between the two tracks
344             dcaCasc = trackPi->GetDCA(trackD0,fBzkG,xdummy,ydummy);
345             // Vertexing
346             twoTrackArrayCasc->AddAt(trackPi,0);
347             twoTrackArrayCasc->AddAt(trackD0,1);
348             AliAODVertex *vertexCasc = ReconstructSecondaryVertex(twoTrackArrayCasc,dispersion);
349             if(!vertexCasc) { 
350               twoTrackArrayCasc->Clear();
351               trackPi=0; 
352               continue; 
353             }
354             ioCascade = MakeCascade(twoTrackArrayCasc,event,vertexCasc,io2Prong,dcaCasc,okDstar);
355             if(okDstar) {
356               if(!okD0) {
357                 // add the D0 to the AOD (if not already done)
358                 AliAODVertex *v2Prong = new(verticesHFRef[iVerticesHF++])AliAODVertex(*vertexp1n1);
359                 if(fInputAOD) AddDaughterRefs(v2Prong,event,twoTrackArray1);
360                 rd = new(aodD0toKpiRef[iD0toKpi++])AliAODRecoDecayHF2Prong(*io2Prong);
361                 rd->SetSecondaryVtx(v2Prong);
362                 v2Prong->SetParent(rd);
363                 okD0=kTRUE; // this is done to add it only once
364               }
365               // add the vertex and the cascade to the AOD
366               AliAODVertex *vCasc = new(verticesHFRef[iVerticesHF++])AliAODVertex(*vertexCasc); 
367               if(fInputAOD) AddDaughterRefs(vCasc,event,twoTrackArrayCasc); // add the pion
368               //vCasc->AddDaughter(rd); // doesn't work, because the D0 and the bachelor have different ROOT process ID
369               rc = new(aodDstarRef[iDstar++])AliAODRecoCascadeHF(*ioCascade);
370               rc->Set2Prong(rd); // add the D0
371               rc->SetSecondaryVtx(vCasc);
372               vCasc->SetParent(rc);
373             }
374             twoTrackArrayCasc->Clear();
375             trackPi=0; 
376             ioCascade=NULL;
377             delete vertexCasc;
378           } // end loop on soft pi tracks
379           if(trackD0) {delete trackD0; trackD0=NULL;}
380         }
381         io2Prong=NULL;
382       }      
383       twoTrackArray1->Clear(); 
384       if(!f3Prong && !f4Prong)  { 
385         negtrack1=0; 
386         delete vertexp1n1; 
387         continue; 
388       }
389
390         
391       // 2nd LOOP  ON  POSITIVE  TRACKS 
392       for(iTrkP2=iTrkP1+1; iTrkP2<nSeleTrks; iTrkP2++) {
393         if(iTrkP2==iTrkP1 || iTrkP2==iTrkN1) continue;
394         if(iTrkP2%1==0) AliDebug(1,Form("    2nd loop on pos: track number %d of %d",iTrkP2,nSeleTrks));  
395         // get track from tracks array
396         postrack2 = (AliESDtrack*)seleTrksArray.UncheckedAt(iTrkP2);
397         if(postrack2->Charge()<0 || !TESTBIT(seleFlags[iTrkP2],kBitDispl)) continue;
398         // back to primary vertex
399         postrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
400         postrack2->PropagateToDCA(fV1,fBzkG,kVeryBig);
401         negtrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
402         //printf("********** %d %d %d\n",postrack1->GetID(),postrack2->GetID(),negtrack1->GetID());
403         dcap2n1 = postrack2->GetDCA(negtrack1,fBzkG,xdummy,ydummy);
404         if(dcap2n1>dcaMax) { postrack2=0; continue; }
405         dcap1p2 = postrack2->GetDCA(postrack1,fBzkG,xdummy,ydummy);
406         if(dcap1p2>dcaMax) { postrack2=0; continue; }
407         
408         // Vertexing
409         twoTrackArray2->AddAt(postrack2,0);
410         twoTrackArray2->AddAt(negtrack1,1);
411         AliAODVertex *vertexp2n1 = ReconstructSecondaryVertex(twoTrackArray2,dispersion);
412         if(!vertexp2n1) { 
413           twoTrackArray2->Clear();
414           postrack2=0; 
415           continue; 
416         }
417         if(f3Prong) { 
418           threeTrackArray->AddAt(postrack1,0);
419           threeTrackArray->AddAt(negtrack1,1);
420           threeTrackArray->AddAt(postrack2,2);
421           AliAODVertex* secVert3PrAOD = ReconstructSecondaryVertex(threeTrackArray,dispersion);
422           io3Prong = Make3Prong(threeTrackArray,event,secVert3PrAOD,dispersion,vertexp1n1,vertexp2n1,dcap1n1,dcap2n1,dcap1p2,ok3Prong);
423           if(ok3Prong) {
424             AliAODVertex *v3Prong = new(verticesHFRef[iVerticesHF++])AliAODVertex(*secVert3PrAOD);
425             if(fInputAOD) AddDaughterRefs(v3Prong,event,threeTrackArray);
426             rd = new(aodCharm3ProngRef[i3Prong++])AliAODRecoDecayHF3Prong(*io3Prong);
427             rd->SetSecondaryVtx(v3Prong);
428             v3Prong->SetParent(rd);
429           }
430           if(io3Prong) io3Prong=NULL; 
431         }
432         if(f4Prong) {
433           // 3rd LOOP  ON  NEGATIVE  TRACKS (for 4 prong) 
434           for(iTrkN2=iTrkN1+1; iTrkN2<nSeleTrks; iTrkN2++) {
435             if(iTrkN2==iTrkP1 || iTrkN2==iTrkP2 || iTrkN2==iTrkN1) continue;
436             if(iTrkN2%1==0) AliDebug(1,Form("    3rd loop on neg: track number %d of %d",iTrkN2,nSeleTrks));  
437             // get track from tracks array
438             negtrack2 = (AliESDtrack*)seleTrksArray.UncheckedAt(iTrkN2);
439             if(negtrack2->Charge()>0 || !TESTBIT(seleFlags[iTrkN2],kBitDispl)) continue;
440             // back to primary vertex
441             postrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
442             postrack2->PropagateToDCA(fV1,fBzkG,kVeryBig);
443             negtrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
444             negtrack2->PropagateToDCA(fV1,fBzkG,kVeryBig);
445             dcap1n2 = postrack1->GetDCA(negtrack2,fBzkG,xdummy,ydummy);
446             if(dcap1n2>dcaMax) { negtrack2=0; continue; }
447             // Vertexing
448             fourTrackArray->AddAt(postrack1,0);
449             fourTrackArray->AddAt(negtrack1,1);
450             fourTrackArray->AddAt(postrack2,2);
451             fourTrackArray->AddAt(negtrack2,3);
452             AliAODVertex* secVert4PrAOD = ReconstructSecondaryVertex(fourTrackArray,dispersion);
453             io4Prong = Make4Prong(fourTrackArray,event,secVert4PrAOD,vertexp1n1,vertexp2n1,dcap1n1,dcap1n2,dcap2n1,ok4Prong);
454             if(ok4Prong) {
455               AliAODVertex *v4Prong = new(verticesHFRef[iVerticesHF++])AliAODVertex(*secVert4PrAOD);
456               if(fInputAOD) AddDaughterRefs(v4Prong,event,fourTrackArray);
457               rd = new(aodCharm4ProngRef[i4Prong++])AliAODRecoDecayHF4Prong(*io4Prong);
458               rd->SetSecondaryVtx(v4Prong);
459               v4Prong->SetParent(rd);
460             }
461             if(io4Prong) io4Prong=NULL; 
462             fourTrackArray->Clear();
463             negtrack2 = 0;
464           } // end loop on negative tracks
465         }
466         postrack2 = 0;
467         delete vertexp2n1;
468       } // end 2nd loop on positive tracks
469       twoTrackArray2->Clear();
470       
471       // 2nd LOOP  ON  NEGATIVE  TRACKS 
472       for(iTrkN2=iTrkN1+1; iTrkN2<nSeleTrks; iTrkN2++) {
473         if(iTrkN2==iTrkP1 || iTrkN2==iTrkP2 || iTrkN2==iTrkN1) continue;
474         if(iTrkN2%1==0) AliDebug(1,Form("    2nd loop on neg: track number %d of %d",iTrkN2,nSeleTrks));  
475         // get track from tracks array
476         negtrack2 = (AliESDtrack*)seleTrksArray.UncheckedAt(iTrkN2);
477         if(negtrack2->Charge()>0 || !TESTBIT(seleFlags[iTrkN2],kBitDispl)) continue;
478         // back to primary vertex
479         postrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
480         negtrack1->PropagateToDCA(fV1,fBzkG,kVeryBig);
481         negtrack2->PropagateToDCA(fV1,fBzkG,kVeryBig);
482         //printf("********** %d %d %d\n",postrack1->GetID(),negtrack1->GetID(),negtrack2->GetID());
483         dcap1n2 = postrack1->GetDCA(negtrack2,fBzkG,xdummy,ydummy);
484         if(dcap1n2>dcaMax) { negtrack2=0; continue; }
485         dcan1n2 = negtrack1->GetDCA(negtrack2,fBzkG,xdummy,ydummy);
486         if(dcan1n2>dcaMax) { negtrack2=0; continue; }
487         
488         // Vertexing
489         twoTrackArray2->AddAt(postrack1,0);
490         twoTrackArray2->AddAt(negtrack2,1);
491
492         AliAODVertex *vertexp1n2 = ReconstructSecondaryVertex(twoTrackArray2,dispersion);
493         if(!vertexp1n2) { 
494           twoTrackArray2->Clear();
495           negtrack2=0; 
496           continue; 
497         }
498         if(f3Prong) { 
499           threeTrackArray->AddAt(negtrack1,0);
500           threeTrackArray->AddAt(postrack1,1);
501           threeTrackArray->AddAt(negtrack2,2);
502           AliAODVertex* secVert3PrAOD = ReconstructSecondaryVertex(threeTrackArray,dispersion);
503           io3Prong = Make3Prong(threeTrackArray,event,secVert3PrAOD,dispersion,vertexp1n1,vertexp1n2,dcap1n1,dcap1n2,dcan1n2,ok3Prong);
504           if(ok3Prong) {
505             AliAODVertex *v3Prong = new(verticesHFRef[iVerticesHF++])AliAODVertex(*secVert3PrAOD);
506             if(fInputAOD) AddDaughterRefs(v3Prong,event,threeTrackArray);
507             rd = new(aodCharm3ProngRef[i3Prong++])AliAODRecoDecayHF3Prong(*io3Prong);
508             rd->SetSecondaryVtx(v3Prong);
509             v3Prong->SetParent(rd);
510           }
511           if(io3Prong) io3Prong=NULL; 
512         }
513         negtrack2 = 0;
514         delete vertexp1n2;
515       } // end 2nd loop on negative tracks
516       twoTrackArray2->Clear();
517       
518       negtrack1 = 0;
519       delete vertexp1n1; 
520     } // end 1st loop on negative tracks
521     
522     postrack1 = 0;
523   }  // end 1st loop on positive tracks
524
525
526   AliDebug(1,Form(" Total HF vertices in event = %d;",
527                   (Int_t)aodVerticesHFTClArr->GetEntriesFast()));
528   if(fD0toKpi) {
529     AliDebug(1,Form(" D0->Kpi in event = %d;",
530                     (Int_t)aodD0toKpiTClArr->GetEntriesFast()));
531   }
532   if(fJPSItoEle) {
533     AliDebug(1,Form(" JPSI->ee in event = %d;",
534                     (Int_t)aodJPSItoEleTClArr->GetEntriesFast()));
535   }
536   if(f3Prong) {
537     AliDebug(1,Form(" Charm->3Prong in event = %d;",
538                     (Int_t)aodCharm3ProngTClArr->GetEntriesFast()));
539   }
540   if(f4Prong) {
541     AliDebug(1,Form(" Charm->4Prong in event = %d;\n",
542                     (Int_t)aodCharm4ProngTClArr->GetEntriesFast()));
543   }
544   if(fDstar) {
545     AliDebug(1,Form(" D*->D0pi in event = %d;\n",
546                     (Int_t)aodDstarTClArr->GetEntriesFast()));
547   }
548     
549
550   twoTrackArray1->Delete();  delete twoTrackArray1;
551   twoTrackArray2->Delete();  delete twoTrackArray2;
552   twoTrackArrayCasc->Delete();  delete twoTrackArrayCasc;
553   threeTrackArray->Clear(); 
554   threeTrackArray->Delete(); delete threeTrackArray;
555   fourTrackArray->Delete();  delete fourTrackArray;
556   delete [] seleFlags; seleFlags=NULL;
557
558   if(fInputAOD) {
559     seleTrksArray.Delete();
560   }
561
562   return;
563 }
564 //----------------------------------------------------------------------------
565 void AliAnalysisVertexingHF::AddDaughterRefs(AliAODVertex *v,AliVEvent *event,
566                                              TObjArray *trkArray) const
567 {
568   // Add the AOD tracks as daughters of the vertex (TRef)
569
570   Int_t nTrks = trkArray->GetEntriesFast();
571
572   AliESDtrack *esdTrack = 0;
573   AliAODTrack *aodTrack = 0;
574   Int_t id;
575
576   for(Int_t i=0; i<nTrks; i++) {
577     esdTrack = (AliESDtrack*)trkArray->UncheckedAt(i);
578     id = (Int_t)esdTrack->GetID();
579     //printf("---> %d\n",id);
580     if(id<0) continue; // this track is a D0
581     aodTrack = (AliAODTrack*)event->GetTrack(fAODMap[id]);
582     v->AddDaughter(aodTrack);
583   }
584
585   return;
586 }       
587 //----------------------------------------------------------------------------
588 AliAODRecoCascadeHF* AliAnalysisVertexingHF::MakeCascade(
589                                    TObjArray *twoTrackArray,AliVEvent *event,
590                                    AliAODVertex *secVert,
591                                    AliAODRecoDecayHF2Prong *rd2Prong,
592                                    Double_t dca,
593                                    Bool_t &okDstar) const
594 {
595   // Make the cascade as a 2Prong decay and check if it passes Dstar
596   // reconstruction cuts
597
598   okDstar = kFALSE;
599
600   Bool_t dummy1,dummy2,dummy3;
601
602   // We use Make2Prong to construct the AliAODRecoCascadeHF
603   // (which inherits from AliAODRecoDecayHF2Prong) 
604   AliAODRecoCascadeHF *theCascade = 
605     (AliAODRecoCascadeHF*)Make2Prong(twoTrackArray,event,secVert,dca,
606                                      dummy1,dummy2,dummy3);
607   if(!theCascade) return 0x0;
608
609   // charge
610   AliESDtrack *trackPi = (AliESDtrack*)twoTrackArray->UncheckedAt(0);
611   theCascade->SetCharge(trackPi->Charge());
612
613   //--- selection cuts
614   //
615   AliAODRecoCascadeHF *tmpCascade = new AliAODRecoCascadeHF(*theCascade);
616   tmpCascade->Set2Prong(rd2Prong); // this is needed for the selection
617   AliAODVertex *primVertexAOD=0;
618   if(!fRecoPrimVtxSkippingTrks && !fRmTrksFromPrimVtx) {
619     // take event primary vertex
620     primVertexAOD = PrimaryVertex(); 
621     tmpCascade->SetOwnPrimaryVtx(primVertexAOD);
622     rd2Prong->SetOwnPrimaryVtx(primVertexAOD);
623   }
624   // select D*->D0pi
625   if(fDstar) {
626     Bool_t testD0=kTRUE;
627     okDstar = tmpCascade->SelectDstar(fDstarCuts,fD0fromDstarCuts,testD0);
628   }
629   tmpCascade->UnsetOwnPrimaryVtx(); 
630   tmpCascade->Unset2Prong(); 
631   delete tmpCascade; tmpCascade=NULL;
632   if(!fRecoPrimVtxSkippingTrks && !fRmTrksFromPrimVtx) {
633     rd2Prong->UnsetOwnPrimaryVtx();
634     delete primVertexAOD; primVertexAOD=NULL;
635   }
636   //---
637   
638   return theCascade;
639 }
640 //-----------------------------------------------------------------------------
641 AliAODRecoDecayHF2Prong *AliAnalysisVertexingHF::Make2Prong(
642                                    TObjArray *twoTrackArray,AliVEvent *event,
643                                    AliAODVertex *secVert,Double_t dca,
644                                    Bool_t &okD0,Bool_t &okJPSI,
645                                    Bool_t &okD0fromDstar) const
646 {
647   // Make 2Prong candidates and check if they pass D0toKpi or BtoJPSI
648   // reconstruction cuts
649   // G.E.Bruno (J/psi), A.Dainese (D0->Kpi)
650
651   okD0=kFALSE; okJPSI=kFALSE; okD0fromDstar=kFALSE;
652
653   Double_t px[2],py[2],pz[2],d0[2],d0err[2];
654
655   AliESDtrack *postrack = (AliESDtrack*)twoTrackArray->UncheckedAt(0);
656   AliESDtrack *negtrack = (AliESDtrack*)twoTrackArray->UncheckedAt(1);
657
658   // propagate tracks to secondary vertex, to compute inv. mass
659   postrack->PropagateToDCA(secVert,fBzkG,kVeryBig);
660   negtrack->PropagateToDCA(secVert,fBzkG,kVeryBig);
661
662   Double_t momentum[3];
663   postrack->GetPxPyPz(momentum);
664   px[0] = momentum[0]; py[0] = momentum[1]; pz[0] = momentum[2]; 
665   negtrack->GetPxPyPz(momentum);
666   px[1] = momentum[0]; py[1] = momentum[1]; pz[1] = momentum[2]; 
667
668
669   // invariant mass cut (try to improve coding here..)
670   Bool_t okMassCut=kFALSE;
671   if(!okMassCut && fD0toKpi)   if(SelectInvMass(0,2,px,py,pz)) okMassCut=kTRUE;
672   if(!okMassCut && fJPSItoEle) if(SelectInvMass(1,2,px,py,pz)) okMassCut=kTRUE;
673   if(!okMassCut && fDstar)     if(SelectInvMass(3,2,px,py,pz)) okMassCut=kTRUE;
674   if(!okMassCut) {
675     AliDebug(2," candidate didn't pass mass cut");
676     return 0x0;    
677   }
678
679   // primary vertex to be used by this candidate
680   AliAODVertex *primVertexAOD  = PrimaryVertex(twoTrackArray,event);
681   if(!primVertexAOD) return 0x0;
682
683   Double_t d0z0[2],covd0z0[3];
684   postrack->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
685   d0[0] = d0z0[0];
686   d0err[0] = TMath::Sqrt(covd0z0[0]);
687   negtrack->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
688   d0[1] = d0z0[0];
689   d0err[1] = TMath::Sqrt(covd0z0[0]);
690
691   // create the object AliAODRecoDecayHF2Prong
692   AliAODRecoDecayHF2Prong *the2Prong = new AliAODRecoDecayHF2Prong(secVert,px,py,pz,d0,d0err,dca);
693   the2Prong->SetOwnPrimaryVtx(primVertexAOD);
694   UShort_t id[2]={(UShort_t)postrack->GetID(),(UShort_t)negtrack->GetID()};
695   the2Prong->SetProngIDs(2,id);
696
697
698   // select D0->Kpi
699   Int_t checkD0,checkD0bar;
700   if(fD0toKpi)   okD0          = the2Prong->SelectD0(fD0toKpiCuts,checkD0,checkD0bar);
701   //if(fDebug && fD0toKpi) printf("   %d\n",(Int_t)okD0);
702   // select J/psi from B
703   Int_t checkJPSI;
704   if(fJPSItoEle) okJPSI        = the2Prong->SelectBtoJPSI(fBtoJPSICuts,checkJPSI);
705   //if(fDebug && fJPSItoEle) printf("   %d\n",(Int_t)okJPSI);
706   // select D0->Kpi from Dstar
707   if(fDstar)     okD0fromDstar = the2Prong->SelectD0(fD0fromDstarCuts,checkD0,checkD0bar);
708   //if(fDebug && fDstar) printf("   %d\n",(Int_t)okD0fromDstar);
709
710   // remove the primary vertex (was used only for selection)
711   if(!fRecoPrimVtxSkippingTrks && !fRmTrksFromPrimVtx) the2Prong->UnsetOwnPrimaryVtx();
712
713   // get PID info from ESD
714   Double_t esdpid0[5];
715   postrack->GetESDpid(esdpid0);
716   Double_t esdpid1[5];
717   negtrack->GetESDpid(esdpid1);
718   Double_t esdpid[10];
719   for(Int_t i=0;i<5;i++) {
720     esdpid[i]   = esdpid0[i];
721     esdpid[5+i] = esdpid1[i];
722   }
723   the2Prong->SetPID(2,esdpid);
724
725   return the2Prong;  
726 }
727 //----------------------------------------------------------------------------
728 AliAODRecoDecayHF3Prong* AliAnalysisVertexingHF::Make3Prong(
729                              TObjArray *threeTrackArray,AliVEvent *event,
730                              AliAODVertex *secVert,Double_t dispersion,
731                              AliAODVertex *vertexp1n1,AliAODVertex *vertexp2n1,
732                              Double_t dcap1n1,Double_t dcap2n1,Double_t dcap1p2,
733                              Bool_t &ok3Prong) const
734 {
735   // Make 3Prong candidates and check if they pass Dplus or Ds or Lambdac
736   // reconstruction cuts 
737   // E.Bruna, F.Prino
738
739
740   ok3Prong=kFALSE;
741   if(!secVert) return 0x0; 
742
743   Double_t px[3],py[3],pz[3],d0[3],d0err[3];
744   Double_t momentum[3];
745
746
747   AliESDtrack *postrack1 = (AliESDtrack*)threeTrackArray->UncheckedAt(0);
748   AliESDtrack *negtrack  = (AliESDtrack*)threeTrackArray->UncheckedAt(1);
749   AliESDtrack *postrack2 = (AliESDtrack*)threeTrackArray->UncheckedAt(2);
750
751   postrack1->PropagateToDCA(secVert,fBzkG,kVeryBig);
752   negtrack->PropagateToDCA(secVert,fBzkG,kVeryBig);
753   postrack2->PropagateToDCA(secVert,fBzkG,kVeryBig);
754   postrack1->GetPxPyPz(momentum);
755   px[0] = momentum[0]; py[0] = momentum[1]; pz[0] = momentum[2]; 
756   negtrack->GetPxPyPz(momentum);
757   px[1] = momentum[0]; py[1] = momentum[1]; pz[1] = momentum[2]; 
758   postrack2->GetPxPyPz(momentum);
759   px[2] = momentum[0]; py[2] = momentum[1]; pz[2] = momentum[2]; 
760
761   // invariant mass cut for D+, Ds, Lc
762   Bool_t okMassCut=kFALSE;
763   if(!okMassCut && f3Prong) if(SelectInvMass(2,3,px,py,pz)) okMassCut=kTRUE;
764   if(!okMassCut) {
765     AliDebug(2," candidate didn't pass mass cut");
766     return 0x0;    
767   }
768
769   // primary vertex to be used by this candidate
770   AliAODVertex *primVertexAOD  = PrimaryVertex(threeTrackArray,event);
771   if(!primVertexAOD) return 0x0;
772
773   Double_t d0z0[2],covd0z0[3];
774   postrack1->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
775   d0[0]=d0z0[0];
776   d0err[0] = TMath::Sqrt(covd0z0[0]);
777   negtrack->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
778   d0[1]=d0z0[0];
779   d0err[1] = TMath::Sqrt(covd0z0[0]);
780   postrack2->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
781   d0[2]=d0z0[0];
782   d0err[2] = TMath::Sqrt(covd0z0[0]);
783
784
785   // create the object AliAODRecoDecayHF3Prong
786   Double_t pos[3]; primVertexAOD->GetXYZ(pos);
787   Double_t dca[3]={dcap1n1,dcap2n1,dcap1p2};
788   Double_t dist12=TMath::Sqrt((vertexp1n1->GetX()-pos[0])*(vertexp1n1->GetX()-pos[0])+(vertexp1n1->GetY()-pos[1])*(vertexp1n1->GetY()-pos[1])+(vertexp1n1->GetZ()-pos[2])*(vertexp1n1->GetZ()-pos[2]));
789   Double_t dist23=TMath::Sqrt((vertexp2n1->GetX()-pos[0])*(vertexp2n1->GetX()-pos[0])+(vertexp2n1->GetY()-pos[1])*(vertexp2n1->GetY()-pos[1])+(vertexp2n1->GetZ()-pos[2])*(vertexp2n1->GetZ()-pos[2]));
790   Short_t charge=(Short_t)(postrack1->Charge()*postrack2->Charge()*negtrack->Charge());
791   AliAODRecoDecayHF3Prong *the3Prong = new AliAODRecoDecayHF3Prong(secVert,px,py,pz,d0,d0err,dca,dispersion,dist12,dist23,charge);
792   the3Prong->SetOwnPrimaryVtx(primVertexAOD);
793   UShort_t id[3]={(UShort_t)postrack1->GetID(),(UShort_t)negtrack->GetID(),(UShort_t)postrack2->GetID()};
794   the3Prong->SetProngIDs(3,id);
795
796
797   // select D+->Kpipi, Ds->KKpi, Lc->pKpi
798   if(f3Prong) {
799     ok3Prong = kFALSE;
800     Int_t ok1,ok2;
801     if(the3Prong->SelectDplus(fDplusCuts))   ok3Prong = kTRUE;
802     if(the3Prong->SelectDs(fDsCuts,ok1,ok2)) ok3Prong = kTRUE;
803     if(the3Prong->SelectLc(fLcCuts,ok1,ok2)) ok3Prong = kTRUE;
804   }
805   //if(fDebug) printf("ok3Prong: %d\n",(Int_t)ok3Prong);
806
807   if(!fRecoPrimVtxSkippingTrks && !fRmTrksFromPrimVtx) the3Prong->UnsetOwnPrimaryVtx();
808
809   // get PID info from ESD
810   Double_t esdpid0[5];
811   postrack1->GetESDpid(esdpid0);
812   Double_t esdpid1[5];
813   negtrack->GetESDpid(esdpid1);
814   Double_t esdpid2[5];
815   postrack2->GetESDpid(esdpid2);
816   
817   Double_t esdpid[15];
818   for(Int_t i=0;i<5;i++) {
819     esdpid[i]    = esdpid0[i];
820     esdpid[5+i]  = esdpid1[i];
821     esdpid[10+i] = esdpid2[i];
822   }
823   the3Prong->SetPID(3,esdpid);
824
825   return the3Prong;
826 }
827 //----------------------------------------------------------------------------
828 AliAODRecoDecayHF4Prong* AliAnalysisVertexingHF::Make4Prong(
829                              TObjArray *fourTrackArray,AliVEvent *event,
830                              AliAODVertex *secVert,
831                              AliAODVertex *vertexp1n1,AliAODVertex *vertexp2n1,
832                              Double_t dcap1n1,Double_t dcap1n2,Double_t dcap2n1,
833                              Bool_t &ok4Prong) const
834 {
835   // Make 4Prong candidates and check if they pass D0toKpipipi
836   // reconstruction cuts
837   // G.E.Bruno, R.Romita
838
839   ok4Prong=kFALSE;
840   if(!secVert) return 0x0; 
841
842   Double_t px[4],py[4],pz[4],d0[4],d0err[4];//d0z[3];
843
844   px[0]=dcap1n1*dcap1n2*dcap2n1; // TO BE CHANGED (done just to removed compilation warning about dca... not used)
845
846   AliESDtrack *postrack1 = (AliESDtrack*)fourTrackArray->UncheckedAt(0);
847   AliESDtrack *negtrack1 = (AliESDtrack*)fourTrackArray->UncheckedAt(1);
848   AliESDtrack *postrack2 = (AliESDtrack*)fourTrackArray->UncheckedAt(2);
849   AliESDtrack *negtrack2 = (AliESDtrack*)fourTrackArray->UncheckedAt(3);
850
851   postrack1->PropagateToDCA(secVert,fBzkG,kVeryBig);
852   negtrack1->PropagateToDCA(secVert,fBzkG,kVeryBig);
853   postrack2->PropagateToDCA(secVert,fBzkG,kVeryBig);
854   negtrack2->PropagateToDCA(secVert,fBzkG,kVeryBig);
855
856   Double_t momentum[3];
857   postrack1->GetPxPyPz(momentum);
858   px[0] = momentum[0]; py[0] = momentum[1]; pz[0] = momentum[2];
859   negtrack1->GetPxPyPz(momentum);
860   px[1] = momentum[0]; py[1] = momentum[1]; pz[1] = momentum[2];
861   postrack2->GetPxPyPz(momentum);
862   px[2] = momentum[0]; py[2] = momentum[1]; pz[2] = momentum[2];
863   negtrack2->GetPxPyPz(momentum);
864   px[3] = momentum[0]; py[3] = momentum[1]; pz[3] = momentum[2];
865
866   // invariant mass cut for rho or D0 (try to improve coding here..)
867   //Bool_t okMassCut=kFALSE;
868   //if(!okMassCut) if(SelectInvMass(2,3,px,py,pz)) okMassCut=kTRUE;
869   //if(!okMassCut) {
870   //  if(fDebug) printf(" candidate didn't pass mass cut\n");
871   //  return 0x0;
872   //}
873
874   // primary vertex to be used by this candidate
875   AliAODVertex *primVertexAOD  = PrimaryVertex(fourTrackArray,event);
876   if(!primVertexAOD) return 0x0;
877
878   /*
879     Double_t d0z0[2],covd0z0[3];
880     postrack1->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
881     negtrack1->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
882     postrack2->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
883     negtrack2->PropagateToDCA(primVertexAOD,fBzkG,kVeryBig,d0z0,covd0z0);
884   */
885
886   // create the object AliAODRecoDecayHF4Prong
887   Double_t pos[3]; primVertexAOD->GetXYZ(pos);
888   Double_t dca[6]={0.,0.,0.,0.,0.,0.}; //  modify it
889   Double_t dist12=TMath::Sqrt((vertexp1n1->GetX()-pos[0])*(vertexp1n1->GetX()-pos[0])+(vertexp1n1->GetY()-pos[1])*(vertexp1n1->GetY()-pos[1])+(vertexp1n1->GetZ()-pos[2])*(vertexp1n1->GetZ()-pos[2]));
890   Double_t dist23=TMath::Sqrt((vertexp2n1->GetX()-pos[0])*(vertexp2n1->GetX()-pos[0])+(vertexp2n1->GetY()-pos[1])*(vertexp2n1->GetY()-pos[1])+(vertexp2n1->GetZ()-pos[2])*(vertexp2n1->GetZ()-pos[2]));
891   Double_t dist14=0.; // to be implemented
892   Double_t dist34=0.; // to be implemented
893   Short_t charge=0;
894   AliAODRecoDecayHF4Prong *the4Prong = new AliAODRecoDecayHF4Prong(secVert,px,py,pz,d0,d0err,dca,dist12,dist23,dist14,dist34,charge);
895   the4Prong->SetOwnPrimaryVtx(primVertexAOD);
896   UShort_t id[4]={(UShort_t)postrack1->GetID(),(UShort_t)negtrack1->GetID(),(UShort_t)postrack2->GetID(),(UShort_t)negtrack2->GetID()};
897   the4Prong->SetProngIDs(4,id);
898
899
900   // use the following two lines once AliAODRecoDecayHF4Prong::SelectD0 is available
901   // select D0->Kpipipi
902   //Int_t checkD0,checkD0bar;   
903   // ok4Prong=the4Prong->SelectD0(fD04pCuts,checkD0,checkD0bar); 
904   ok4Prong=kFALSE;  //for the time being ...
905
906   if(!fRecoPrimVtxSkippingTrks && !fRmTrksFromPrimVtx) the4Prong->UnsetOwnPrimaryVtx();
907
908   // get PID info from ESD
909   Double_t esdpid0[5];
910   postrack1->GetESDpid(esdpid0);
911   Double_t esdpid1[5];
912   negtrack1->GetESDpid(esdpid1);
913   Double_t esdpid2[5];
914   postrack2->GetESDpid(esdpid2);
915   Double_t esdpid3[5];
916   negtrack2->GetESDpid(esdpid3);
917
918   Double_t esdpid[20];
919   for(Int_t i=0;i<5;i++) {
920     esdpid[i]    = esdpid0[i];
921     esdpid[5+i]  = esdpid1[i];
922     esdpid[10+i] = esdpid2[i];
923     esdpid[15+i] = esdpid3[i];
924   }
925   the4Prong->SetPID(4,esdpid);
926
927   return the4Prong;
928 }
929 //-----------------------------------------------------------------------------
930 AliAODVertex* AliAnalysisVertexingHF::PrimaryVertex(TObjArray *trkArray,
931                                                     AliVEvent *event) const
932 {
933   // Returns primary vertex to be used for this candidate
934  
935   AliESDVertex *vertexESD = 0;
936   AliAODVertex *vertexAOD = 0;
937
938
939   if(!fRecoPrimVtxSkippingTrks && !fRmTrksFromPrimVtx) { 
940     // primary vertex from the input event
941     
942     vertexESD = new AliESDVertex(*fV1);
943
944   } else {
945     // primary vertex specific to this candidate
946
947     Int_t nTrks = trkArray->GetEntriesFast();
948     AliVertexerTracks *vertexer = new AliVertexerTracks(event->GetMagneticField());
949
950     if(fRecoPrimVtxSkippingTrks) { 
951       // recalculating the vertex
952       
953       if(strstr(fV1->GetTitle(),"VertexerTracksWithConstraint")) {
954         Float_t diamondcovxy[3];
955         event->GetDiamondCovXY(diamondcovxy);
956         Double_t pos[3]={event->GetDiamondX(),event->GetDiamondY(),0.};
957         Double_t cov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,10.};
958         AliESDVertex *diamond = new AliESDVertex(pos,cov,1.,1);
959         vertexer->SetVtxStart(diamond);
960         delete diamond; diamond=NULL;
961         if(strstr(fV1->GetTitle(),"VertexerTracksWithConstraintOnlyFitter")) 
962           vertexer->SetOnlyFitter();
963       }
964       Int_t skipped[10];
965       AliESDtrack *t = 0;
966       for(Int_t i=0; i<nTrks; i++) {
967         t = (AliESDtrack*)trkArray->UncheckedAt(i);
968         skipped[i] = (Int_t)t->GetID();
969       }
970       vertexer->SetSkipTracks(nTrks,skipped);
971       vertexESD = (AliESDVertex*)vertexer->FindPrimaryVertex(event); 
972       
973     } else if(fRmTrksFromPrimVtx) { 
974       // removing the prongs tracks
975       
976       TObjArray rmArray(nTrks);
977       UShort_t *rmId = new UShort_t[nTrks];
978       AliESDtrack *esdTrack = 0;
979       AliESDtrack *t = 0;
980       for(Int_t i=0; i<nTrks; i++) {
981         t = (AliESDtrack*)trkArray->UncheckedAt(i);
982         esdTrack = new AliESDtrack(*t);
983         rmArray.AddLast(esdTrack);
984         if(esdTrack->GetID()>=0) {
985           rmId[i]=(UShort_t)esdTrack->GetID();
986         } else {
987           rmId[i]=9999;
988         }
989       }
990       Float_t diamondxy[2]={event->GetDiamondX(),event->GetDiamondY()};
991       vertexESD = vertexer->RemoveTracksFromVertex(fV1,&rmArray,rmId,diamondxy);
992       delete [] rmId; rmId=NULL;
993       rmArray.Delete();
994       
995     }
996
997     if(!vertexESD) return vertexAOD;
998     if(vertexESD->GetNContributors()<=0) { 
999       AliDebug(2,"vertexing failed"); 
1000       delete vertexESD; vertexESD=NULL;
1001       return vertexAOD;
1002     }
1003
1004     delete vertexer; vertexer=NULL;
1005
1006   }
1007
1008   // convert to AliAODVertex
1009   Double_t pos[3],cov[6],chi2perNDF;
1010   vertexESD->GetXYZ(pos); // position
1011   vertexESD->GetCovMatrix(cov); //covariance matrix
1012   chi2perNDF = vertexESD->GetChi2toNDF();
1013   delete vertexESD; vertexESD=NULL;
1014
1015   vertexAOD = new AliAODVertex(pos,cov,chi2perNDF);
1016
1017   return vertexAOD;
1018 }
1019 //-----------------------------------------------------------------------------
1020 void AliAnalysisVertexingHF::PrintStatus() const {
1021   // Print parameters being used
1022
1023   //printf("Preselections:\n");
1024   //   fTrackFilter->Dump();
1025   if(fSecVtxWithKF) {
1026     printf("Secondary vertex with Kalman filter package (AliKFParticle)\n");
1027   } else {
1028     printf("Secondary vertex with AliVertexerTracks\n");
1029   }
1030   if(fRecoPrimVtxSkippingTrks) printf("RecoPrimVtxSkippingTrks\n");
1031   if(fRmTrksFromPrimVtx) printf("RmTrksFromPrimVtx\n");
1032   if(fD0toKpi) {
1033     printf("Reconstruct D0->Kpi candidates with cuts:\n");
1034     printf("    |M-MD0| [GeV]    < %f\n",fD0toKpiCuts[0]);
1035     printf("    dca    [cm]  < %f\n",fD0toKpiCuts[1]);
1036     printf("    cosThetaStar     < %f\n",fD0toKpiCuts[2]);
1037     printf("    pTK     [GeV/c]    > %f\n",fD0toKpiCuts[3]);
1038     printf("    pTpi    [GeV/c]    > %f\n",fD0toKpiCuts[4]);
1039     printf("    |d0K|  [cm]  < %f\n",fD0toKpiCuts[5]);
1040     printf("    |d0pi| [cm]  < %f\n",fD0toKpiCuts[6]);
1041     printf("    d0d0  [cm^2] < %f\n",fD0toKpiCuts[7]);
1042     printf("    cosThetaPoint    > %f\n",fD0toKpiCuts[8]);
1043   }
1044   if(fDstar) {
1045     printf("    |M-MD*| [GeV]    < %f\n",fDstarCuts[0]);
1046     printf("    |M_Kpipi-M_Kpi-(MD*-MD0)| [GeV]  < %f\n",fDstarCuts[1]);
1047     printf("    pTpisoft [GeV/c]    > %f\n",fDstarCuts[2]);
1048     printf("    pTpisoft [GeV/c]    < %f\n",fDstarCuts[3]);
1049     printf("    Theta(pisoft,D0plane) < %f\n",fDstarCuts[4]);
1050     printf("Reconstruct D*->D0pi candidates with cuts:\n");
1051     printf("   D0 from D* cuts:\n");
1052     printf("    |M-MD0| [GeV]    < %f\n",fD0fromDstarCuts[0]);
1053     printf("    dca    [cm]  < %f\n",fD0fromDstarCuts[1]);
1054     printf("    cosThetaStar     < %f\n",fD0fromDstarCuts[2]);
1055     printf("    pTK     [GeV/c]    > %f\n",fD0fromDstarCuts[3]);
1056     printf("    pTpi    [GeV/c]    > %f\n",fD0fromDstarCuts[4]);
1057     printf("    |d0K|  [cm]  < %f\n",fD0fromDstarCuts[5]);
1058     printf("    |d0pi| [cm]  < %f\n",fD0fromDstarCuts[6]);
1059     printf("    d0d0  [cm^2] < %f\n",fD0fromDstarCuts[7]);
1060     printf("    cosThetaPoint    > %f\n",fD0fromDstarCuts[8]);
1061   }
1062   if(fJPSItoEle) {
1063     printf("Reconstruct J/psi from B candidates with cuts:\n");
1064     printf("    |M-MJPSI| [GeV]    < %f\n",fBtoJPSICuts[0]);
1065     printf("    dca    [cm]  < %f\n",fBtoJPSICuts[1]);
1066     printf("    cosThetaStar     < %f\n",fBtoJPSICuts[2]);
1067     printf("    pTP     [GeV/c]    > %f\n",fBtoJPSICuts[3]);
1068     printf("    pTN    [GeV/c]    > %f\n",fBtoJPSICuts[4]);
1069     printf("    |d0P|  [cm]  < %f\n",fBtoJPSICuts[5]);
1070     printf("    |d0N| [cm]  < %f\n",fBtoJPSICuts[6]);
1071     printf("    d0d0  [cm^2] < %f\n",fBtoJPSICuts[7]);
1072     printf("    cosThetaPoint    > %f\n",fBtoJPSICuts[8]);
1073   }
1074   if(f3Prong) {
1075     printf("Reconstruct 3 prong candidates.\n");
1076     printf("  D+->Kpipi cuts:\n");
1077     printf("    |M-MD+| [GeV]    < %f\n",fDplusCuts[0]);
1078     printf("    pTK     [GeV/c]    > %f\n",fDplusCuts[1]);
1079     printf("    pTPi    [GeV/c]    > %f\n",fDplusCuts[2]);
1080     printf("    |d0K|  [cm]  > %f\n",fDplusCuts[3]);
1081     printf("    |d0Pi| [cm]  > %f\n",fDplusCuts[4]);
1082     printf("    dist12    [cm]  < %f\n",fDplusCuts[5]);
1083     printf("    sigmavert [cm]   < %f\n",fDplusCuts[6]);
1084     printf("    dist prim-sec [cm] > %f\n",fDplusCuts[7]);
1085     printf("    pM=Max{pT1,pT2,pT3} [GeV/c] > %f\n",fDplusCuts[8]);
1086     printf("    cosThetaPoint    > %f\n",fDplusCuts[9]);
1087     printf("    Sum d0^2 [cm^2]  > %f\n",fDplusCuts[10]);
1088     printf("    dca cut [cm]  < %f\n",fDplusCuts[11]);
1089     printf("  Ds->KKpi cuts:\n");
1090     printf("    |M-MDs| [GeV]    < %f\n",fDsCuts[0]);
1091     printf("    pTK     [GeV/c]    > %f\n",fDsCuts[1]);
1092     printf("    pTPi    [GeV/c]    > %f\n",fDsCuts[2]);
1093     printf("    |d0K|  [cm]  > %f\n",fDsCuts[3]);
1094     printf("    |d0Pi| [cm]  > %f\n",fDsCuts[4]);
1095     printf("    dist12    [cm]  < %f\n",fDsCuts[5]);
1096     printf("    sigmavert [cm]   < %f\n",fDsCuts[6]);
1097     printf("    dist prim-sec [cm] > %f\n",fDsCuts[7]);
1098     printf("    pM=Max{pT1,pT2,pT3} [GeV/c] > %f\n",fDsCuts[8]);
1099     printf("    cosThetaPoint    > %f\n",fDsCuts[9]);
1100     printf("    Sum d0^2 [cm^2]  > %f\n",fDsCuts[10]);
1101     printf("    dca cut [cm]  < %f\n",fDsCuts[11]);
1102     printf("    Inv. Mass  phi/K0* [GeV]  < %f\n",fDsCuts[12]);
1103     printf("  Lc->pKpi cuts:\n");
1104     printf("    |M-MLc| [GeV]    < %f\n",fLcCuts[0]);
1105     printf("    pTP     [GeV/c]    > %f\n",fLcCuts[1]);
1106     printf("    pTPi and pTK [GeV/c]    > %f\n",fLcCuts[2]);
1107     printf("    |d0P|  [cm]  > %f\n",fLcCuts[3]);
1108     printf("    |d0Pi| and |d0K| [cm]  > %f\n",fLcCuts[4]);
1109     printf("    dist12    [cm]  < %f\n",fLcCuts[5]);
1110     printf("    sigmavert [cm]   < %f\n",fLcCuts[6]);
1111     printf("    dist prim-sec [cm] > %f\n",fLcCuts[7]);
1112     printf("    pM=Max{pT1,pT2,pT3} [GeV/c] > %f\n",fLcCuts[8]);
1113     printf("    cosThetaPoint    > %f\n",fLcCuts[9]);
1114     printf("    Sum d0^2 [cm^2]  > %f\n",fLcCuts[10]);
1115     printf("    dca cut [cm]  < %f\n",fLcCuts[11]);
1116     printf("  Ds->KKpi cuts:\n");
1117   }
1118
1119   return;
1120 }
1121 //-----------------------------------------------------------------------------
1122 AliAODVertex* AliAnalysisVertexingHF::ReconstructSecondaryVertex(TObjArray *trkArray,
1123                                                                  Double_t &dispersion) const
1124 {
1125   // Secondary vertex reconstruction with AliVertexerTracks or AliKFParticle
1126
1127   AliESDVertex *vertexESD = 0;
1128   AliAODVertex *vertexAOD = 0;
1129
1130   if(!fSecVtxWithKF) { // AliVertexerTracks
1131
1132     AliVertexerTracks *vertexer = new AliVertexerTracks(fBzkG);
1133     vertexer->SetVtxStart(fV1);
1134     vertexESD = (AliESDVertex*)vertexer->VertexForSelectedESDTracks(trkArray);
1135     delete vertexer; vertexer=NULL;
1136
1137     if(!vertexESD) return vertexAOD;
1138
1139     if(vertexESD->GetNContributors()!=trkArray->GetEntriesFast()) { 
1140       AliDebug(2,"vertexing failed"); 
1141       delete vertexESD; vertexESD=NULL;
1142       return vertexAOD;
1143     }
1144
1145   } else { // Kalman Filter vertexer (AliKFParticle)
1146
1147     AliKFParticle::SetField(fBzkG);
1148
1149     AliKFVertex vertexKF;
1150
1151     Int_t nTrks = trkArray->GetEntriesFast();
1152     for(Int_t i=0; i<nTrks; i++) {
1153       AliESDtrack *esdTrack = (AliESDtrack*)trkArray->At(i);
1154       AliKFParticle daughterKF(*esdTrack,211);
1155       vertexKF.AddDaughter(daughterKF);
1156     }
1157     vertexESD = new AliESDVertex(vertexKF.Parameters(),
1158                                  vertexKF.CovarianceMatrix(),
1159                                  vertexKF.GetChi2(),
1160                                  vertexKF.GetNContributors());
1161
1162   }
1163
1164   // convert to AliAODVertex
1165   Double_t pos[3],cov[6],chi2perNDF;
1166   vertexESD->GetXYZ(pos); // position
1167   vertexESD->GetCovMatrix(cov); //covariance matrix
1168   chi2perNDF = vertexESD->GetChi2toNDF();
1169   dispersion = vertexESD->GetDispersion();
1170   delete vertexESD; vertexESD=NULL;
1171
1172   vertexAOD = new AliAODVertex(pos,cov,chi2perNDF);
1173   vertexAOD->SetType(AliAODVertex::kUndef); // to be changed
1174
1175   return vertexAOD;
1176 }
1177 //-----------------------------------------------------------------------------
1178 Bool_t AliAnalysisVertexingHF::SelectInvMass(Int_t decay,
1179                                              Int_t nprongs,
1180                                              Double_t *px,
1181                                              Double_t *py,
1182                                              Double_t *pz) const {
1183   // Check invariant mass cut
1184
1185   Short_t dummycharge=0;
1186   Double_t *dummyd0 = new Double_t[nprongs];
1187   for(Int_t ip=0;ip<nprongs;ip++) dummyd0[ip]=0.;
1188   AliAODRecoDecay *rd = new AliAODRecoDecay(0x0,nprongs,dummycharge,px,py,pz,dummyd0);
1189   delete [] dummyd0;
1190
1191   UInt_t pdg2[2],pdg3[3];
1192   Double_t mPDG,minv;
1193
1194   Bool_t retval=kFALSE;
1195   switch (decay) 
1196     { 
1197     case 0:                  // D0->Kpi
1198       pdg2[0]=211; pdg2[1]=321;
1199       mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
1200       minv = rd->InvMass(nprongs,pdg2);
1201       if(TMath::Abs(minv-mPDG)<fD0toKpiCuts[0]) retval=kTRUE;
1202       pdg2[0]=321; pdg2[1]=211;
1203       minv = rd->InvMass(nprongs,pdg2);
1204       if(TMath::Abs(minv-mPDG)<fD0toKpiCuts[0]) retval=kTRUE;
1205       break;
1206     case 1:                  // JPSI->ee
1207       pdg2[0]=11; pdg2[1]=11;
1208       mPDG=TDatabasePDG::Instance()->GetParticle(443)->Mass();
1209       minv = rd->InvMass(nprongs,pdg2);
1210       if(TMath::Abs(minv-mPDG)<fBtoJPSICuts[0]) retval=kTRUE;
1211       break;
1212     case 2:                  // D+->Kpipi
1213       pdg3[0]=211; pdg3[1]=321; pdg3[2]=211;
1214       mPDG=TDatabasePDG::Instance()->GetParticle(411)->Mass();
1215       minv = rd->InvMass(nprongs,pdg3);
1216       if(TMath::Abs(minv-mPDG)<fDplusCuts[0]) retval=kTRUE;
1217                             // Ds+->KKpi
1218       pdg3[0]=321; pdg3[1]=321; pdg3[2]=211;
1219       mPDG=TDatabasePDG::Instance()->GetParticle(431)->Mass();
1220       minv = rd->InvMass(nprongs,pdg3);
1221       if(TMath::Abs(minv-mPDG)<fDsCuts[0]) retval=kTRUE;
1222       pdg3[0]=211; pdg3[1]=321; pdg3[2]=321;
1223       minv = rd->InvMass(nprongs,pdg3);
1224       if(TMath::Abs(minv-mPDG)<fDsCuts[0]) retval=kTRUE;
1225                             // Lc->pKpi
1226       pdg3[0]=2212; pdg3[1]=321; pdg3[2]=211;
1227       mPDG=TDatabasePDG::Instance()->GetParticle(4122)->Mass();
1228       minv = rd->InvMass(nprongs,pdg3);
1229       if(TMath::Abs(minv-mPDG)<fLcCuts[0]) retval=kTRUE;
1230       pdg3[0]=211; pdg3[1]=321; pdg3[2]=2212;
1231       minv = rd->InvMass(nprongs,pdg3);
1232       if(TMath::Abs(minv-mPDG)<fLcCuts[0]) retval=kTRUE; 
1233       break;
1234     case 3:                  // D*->D0pi
1235       pdg2[0]=421; pdg2[1]=211;
1236       mPDG=TDatabasePDG::Instance()->GetParticle(413)->Mass();
1237       minv = rd->InvMass(nprongs,pdg2);
1238       if(TMath::Abs(minv-mPDG)<fDstarCuts[0]) retval=kTRUE;
1239       break;
1240     default:
1241       printf("SelectInvMass(): wrong decay selection\n");
1242       break;
1243     }
1244
1245   delete rd;
1246
1247   return retval;
1248 }
1249 //-----------------------------------------------------------------------------
1250 void AliAnalysisVertexingHF::SelectTracksAndCopyVertex(AliVEvent *event,
1251                                    TObjArray &seleTrksArray,Int_t &nSeleTrks,
1252                                                        UChar_t *seleFlags)
1253 {
1254   // Apply single-track preselection.
1255   // Fill a TObjArray with selected tracks (for displaced vertices or
1256   // soft pion from D*). Selection flag stored in seleFlags.
1257   // Create the AliESDVertex object (convert from AliAODVertex if necessary)
1258   // In case of AOD input, also fill fAODMap for track index<->ID
1259
1260   const AliVVertex *vprimary = event->GetPrimaryVertex();
1261
1262   if(fV1) { delete fV1; fV1=NULL; }
1263   if(fAODMap) { delete [] fAODMap; fAODMap=NULL; }
1264
1265   Int_t nindices=0;
1266   UShort_t *indices = 0;
1267   Double_t pos[3],cov[6];
1268
1269   if(!fInputAOD) { // ESD
1270     fV1 = new AliESDVertex(*((AliESDVertex*)vprimary));
1271   } else {         // AOD
1272     vprimary->GetXYZ(pos);
1273     vprimary->GetCovarianceMatrix(cov);
1274     fV1 = new AliESDVertex(pos,cov,100.,100,vprimary->GetName());
1275     indices = new UShort_t[event->GetNumberOfTracks()];
1276     fAODMap = new Int_t[100000];
1277   }
1278
1279
1280   Int_t entries = (Int_t)event->GetNumberOfTracks();
1281   Bool_t okDisplaced=kFALSE,okSoftPi=kFALSE;
1282   nSeleTrks=0;
1283  
1284   // transfer ITS tracks from event to arrays
1285   for(Int_t i=0; i<entries; i++) {
1286     AliVTrack *track = (AliVTrack*)event->GetTrack(i);
1287
1288     if(fInputAOD) {
1289       AliAODTrack *aodt = (AliAODTrack*)track;
1290       if(aodt->GetUsedForPrimVtxFit()) { 
1291         indices[nindices]=aodt->GetID(); nindices++; 
1292       }
1293       fAODMap[(Int_t)aodt->GetID()] = i;
1294     }
1295
1296     AliESDtrack *esdt = 0;
1297     if(!fInputAOD) {
1298       esdt = (AliESDtrack*)track;
1299     } else {
1300       esdt = new AliESDtrack(track);
1301     }
1302
1303     // single track selection
1304     okDisplaced=kFALSE; okSoftPi=kFALSE;
1305     if(SingleTrkCuts(esdt,okDisplaced,okSoftPi)) {
1306       seleTrksArray.AddLast(esdt);
1307       seleFlags[nSeleTrks]=0;
1308       if(okDisplaced) SETBIT(seleFlags[nSeleTrks],kBitDispl);
1309       if(okSoftPi)    SETBIT(seleFlags[nSeleTrks],kBitSoftPi);
1310       nSeleTrks++;
1311     } else {
1312       if(fInputAOD) delete esdt; 
1313       esdt = NULL;
1314       continue;
1315     } 
1316
1317   } // end loop on tracks
1318
1319   // primary vertex from AOD
1320   if(fInputAOD) {
1321     delete fV1; fV1=NULL;
1322     vprimary->GetXYZ(pos);
1323     vprimary->GetCovarianceMatrix(cov);
1324     Double_t chi2toNDF = vprimary->GetChi2perNDF();
1325     Int_t ncontr=nindices;
1326     if(!strcmp(vprimary->GetTitle(),"VertexerTracksWithContraint")) ncontr += 1;
1327     Double_t chi2=chi2toNDF*(2.*(Double_t)ncontr-3.); 
1328     fV1 = new AliESDVertex(pos,cov,chi2,ncontr,vprimary->GetName());
1329     fV1->SetTitle(vprimary->GetTitle());
1330     fV1->SetIndices(nindices,indices);
1331   }
1332   if(indices) { delete [] indices; indices=NULL; }
1333
1334
1335   return;
1336 }
1337 //-----------------------------------------------------------------------------
1338 void AliAnalysisVertexingHF::SetD0toKpiCuts(Double_t cut0,Double_t cut1,
1339                                    Double_t cut2,Double_t cut3,Double_t cut4,
1340                                    Double_t cut5,Double_t cut6,
1341                                    Double_t cut7,Double_t cut8) 
1342 {
1343   // Set the cuts for D0 selection
1344   fD0toKpiCuts[0] = cut0;
1345   fD0toKpiCuts[1] = cut1;
1346   fD0toKpiCuts[2] = cut2;
1347   fD0toKpiCuts[3] = cut3;
1348   fD0toKpiCuts[4] = cut4;
1349   fD0toKpiCuts[5] = cut5;
1350   fD0toKpiCuts[6] = cut6;
1351   fD0toKpiCuts[7] = cut7;
1352   fD0toKpiCuts[8] = cut8;
1353
1354   return;
1355 }
1356 //-----------------------------------------------------------------------------
1357 void AliAnalysisVertexingHF::SetD0toKpiCuts(const Double_t cuts[9]) 
1358 {
1359   // Set the cuts for D0 selection
1360
1361   for(Int_t i=0; i<9; i++) fD0toKpiCuts[i] = cuts[i];
1362
1363   return;
1364 }
1365 //-----------------------------------------------------------------------------
1366 void AliAnalysisVertexingHF::SetD0fromDstarCuts(Double_t cut0,Double_t cut1,
1367                                    Double_t cut2,Double_t cut3,Double_t cut4,
1368                                    Double_t cut5,Double_t cut6,
1369                                    Double_t cut7,Double_t cut8) 
1370 {
1371   // Set the cuts for D0 from D* selection
1372   fD0fromDstarCuts[0] = cut0;
1373   fD0fromDstarCuts[1] = cut1;
1374   fD0fromDstarCuts[2] = cut2;
1375   fD0fromDstarCuts[3] = cut3;
1376   fD0fromDstarCuts[4] = cut4;
1377   fD0fromDstarCuts[5] = cut5;
1378   fD0fromDstarCuts[6] = cut6;
1379   fD0fromDstarCuts[7] = cut7;
1380   fD0fromDstarCuts[8] = cut8;
1381
1382   return;
1383 }
1384 //-----------------------------------------------------------------------------
1385 void AliAnalysisVertexingHF::SetD0fromDstarCuts(const Double_t cuts[9]) 
1386 {
1387   // Set the cuts for D0 from D* selection
1388
1389   for(Int_t i=0; i<9; i++) fD0fromDstarCuts[i] = cuts[i];
1390
1391   return;
1392 }
1393 //-----------------------------------------------------------------------------
1394 void AliAnalysisVertexingHF::SetDstarCuts(Double_t cut0,Double_t cut1,
1395                                           Double_t cut2,Double_t cut3,
1396                                           Double_t cut4)
1397 {
1398   // Set the cuts for D* selection
1399   fDstarCuts[0] = cut0;
1400   fDstarCuts[1] = cut1;
1401   fDstarCuts[2] = cut2;
1402   fDstarCuts[3] = cut3;
1403   fDstarCuts[4] = cut4;
1404
1405   return;
1406 }
1407 //-----------------------------------------------------------------------------
1408 void AliAnalysisVertexingHF::SetDstarCuts(const Double_t cuts[5])
1409 {
1410   // Set the cuts for D* selection
1411
1412   for(Int_t i=0; i<5; i++) fDstarCuts[i] = cuts[i];
1413
1414   return;
1415 }
1416 //-----------------------------------------------------------------------------
1417 void AliAnalysisVertexingHF::SetBtoJPSICuts(Double_t cut0,Double_t cut1,
1418                                    Double_t cut2,Double_t cut3,Double_t cut4,
1419                                    Double_t cut5,Double_t cut6,
1420                                    Double_t cut7,Double_t cut8) 
1421 {
1422   // Set the cuts for J/psi from B selection
1423   fBtoJPSICuts[0] = cut0;
1424   fBtoJPSICuts[1] = cut1;
1425   fBtoJPSICuts[2] = cut2;
1426   fBtoJPSICuts[3] = cut3;
1427   fBtoJPSICuts[4] = cut4;
1428   fBtoJPSICuts[5] = cut5;
1429   fBtoJPSICuts[6] = cut6;
1430   fBtoJPSICuts[7] = cut7;
1431   fBtoJPSICuts[8] = cut8;
1432
1433   return;
1434 }
1435 //-----------------------------------------------------------------------------
1436 void AliAnalysisVertexingHF::SetBtoJPSICuts(const Double_t cuts[9]) 
1437 {
1438   // Set the cuts for J/psi from B selection
1439
1440   for(Int_t i=0; i<9; i++) fBtoJPSICuts[i] = cuts[i];
1441
1442   return;
1443 }
1444 //-----------------------------------------------------------------------------
1445 void AliAnalysisVertexingHF::SetDplusCuts(Double_t cut0,Double_t cut1,
1446                                    Double_t cut2,Double_t cut3,Double_t cut4,
1447                                    Double_t cut5,Double_t cut6,
1448                                    Double_t cut7,Double_t cut8,
1449                                    Double_t cut9,Double_t cut10,Double_t cut11)
1450 {
1451   // Set the cuts for Dplus->Kpipi selection
1452   fDplusCuts[0] = cut0;
1453   fDplusCuts[1] = cut1;
1454   fDplusCuts[2] = cut2;
1455   fDplusCuts[3] = cut3;
1456   fDplusCuts[4] = cut4;
1457   fDplusCuts[5] = cut5;
1458   fDplusCuts[6] = cut6;
1459   fDplusCuts[7] = cut7;
1460   fDplusCuts[8] = cut8;
1461   fDplusCuts[9] = cut9;
1462   fDplusCuts[10] = cut10;
1463   fDplusCuts[11] = cut11;
1464
1465   return;
1466 }
1467 //-----------------------------------------------------------------------------
1468 void AliAnalysisVertexingHF::SetDplusCuts(const Double_t cuts[12]) 
1469 {
1470   // Set the cuts for Dplus->Kpipi selection
1471
1472   for(Int_t i=0; i<12; i++) fDplusCuts[i] = cuts[i];
1473
1474   return;
1475 }
1476 //-----------------------------------------------------------------------------
1477 void AliAnalysisVertexingHF::SetDsCuts(Double_t cut0,Double_t cut1,
1478                                    Double_t cut2,Double_t cut3,Double_t cut4,
1479                                    Double_t cut5,Double_t cut6,
1480                                    Double_t cut7,Double_t cut8,
1481                                    Double_t cut9,Double_t cut10,
1482                                    Double_t cut11,Double_t cut12)
1483 {
1484   // Set the cuts for Ds->KKpi selection
1485   fDsCuts[0] = cut0;
1486   fDsCuts[1] = cut1;
1487   fDsCuts[2] = cut2;
1488   fDsCuts[3] = cut3;
1489   fDsCuts[4] = cut4;
1490   fDsCuts[5] = cut5;
1491   fDsCuts[6] = cut6;
1492   fDsCuts[7] = cut7;
1493   fDsCuts[8] = cut8;
1494   fDsCuts[9] = cut9;
1495   fDsCuts[10] = cut10;
1496   fDsCuts[11] = cut11;
1497   fDsCuts[12] = cut12;
1498
1499   return;
1500 }
1501 //-----------------------------------------------------------------------------
1502 void AliAnalysisVertexingHF::SetDsCuts(const Double_t cuts[13]) 
1503 {
1504   // Set the cuts for Ds->KKpi selection
1505
1506   for(Int_t i=0; i<13; i++) fDsCuts[i] = cuts[i];
1507
1508   return;
1509 }
1510 //-----------------------------------------------------------------------------
1511 void AliAnalysisVertexingHF::SetLcCuts(Double_t cut0,Double_t cut1,
1512                                    Double_t cut2,Double_t cut3,Double_t cut4,
1513                                    Double_t cut5,Double_t cut6,
1514                                    Double_t cut7,Double_t cut8,
1515                                    Double_t cut9,Double_t cut10,Double_t cut11)
1516 {
1517   // Set the cuts for Lc->pKpi selection
1518   fLcCuts[0] = cut0;
1519   fLcCuts[1] = cut1;
1520   fLcCuts[2] = cut2;
1521   fLcCuts[3] = cut3;
1522   fLcCuts[4] = cut4;
1523   fLcCuts[5] = cut5;
1524   fLcCuts[6] = cut6;
1525   fLcCuts[7] = cut7;
1526   fLcCuts[8] = cut8;
1527   fLcCuts[9] = cut9;
1528   fLcCuts[10] = cut10;
1529   fLcCuts[11] = cut11;
1530
1531   return;
1532 }
1533 //-----------------------------------------------------------------------------
1534 void AliAnalysisVertexingHF::SetLcCuts(const Double_t cuts[12]) 
1535 {
1536   // Set the cuts for Lc->pKpi selection
1537
1538   for(Int_t i=0; i<12; i++) fLcCuts[i] = cuts[i];
1539
1540   return;
1541 }
1542 //-----------------------------------------------------------------------------
1543 Bool_t AliAnalysisVertexingHF::SingleTrkCuts(AliESDtrack *trk,
1544                                              Bool_t &okDisplaced,Bool_t &okSoftPi) const 
1545 {
1546   // Check if track passes some kinematical cuts  
1547
1548   // this is needed to store the impact parameters
1549   trk->RelateToVertex(fV1,fBzkG,kVeryBig);
1550
1551   UInt_t selectInfo;
1552   //
1553   // Track selection, displaced tracks
1554   selectInfo = 0; 
1555   if(fTrackFilter) {
1556     selectInfo = fTrackFilter->IsSelected(trk);
1557   }
1558   if(selectInfo) okDisplaced=kTRUE;
1559   // Track selection, soft pions
1560   selectInfo = 0; 
1561   if(fDstar && fTrackFilterSoftPi) {
1562     selectInfo = fTrackFilterSoftPi->IsSelected(trk);
1563   }
1564   if(selectInfo) okSoftPi=kTRUE;
1565
1566   if(okDisplaced || okSoftPi) return kTRUE;
1567
1568   return kFALSE;
1569 }
1570 //-----------------------------------------------------------------------------