ROOT logo
#ifndef ALIANALYSISVERTEXINGHF_H
#define ALIANALYSISVERTEXINGHF_H
/* Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice                               */

/* $Id$ */ 

//-------------------------------------------------------------------------
//                      Class AliAnalysisVertexingHF
//            Reconstruction of heavy-flavour decay candidates
//      
//  Origin: E.Bruna, G.E.Bruno, A.Dainese, F.Prino, R.Romita, X.M.Zhang
//  Contact: andrea.dainese@pd.infn.it
//-------------------------------------------------------------------------

#include <TNamed.h>
#include <TList.h>

#include "AliAnalysisFilter.h"
#include "AliESDtrackCuts.h"

class AliPIDResponse;
class AliESDVertex;
class AliAODRecoDecay;
class AliAODRecoDecayHF;
class AliAODRecoDecayHF2Prong;
class AliAODRecoDecayHF3Prong;
class AliAODRecoDecayHF4Prong;
class AliAODRecoCascadeHF;
class AliAnalysisFilter;
class AliRDHFCuts;
class AliRDHFCutsD0toKpi;
class AliRDHFCutsJpsitoee;
class AliRDHFCutsDplustoKpipi;
class AliRDHFCutsDstoKKpi;
class AliRDHFCutsLctopKpi;
class AliRDHFCutsLctoV0;
class AliRDHFCutsD0toKpipipi;
class AliRDHFCutsDStartoKpipi;
class AliESDtrack;
class AliVEvent;
class AliAODVertex;
class AliVertexerTracks;
class AliESDv0; 
class AliAODv0; 

//-----------------------------------------------------------------------------
class AliAnalysisVertexingHF : public TNamed {
 public:
  //
  AliAnalysisVertexingHF();
  AliAnalysisVertexingHF(const AliAnalysisVertexingHF& source);
  AliAnalysisVertexingHF& operator=(const AliAnalysisVertexingHF& source); 
  virtual ~AliAnalysisVertexingHF();

  void FindCandidates(AliVEvent *event,
		      TClonesArray *aodVerticesHFTClArr,
		      TClonesArray *aodD0toKpiTClArr,
		      TClonesArray *aodJPSItoEleTClArr,
		      TClonesArray *aodCharm3ProngTClArr,
		      TClonesArray *aodCharm4ProngTClArr,
		      TClonesArray *aodDstarTClArr,
		      TClonesArray *aodCascadesTClArr,
		      TClonesArray *aodLikeSign2ProngTClArr,
		      TClonesArray *aodLikeSign3ProngTClArr);

  TList* FillListOfCuts();
  void FixReferences(AliAODEvent *aod);  
  void PrintStatus() const;
  void SetSecVtxWithKF() { fSecVtxWithKF=kTRUE; }
  void SetD0toKpiOn() { fD0toKpi=kTRUE; }
  void SetD0toKpiOff() { fD0toKpi=kFALSE; }
  void SetJPSItoEleOn() { fJPSItoEle=kTRUE; }
  void SetJPSItoEleOff() { fJPSItoEle=kFALSE; }
  void Set3ProngOn() { f3Prong=kTRUE; }
  void Set3ProngOff() { f3Prong=kFALSE; }
  void Set4ProngOn() { f4Prong=kTRUE; }
  void Set4ProngOff() { f4Prong=kFALSE; }
  void SetDstarOn() { fDstar=kTRUE; }
  void SetDstarOff() { fDstar=kFALSE; }
  void SetCascadesOn() { fCascades=kTRUE; }
  void SetCascadesOff() { fCascades=kFALSE; }
  void SetLikeSignOn() { fLikeSign=kTRUE; }
  void SetLikeSignOff() {fLikeSign=kFALSE; fLikeSign3prong=kFALSE;}
  void SetLikeSign3prongOn() { fLikeSign=kTRUE; fLikeSign3prong=kTRUE; }
  void SetLikeSign3prongOff() { fLikeSign3prong=kFALSE; }
  void SetMixEventOn() { fMixEvent=kTRUE; }
  void SetMixEventOff() { fMixEvent=kFALSE; }
  void SetInputAOD() { fInputAOD=kTRUE; }
  Bool_t GetD0toKpi() const { return fD0toKpi; }
  Bool_t GetJPSItoEle() const { return fJPSItoEle; }
  Bool_t Get3Prong() const { return f3Prong; }
  Bool_t Get4Prong() const { return f4Prong; }
  Bool_t GetDstar()  const { return fDstar; }
  Bool_t GetCascades() const { return fCascades; }
  Bool_t GetLikeSign() const { return fLikeSign; }
  Bool_t GetLikeSign3prong() const { return fLikeSign3prong; }
  Bool_t GetMixEvent() const { return fMixEvent; }
  Bool_t GetInputAOD() const { return fInputAOD; }
  Bool_t GetRecoPrimVtxSkippingTrks() const {return fRecoPrimVtxSkippingTrks;}
  Bool_t GetRmTrksFromPrimVtx() const {return fRmTrksFromPrimVtx;}
  void SetFindVertexForDstar(Bool_t vtx=kTRUE) { fFindVertexForDstar=vtx; }
  void SetFindVertexForCascades(Bool_t vtx=kTRUE) { fFindVertexForCascades=vtx; }

  void  SetV0TypeForCascadeVertex(Int_t type) {fV0TypeForCascadeVertex = type;}
  Int_t GetV0TypeForCascadeVertex()           { return fV0TypeForCascadeVertex;}
  
  void SetRecoPrimVtxSkippingTrks() 
    { fRecoPrimVtxSkippingTrks=kTRUE; fRmTrksFromPrimVtx=kFALSE;}
  void UnsetRecoPrimVtxSkippingTrks()
    { fRecoPrimVtxSkippingTrks=kFALSE; fRmTrksFromPrimVtx=kFALSE;}
  void SetRmTrksFromPrimVtx() 
    {fRmTrksFromPrimVtx=kTRUE; fRecoPrimVtxSkippingTrks=kFALSE; }
  void SetTrackFilter(AliAnalysisFilter* trackF) {
    // switch off the TOF selection that cannot be applied with AODTracks 
    TList *l = (TList*)trackF->GetCuts();
    AliESDtrackCuts *tcuts = (AliESDtrackCuts*)l->FindObject("AliESDtrackCuts");
    if(tcuts->GetFlagCutTOFdistance()) tcuts->SetFlagCutTOFdistance(kFALSE);
    fTrackFilter = trackF; 
  }
  void SetTrackFilter2prongPbCentral(Float_t maxPercentile, AliAnalysisFilter* trackF) {
    // switch off the TOF selection that cannot be applied with AODTracks 
    TList *l = (TList*)trackF->GetCuts();
    AliESDtrackCuts *tcuts = (AliESDtrackCuts*)l->FindObject("AliESDtrackCuts");
    if(tcuts->GetFlagCutTOFdistance()) tcuts->SetFlagCutTOFdistance(kFALSE);
    fTrackFilter2prongCentral = trackF; 
    fMaxCentPercentileForTightCuts=maxPercentile;
  }
  void SetTrackFilter3prongPbCentral(Float_t maxPercentile, AliAnalysisFilter* trackF) {
    // switch off the TOF selection that cannot be applied with AODTracks 
    TList *l = (TList*)trackF->GetCuts();
    AliESDtrackCuts *tcuts = (AliESDtrackCuts*)l->FindObject("AliESDtrackCuts");
    if(tcuts->GetFlagCutTOFdistance()) tcuts->SetFlagCutTOFdistance(kFALSE);
    fTrackFilter3prongCentral = trackF; 
    fMaxCentPercentileForTightCuts=maxPercentile;
  }
  void SetTrackFilterSoftPi(AliAnalysisFilter* trackF) { 
    // switch off the TOF selection that cannot be applied with AODTracks 
    TList *l = (TList*)trackF->GetCuts();
    AliESDtrackCuts *tcuts = (AliESDtrackCuts*)l->FindObject("AliESDtrackCuts");
    if(tcuts->GetFlagCutTOFdistance()) tcuts->SetFlagCutTOFdistance(kFALSE);
    fTrackFilterSoftPi = trackF; 
  }
  AliAnalysisFilter* GetTrackFilter() const { return fTrackFilter; }
  AliAnalysisFilter* GetTrackFilterSoftPi() const { return fTrackFilterSoftPi; }
  void SetCutsD0toKpi(AliRDHFCutsD0toKpi* cuts) { fCutsD0toKpi = cuts; }
  AliRDHFCutsD0toKpi* GetCutsD0toKpi() const { return fCutsD0toKpi; }
  void SetCutsJpsitoee(AliRDHFCutsJpsitoee* cuts) { fCutsJpsitoee = cuts; }
  AliRDHFCutsJpsitoee* GetCutsJpsitoee() const { return fCutsJpsitoee; }
  void SetCutsDplustoKpipi(AliRDHFCutsDplustoKpipi* cuts) { fCutsDplustoKpipi = cuts; }
  AliRDHFCutsDplustoKpipi* GetCutsDplustoKpipi() const { return fCutsDplustoKpipi; }
  void SetCutsDstoKKpi(AliRDHFCutsDstoKKpi* cuts) { fCutsDstoKKpi = cuts; }
  AliRDHFCutsDstoKKpi* GetCutsDstoKKpi() const { return fCutsDstoKKpi; }
  void SetCutsLctopKpi(AliRDHFCutsLctopKpi* cuts) { fCutsLctopKpi = cuts; }
  AliRDHFCutsLctopKpi* GetCutsLctopKpi() const { return fCutsLctopKpi; }
  void SetCutsLctoV0(AliRDHFCutsLctoV0* cuts) { fCutsLctoV0 = cuts; }
  AliRDHFCutsLctoV0* GetCutsLctoV0() const { return fCutsLctoV0; }
  void SetCutsD0toKpipipi(AliRDHFCutsD0toKpipipi* cuts) { fCutsD0toKpipipi = cuts; }
  AliRDHFCutsD0toKpipipi* GetCutsD0toKpipipi() const { return fCutsD0toKpipipi; }
  void SetCutsDStartoKpipi(AliRDHFCutsDStartoKpipi* cuts) { fCutsDStartoKpipi = cuts; }
  AliRDHFCutsDStartoKpipi* GetCutsDStartoKpipi() const { return fCutsDStartoKpipi; }
  void SetMassCutBeforeVertexing(Bool_t flag) { fMassCutBeforeVertexing=flag; } 

  void SetMasses();
  Bool_t CheckCutsConsistency();

  void SetUseTPCPID(Bool_t opt=kTRUE){fUseTPCPID=opt;}
  void SetUseTOFPID(Bool_t opt=kTRUE){fUseTOFPID=opt;}
  void SetUseTPCPIDOnlyIfNoTOF(Bool_t opt=kTRUE){fUseTPCPIDOnlyIfNoTOF=opt;}
  void SetMaxMomForTPCPid(Double_t mom){fMaxMomForTPCPid=mom;}
  void SetnSigmaTPCforPionSel(Double_t nsl, Double_t nsh){
    fnSigmaTPCPionLow=nsl; fnSigmaTPCPionHi=nsh;}
  void SetnSigmaTOFforPionSel(Double_t nsl, Double_t nsh){
    fnSigmaTOFPionLow=nsl; fnSigmaTOFPionHi=nsh;}
  void SetnSigmaTPCforKaonSel(Double_t nsl, Double_t nsh){
    fnSigmaTPCKaonLow=nsl; fnSigmaTPCKaonHi=nsh;}
  void SetnSigmaTOFforKaonSel(Double_t nsl, Double_t nsh){
    fnSigmaTOFKaonLow=nsl; fnSigmaTOFKaonHi=nsh;}
  void SetnSigmaTPCforProtonSel(Double_t nsl, Double_t nsh){
    fnSigmaTPCProtonLow=nsl; fnSigmaTPCProtonHi=nsh;}
  void SetnSigmaTOFforProtonSel(Double_t nsl, Double_t nsh){
    fnSigmaTOFProtonLow=nsl; fnSigmaTOFProtonHi=nsh;}

  void SetUseKaonPIDfor3Prong(Bool_t opt=kTRUE){fUseKaonPIDfor3Prong=opt;}
  void SetNotUseProtonPIDforLambdaC(){fUsePIDforLc=0;}
  void SetUseProtonPIDforLambdaC(){fUsePIDforLc=1;}
  void SetUseProtonAndPionPIDforLambdaC(){fUsePIDforLc=2;}
  void SetUseKaonPIDforDs(Bool_t opt=kTRUE){fUseKaonPIDforDs=opt;}

  void SetNotUseProtonPIDforLambdaC2V0(){fUsePIDforLc2V0=kFALSE;} //clm
  void SetUseProtonPIDforLambdaC2V0(){fUsePIDforLc2V0=kTRUE;}     //clm



  void GetnSigmaTOFforPionSel(Double_t& minnsigma, Double_t& maxnsigma) const {
    minnsigma=fnSigmaTOFPionLow;maxnsigma=fnSigmaTOFPionHi;
  }
  void GetnSigmaTPCforPionSel(Double_t& minnsigma, Double_t& maxnsigma) const {
    minnsigma=fnSigmaTPCPionLow;maxnsigma=fnSigmaTPCPionHi;
  }
  void GetnSigmaTOFforKaonSel(Double_t& minnsigma, Double_t& maxnsigma) const {
    minnsigma=fnSigmaTOFKaonLow;maxnsigma=fnSigmaTOFKaonHi;
  }
  void GetnSigmaTPCforKaonSel(Double_t& minnsigma, Double_t& maxnsigma) const {
    minnsigma=fnSigmaTPCKaonLow;maxnsigma=fnSigmaTPCKaonHi; 
  }
  void GetnSigmaTOFforProtonSel(Double_t& minnsigma, Double_t& maxnsigma) const {
    minnsigma=fnSigmaTOFProtonLow;maxnsigma=fnSigmaTOFProtonHi;
  }
  void GetnSigmaTPCforProtonSel(Double_t& minnsigma, Double_t& maxnsigma) const {
    minnsigma=fnSigmaTPCProtonLow;maxnsigma=fnSigmaTPCProtonHi;
  }

  Bool_t GetUseTPCPID() const {return fUseTPCPID;}
  Bool_t GetUseTOFPID() const {return fUseTOFPID;}
  Bool_t GetUseTPCPIDOnlyIfNoTOF() const {return fUseTPCPIDOnlyIfNoTOF;}
  Double_t GetMaxMomForTPCPid() const {return fMaxMomForTPCPid;}

  Bool_t GetUseKaonPIDfor3Prong() const {return fUseKaonPIDfor3Prong;}
  Int_t GetUseProtonPIDforLambdaC() const {return fUsePIDforLc;}
  Bool_t GetUseKaonPIDforDs() const {return fUseKaonPIDforDs;}
  Bool_t GetUseProtonPIDforLambdaC2V0() const {return fUsePIDforLc2V0;}
 
  void SetPidResponse(AliPIDResponse* p){fPidResponse=p;}

  //
 private:
  //
  enum { kBitDispl = 0, kBitSoftPi = 1, kBit3Prong = 2, kBitPionCompat = 3, kBitKaonCompat = 4, kBitProtonCompat = 5};

  Bool_t fInputAOD; // input from AOD (kTRUE) or ESD (kFALSE) 
  Int_t fAODMapSize; // size of fAODMap 
  Int_t *fAODMap; //[fAODMapSize] map between index and ID for AOD tracks

  AliVertexerTracks* fVertexerTracks; // vertexer, to compute secondary vertices
  Double_t fBzkG; // z componenent of field in kG

  Bool_t fSecVtxWithKF; // if kTRUE use KF vertexer, else AliVertexerTracks

  Bool_t fRecoPrimVtxSkippingTrks; // flag for primary vertex reco on the fly
                                   // for each candidate, w/o its daughters
  Bool_t fRmTrksFromPrimVtx; // flag for fast removal of daughters from 
                             // the primary vertex

  AliESDVertex *fV1; // primary vertex

  // flag to enable candidates production
  Bool_t fD0toKpi;   // D0->Kpi 
  Bool_t fJPSItoEle; // Jpsi->ee
  Bool_t f3Prong;    // D+,Ds,Lc
  Bool_t f4Prong;    // D0->Kpipipi
  Bool_t fDstar;     // D*->D0pi
  Bool_t fCascades;  // cascades, Lc --> v0+track
  Bool_t fLikeSign;  // Like-sign pairs
  Bool_t fLikeSign3prong;  // Like-sign triplets
  Bool_t fMixEvent; // event mixing

  AliPIDResponse* fPidResponse; // PID response
  Bool_t fUseKaonPIDfor3Prong;  // Kaon PID usage for 3 prongs
  Int_t  fUsePIDforLc;          // PID for Lambdac: 0=no, 1=proton, 2=p and pi
  Bool_t fUsePIDforLc2V0;       // PID for Lambdac 2 V0: 0=no, 1=proton,
  Bool_t fUseKaonPIDforDs;      // Kaon PID usage for Ds
  Bool_t fUseTPCPID;            // switch use/not use TPC PID
  Bool_t fUseTOFPID;            // switch use/not use TOF PID
  Bool_t fUseTPCPIDOnlyIfNoTOF; // use TPC PID only for tracks that without TOF
  Double_t fMaxMomForTPCPid;    // upper momentum limit to apply TPC PID
  Double_t fnSigmaTPCPionLow;   //Low cut value on n. of sigmas for pi TPC PID
  Double_t fnSigmaTPCPionHi;    //High cut value on n. of sigmas for pi TPC PID
  Double_t fnSigmaTOFPionLow;   //Low cut value on n. of sigmas for pi TOF PID
  Double_t fnSigmaTOFPionHi;    //High cut value on n. of sigmas for pi TOF PID
  Double_t fnSigmaTPCKaonLow;   //Low cut value on n. of sigmas for K TPC PID
  Double_t fnSigmaTPCKaonHi;    //High cut value on n. of sigmas for K TPC PID
  Double_t fnSigmaTOFKaonLow;   //Low cut value on n. of sigmas for K TOF PID
  Double_t fnSigmaTOFKaonHi;    //High cut value on n. of sigmas for K TOF PID
  Double_t fnSigmaTPCProtonLow; //Low cut value on n. of sigmas for p TPC PID
  Double_t fnSigmaTPCProtonHi;  //High cut value on n. of sigmas for p TPC PID
  Double_t fnSigmaTOFProtonLow; //Low cut value on n. of sigmas for p TOF PID
  Double_t fnSigmaTOFProtonHi;  //High cut value on n. of sigmas for p TOF PID

  Float_t fMaxCentPercentileForTightCuts; //max. centrality percentile for using tight cuts

  // single-track cuts
  AliAnalysisFilter *fTrackFilter; //  Track Filter for displaced vertices
  AliAnalysisFilter *fTrackFilter2prongCentral; //  Track Filter for displaced vertices in PbPb central events (tighter cuts) for 2 prong (D0->Kpi)
  AliAnalysisFilter *fTrackFilter3prongCentral; //  Track Filter for displaced vertices in PbPb central events (tighter cuts) for 3 prong (D+, Ds, Lc)
  AliAnalysisFilter *fTrackFilterSoftPi; //  Track Filter for D* soft pion
  // candidates cuts
  AliRDHFCutsD0toKpi *fCutsD0toKpi; // D0->Kpi cuts
  AliRDHFCutsJpsitoee *fCutsJpsitoee; // J/psi->ee cuts
  AliRDHFCutsDplustoKpipi *fCutsDplustoKpipi; // D+->Kpipi cuts
  AliRDHFCutsDstoKKpi *fCutsDstoKKpi; // Ds->KKpi cuts
  AliRDHFCutsLctopKpi *fCutsLctopKpi; // Lc->pKpi cuts
  AliRDHFCutsLctoV0 *fCutsLctoV0; // Lc --> v0 + bachelor cuts
  AliRDHFCutsD0toKpipipi *fCutsD0toKpipipi; // D0->Kpipipi cuts
  AliRDHFCutsDStartoKpipi *fCutsDStartoKpipi; // Dstar->D0pi cuts

  TList *fListOfCuts;    // pointer to list of cuts for output file
  Bool_t fFindVertexForDstar; // reconstruct a secondary vertex or assume it's from the primary vertex
  Bool_t fFindVertexForCascades;  // reconstruct a secondary vertex or assume it's from the primary vertex
  Int_t  fV0TypeForCascadeVertex;  // Select which V0 type we want to use for the cascas
  Bool_t fMassCutBeforeVertexing; // to go faster in PbPb
  // dummies for invariant mass calculation
  AliAODRecoDecay *fMassCalc2; // for 2 prong
  AliAODRecoDecay *fMassCalc3; // for 3 prong
  AliAODRecoDecay *fMassCalc4; // for 4 prong
  Bool_t fOKInvMassD0; // pair fullfilling D0 inv mass selection
  Bool_t fOKInvMassJpsi; // pair fullfilling Jpsi inv mass selection
  Bool_t fOKInvMassDplus; // triplet fullfilling D+ inv mass selection
  Bool_t fOKInvMassDs; // triplet fullfilling Ds inv mass selection
  Bool_t fOKInvMassLc; // triplet fullfilling Lc inv mass selection
  Bool_t fOKInvMassDstar; // combination fullfilling D* inv mass selection
  Bool_t fOKInvMassD0to4p; // 4tracks fullfilling D0 inv mass selection
  Bool_t fOKInvMassLctoV0; // triplet fullfilling Lc inv mass selection

  Int_t fnTrksTotal;
  Int_t fnSeleTrksTotal;

  Double_t fMassDzero;
  Double_t fMassDplus;
  Double_t fMassDs;
  Double_t fMassLambdaC;
  Double_t fMassDstar;
  Double_t fMassJpsi;


  //
  void AddRefs(AliAODVertex *v,AliAODRecoDecayHF *rd,const AliVEvent *event,
	       const TObjArray *trkArray) const;
  void AddDaughterRefs(AliAODVertex *v,const AliVEvent *event,
		       const TObjArray *trkArray) const;
  AliAODRecoDecayHF2Prong* Make2Prong(TObjArray *twoTrackArray1,AliVEvent *event,
				      AliAODVertex *secVert,Double_t dcap1n1,
				      Bool_t &okD0,Bool_t &okJPSI,Bool_t &okD0fromDstar);
  AliAODRecoDecayHF3Prong* Make3Prong(TObjArray *threeTrackArray,AliVEvent *event,
				      AliAODVertex *secVert,
				      Double_t dispersion,
				      const AliAODVertex *vertexp1n1,
				      const AliAODVertex *vertexp2n1,
				      Double_t dcap1n1,Double_t dcap2n1,Double_t dcap1p2, 
				      Bool_t useForLc, Bool_t useForDs,
				      Bool_t &ok3Prong);
  AliAODRecoDecayHF4Prong* Make4Prong(TObjArray *fourTrackArray,AliVEvent *event,
                                      AliAODVertex *secVert,
                                      const AliAODVertex *vertexp1n1,
                                      const AliAODVertex *vertexp1n1p2,
                                      Double_t dcap1n1,Double_t dcap1n2,
                                      Double_t dcap2n1,Double_t dcap2n2,
                                      Bool_t &ok4Prong);
  AliAODRecoCascadeHF* MakeCascade(TObjArray *twoTrackArray,AliVEvent *event,
				   AliAODVertex *secVert,
				   AliAODRecoDecayHF2Prong *rd2Prong,
				   Double_t dca,
				   Bool_t &okDstar);
  AliAODRecoCascadeHF* MakeCascade(TObjArray *twoTrackArray,AliVEvent *event,
				   AliAODVertex *secVert,
				   AliAODv0 *v0,
				   Double_t dca,
				   Bool_t &okCascades);

  AliAODVertex* PrimaryVertex(const TObjArray *trkArray=0x0,AliVEvent *event=0x0) const;
  AliAODVertex* ReconstructSecondaryVertex(TObjArray *trkArray,Double_t &dispersion,Bool_t useTRefArray=kTRUE) const;

  Bool_t SelectInvMassAndPt3prong(Double_t *px,Double_t *py,Double_t *pz, Int_t pidLcStatus=3);
  Bool_t SelectInvMassAndPt4prong(Double_t *px,Double_t *py,Double_t *pz);
  Bool_t SelectInvMassAndPtD0Kpi(Double_t *px,Double_t *py,Double_t *pz);
  Bool_t SelectInvMassAndPtJpsiee(Double_t *px,Double_t *py,Double_t *pz);
  Bool_t SelectInvMassAndPtDstarD0pi(Double_t *px,Double_t *py,Double_t *pz);
  Bool_t SelectInvMassAndPtCascade(Double_t *px,Double_t *py,Double_t *pz);

  Bool_t SelectInvMassAndPt3prong(TObjArray *trkArray);
  Bool_t SelectInvMassAndPt4prong(TObjArray *trkArray);
  Bool_t SelectInvMassAndPtDstarD0pi(TObjArray *trkArray);

  void   SelectTracksAndCopyVertex(const AliVEvent *event,Int_t trkEntries,
				   TObjArray &seleTrksArray,
				   TObjArray &tracksAtVertex,
				   Int_t &nSeleTrks,
				   UChar_t *seleFlags,Int_t *evtNumber);
  void SetParametersAtVertex(AliESDtrack* esdt, const AliExternalTrackParam* extpar) const;

  Bool_t SingleTrkCuts(AliESDtrack *trk,Float_t centralityperc, Bool_t &okDisplaced,Bool_t &okSoftPi, Bool_t &ok3prong) const;

  void   SetSelectionBitForPID(AliRDHFCuts *cuts,AliAODRecoDecayHF *rd,Int_t bit);

  AliAODv0* TransformESDv0toAODv0(AliESDv0 *esdv0, 
				  TObjArray *twoTrackArrayV0);

  //
  ClassDef(AliAnalysisVertexingHF,24);  // Reconstruction of HF decay candidates
};


#endif








 AliAnalysisVertexingHF.h:1
 AliAnalysisVertexingHF.h:2
 AliAnalysisVertexingHF.h:3
 AliAnalysisVertexingHF.h:4
 AliAnalysisVertexingHF.h:5
 AliAnalysisVertexingHF.h:6
 AliAnalysisVertexingHF.h:7
 AliAnalysisVertexingHF.h:8
 AliAnalysisVertexingHF.h:9
 AliAnalysisVertexingHF.h:10
 AliAnalysisVertexingHF.h:11
 AliAnalysisVertexingHF.h:12
 AliAnalysisVertexingHF.h:13
 AliAnalysisVertexingHF.h:14
 AliAnalysisVertexingHF.h:15
 AliAnalysisVertexingHF.h:16
 AliAnalysisVertexingHF.h:17
 AliAnalysisVertexingHF.h:18
 AliAnalysisVertexingHF.h:19
 AliAnalysisVertexingHF.h:20
 AliAnalysisVertexingHF.h:21
 AliAnalysisVertexingHF.h:22
 AliAnalysisVertexingHF.h:23
 AliAnalysisVertexingHF.h:24
 AliAnalysisVertexingHF.h:25
 AliAnalysisVertexingHF.h:26
 AliAnalysisVertexingHF.h:27
 AliAnalysisVertexingHF.h:28
 AliAnalysisVertexingHF.h:29
 AliAnalysisVertexingHF.h:30
 AliAnalysisVertexingHF.h:31
 AliAnalysisVertexingHF.h:32
 AliAnalysisVertexingHF.h:33
 AliAnalysisVertexingHF.h:34
 AliAnalysisVertexingHF.h:35
 AliAnalysisVertexingHF.h:36
 AliAnalysisVertexingHF.h:37
 AliAnalysisVertexingHF.h:38
 AliAnalysisVertexingHF.h:39
 AliAnalysisVertexingHF.h:40
 AliAnalysisVertexingHF.h:41
 AliAnalysisVertexingHF.h:42
 AliAnalysisVertexingHF.h:43
 AliAnalysisVertexingHF.h:44
 AliAnalysisVertexingHF.h:45
 AliAnalysisVertexingHF.h:46
 AliAnalysisVertexingHF.h:47
 AliAnalysisVertexingHF.h:48
 AliAnalysisVertexingHF.h:49
 AliAnalysisVertexingHF.h:50
 AliAnalysisVertexingHF.h:51
 AliAnalysisVertexingHF.h:52
 AliAnalysisVertexingHF.h:53
 AliAnalysisVertexingHF.h:54
 AliAnalysisVertexingHF.h:55
 AliAnalysisVertexingHF.h:56
 AliAnalysisVertexingHF.h:57
 AliAnalysisVertexingHF.h:58
 AliAnalysisVertexingHF.h:59
 AliAnalysisVertexingHF.h:60
 AliAnalysisVertexingHF.h:61
 AliAnalysisVertexingHF.h:62
 AliAnalysisVertexingHF.h:63
 AliAnalysisVertexingHF.h:64
 AliAnalysisVertexingHF.h:65
 AliAnalysisVertexingHF.h:66
 AliAnalysisVertexingHF.h:67
 AliAnalysisVertexingHF.h:68
 AliAnalysisVertexingHF.h:69
 AliAnalysisVertexingHF.h:70
 AliAnalysisVertexingHF.h:71
 AliAnalysisVertexingHF.h:72
 AliAnalysisVertexingHF.h:73
 AliAnalysisVertexingHF.h:74
 AliAnalysisVertexingHF.h:75
 AliAnalysisVertexingHF.h:76
 AliAnalysisVertexingHF.h:77
 AliAnalysisVertexingHF.h:78
 AliAnalysisVertexingHF.h:79
 AliAnalysisVertexingHF.h:80
 AliAnalysisVertexingHF.h:81
 AliAnalysisVertexingHF.h:82
 AliAnalysisVertexingHF.h:83
 AliAnalysisVertexingHF.h:84
 AliAnalysisVertexingHF.h:85
 AliAnalysisVertexingHF.h:86
 AliAnalysisVertexingHF.h:87
 AliAnalysisVertexingHF.h:88
 AliAnalysisVertexingHF.h:89
 AliAnalysisVertexingHF.h:90
 AliAnalysisVertexingHF.h:91
 AliAnalysisVertexingHF.h:92
 AliAnalysisVertexingHF.h:93
 AliAnalysisVertexingHF.h:94
 AliAnalysisVertexingHF.h:95
 AliAnalysisVertexingHF.h:96
 AliAnalysisVertexingHF.h:97
 AliAnalysisVertexingHF.h:98
 AliAnalysisVertexingHF.h:99
 AliAnalysisVertexingHF.h:100
 AliAnalysisVertexingHF.h:101
 AliAnalysisVertexingHF.h:102
 AliAnalysisVertexingHF.h:103
 AliAnalysisVertexingHF.h:104
 AliAnalysisVertexingHF.h:105
 AliAnalysisVertexingHF.h:106
 AliAnalysisVertexingHF.h:107
 AliAnalysisVertexingHF.h:108
 AliAnalysisVertexingHF.h:109
 AliAnalysisVertexingHF.h:110
 AliAnalysisVertexingHF.h:111
 AliAnalysisVertexingHF.h:112
 AliAnalysisVertexingHF.h:113
 AliAnalysisVertexingHF.h:114
 AliAnalysisVertexingHF.h:115
 AliAnalysisVertexingHF.h:116
 AliAnalysisVertexingHF.h:117
 AliAnalysisVertexingHF.h:118
 AliAnalysisVertexingHF.h:119
 AliAnalysisVertexingHF.h:120
 AliAnalysisVertexingHF.h:121
 AliAnalysisVertexingHF.h:122
 AliAnalysisVertexingHF.h:123
 AliAnalysisVertexingHF.h:124
 AliAnalysisVertexingHF.h:125
 AliAnalysisVertexingHF.h:126
 AliAnalysisVertexingHF.h:127
 AliAnalysisVertexingHF.h:128
 AliAnalysisVertexingHF.h:129
 AliAnalysisVertexingHF.h:130
 AliAnalysisVertexingHF.h:131
 AliAnalysisVertexingHF.h:132
 AliAnalysisVertexingHF.h:133
 AliAnalysisVertexingHF.h:134
 AliAnalysisVertexingHF.h:135
 AliAnalysisVertexingHF.h:136
 AliAnalysisVertexingHF.h:137
 AliAnalysisVertexingHF.h:138
 AliAnalysisVertexingHF.h:139
 AliAnalysisVertexingHF.h:140
 AliAnalysisVertexingHF.h:141
 AliAnalysisVertexingHF.h:142
 AliAnalysisVertexingHF.h:143
 AliAnalysisVertexingHF.h:144
 AliAnalysisVertexingHF.h:145
 AliAnalysisVertexingHF.h:146
 AliAnalysisVertexingHF.h:147
 AliAnalysisVertexingHF.h:148
 AliAnalysisVertexingHF.h:149
 AliAnalysisVertexingHF.h:150
 AliAnalysisVertexingHF.h:151
 AliAnalysisVertexingHF.h:152
 AliAnalysisVertexingHF.h:153
 AliAnalysisVertexingHF.h:154
 AliAnalysisVertexingHF.h:155
 AliAnalysisVertexingHF.h:156
 AliAnalysisVertexingHF.h:157
 AliAnalysisVertexingHF.h:158
 AliAnalysisVertexingHF.h:159
 AliAnalysisVertexingHF.h:160
 AliAnalysisVertexingHF.h:161
 AliAnalysisVertexingHF.h:162
 AliAnalysisVertexingHF.h:163
 AliAnalysisVertexingHF.h:164
 AliAnalysisVertexingHF.h:165
 AliAnalysisVertexingHF.h:166
 AliAnalysisVertexingHF.h:167
 AliAnalysisVertexingHF.h:168
 AliAnalysisVertexingHF.h:169
 AliAnalysisVertexingHF.h:170
 AliAnalysisVertexingHF.h:171
 AliAnalysisVertexingHF.h:172
 AliAnalysisVertexingHF.h:173
 AliAnalysisVertexingHF.h:174
 AliAnalysisVertexingHF.h:175
 AliAnalysisVertexingHF.h:176
 AliAnalysisVertexingHF.h:177
 AliAnalysisVertexingHF.h:178
 AliAnalysisVertexingHF.h:179
 AliAnalysisVertexingHF.h:180
 AliAnalysisVertexingHF.h:181
 AliAnalysisVertexingHF.h:182
 AliAnalysisVertexingHF.h:183
 AliAnalysisVertexingHF.h:184
 AliAnalysisVertexingHF.h:185
 AliAnalysisVertexingHF.h:186
 AliAnalysisVertexingHF.h:187
 AliAnalysisVertexingHF.h:188
 AliAnalysisVertexingHF.h:189
 AliAnalysisVertexingHF.h:190
 AliAnalysisVertexingHF.h:191
 AliAnalysisVertexingHF.h:192
 AliAnalysisVertexingHF.h:193
 AliAnalysisVertexingHF.h:194
 AliAnalysisVertexingHF.h:195
 AliAnalysisVertexingHF.h:196
 AliAnalysisVertexingHF.h:197
 AliAnalysisVertexingHF.h:198
 AliAnalysisVertexingHF.h:199
 AliAnalysisVertexingHF.h:200
 AliAnalysisVertexingHF.h:201
 AliAnalysisVertexingHF.h:202
 AliAnalysisVertexingHF.h:203
 AliAnalysisVertexingHF.h:204
 AliAnalysisVertexingHF.h:205
 AliAnalysisVertexingHF.h:206
 AliAnalysisVertexingHF.h:207
 AliAnalysisVertexingHF.h:208
 AliAnalysisVertexingHF.h:209
 AliAnalysisVertexingHF.h:210
 AliAnalysisVertexingHF.h:211
 AliAnalysisVertexingHF.h:212
 AliAnalysisVertexingHF.h:213
 AliAnalysisVertexingHF.h:214
 AliAnalysisVertexingHF.h:215
 AliAnalysisVertexingHF.h:216
 AliAnalysisVertexingHF.h:217
 AliAnalysisVertexingHF.h:218
 AliAnalysisVertexingHF.h:219
 AliAnalysisVertexingHF.h:220
 AliAnalysisVertexingHF.h:221
 AliAnalysisVertexingHF.h:222
 AliAnalysisVertexingHF.h:223
 AliAnalysisVertexingHF.h:224
 AliAnalysisVertexingHF.h:225
 AliAnalysisVertexingHF.h:226
 AliAnalysisVertexingHF.h:227
 AliAnalysisVertexingHF.h:228
 AliAnalysisVertexingHF.h:229
 AliAnalysisVertexingHF.h:230
 AliAnalysisVertexingHF.h:231
 AliAnalysisVertexingHF.h:232
 AliAnalysisVertexingHF.h:233
 AliAnalysisVertexingHF.h:234
 AliAnalysisVertexingHF.h:235
 AliAnalysisVertexingHF.h:236
 AliAnalysisVertexingHF.h:237
 AliAnalysisVertexingHF.h:238
 AliAnalysisVertexingHF.h:239
 AliAnalysisVertexingHF.h:240
 AliAnalysisVertexingHF.h:241
 AliAnalysisVertexingHF.h:242
 AliAnalysisVertexingHF.h:243
 AliAnalysisVertexingHF.h:244
 AliAnalysisVertexingHF.h:245
 AliAnalysisVertexingHF.h:246
 AliAnalysisVertexingHF.h:247
 AliAnalysisVertexingHF.h:248
 AliAnalysisVertexingHF.h:249
 AliAnalysisVertexingHF.h:250
 AliAnalysisVertexingHF.h:251
 AliAnalysisVertexingHF.h:252
 AliAnalysisVertexingHF.h:253
 AliAnalysisVertexingHF.h:254
 AliAnalysisVertexingHF.h:255
 AliAnalysisVertexingHF.h:256
 AliAnalysisVertexingHF.h:257
 AliAnalysisVertexingHF.h:258
 AliAnalysisVertexingHF.h:259
 AliAnalysisVertexingHF.h:260
 AliAnalysisVertexingHF.h:261
 AliAnalysisVertexingHF.h:262
 AliAnalysisVertexingHF.h:263
 AliAnalysisVertexingHF.h:264
 AliAnalysisVertexingHF.h:265
 AliAnalysisVertexingHF.h:266
 AliAnalysisVertexingHF.h:267
 AliAnalysisVertexingHF.h:268
 AliAnalysisVertexingHF.h:269
 AliAnalysisVertexingHF.h:270
 AliAnalysisVertexingHF.h:271
 AliAnalysisVertexingHF.h:272
 AliAnalysisVertexingHF.h:273
 AliAnalysisVertexingHF.h:274
 AliAnalysisVertexingHF.h:275
 AliAnalysisVertexingHF.h:276
 AliAnalysisVertexingHF.h:277
 AliAnalysisVertexingHF.h:278
 AliAnalysisVertexingHF.h:279
 AliAnalysisVertexingHF.h:280
 AliAnalysisVertexingHF.h:281
 AliAnalysisVertexingHF.h:282
 AliAnalysisVertexingHF.h:283
 AliAnalysisVertexingHF.h:284
 AliAnalysisVertexingHF.h:285
 AliAnalysisVertexingHF.h:286
 AliAnalysisVertexingHF.h:287
 AliAnalysisVertexingHF.h:288
 AliAnalysisVertexingHF.h:289
 AliAnalysisVertexingHF.h:290
 AliAnalysisVertexingHF.h:291
 AliAnalysisVertexingHF.h:292
 AliAnalysisVertexingHF.h:293
 AliAnalysisVertexingHF.h:294
 AliAnalysisVertexingHF.h:295
 AliAnalysisVertexingHF.h:296
 AliAnalysisVertexingHF.h:297
 AliAnalysisVertexingHF.h:298
 AliAnalysisVertexingHF.h:299
 AliAnalysisVertexingHF.h:300
 AliAnalysisVertexingHF.h:301
 AliAnalysisVertexingHF.h:302
 AliAnalysisVertexingHF.h:303
 AliAnalysisVertexingHF.h:304
 AliAnalysisVertexingHF.h:305
 AliAnalysisVertexingHF.h:306
 AliAnalysisVertexingHF.h:307
 AliAnalysisVertexingHF.h:308
 AliAnalysisVertexingHF.h:309
 AliAnalysisVertexingHF.h:310
 AliAnalysisVertexingHF.h:311
 AliAnalysisVertexingHF.h:312
 AliAnalysisVertexingHF.h:313
 AliAnalysisVertexingHF.h:314
 AliAnalysisVertexingHF.h:315
 AliAnalysisVertexingHF.h:316
 AliAnalysisVertexingHF.h:317
 AliAnalysisVertexingHF.h:318
 AliAnalysisVertexingHF.h:319
 AliAnalysisVertexingHF.h:320
 AliAnalysisVertexingHF.h:321
 AliAnalysisVertexingHF.h:322
 AliAnalysisVertexingHF.h:323
 AliAnalysisVertexingHF.h:324
 AliAnalysisVertexingHF.h:325
 AliAnalysisVertexingHF.h:326
 AliAnalysisVertexingHF.h:327
 AliAnalysisVertexingHF.h:328
 AliAnalysisVertexingHF.h:329
 AliAnalysisVertexingHF.h:330
 AliAnalysisVertexingHF.h:331
 AliAnalysisVertexingHF.h:332
 AliAnalysisVertexingHF.h:333
 AliAnalysisVertexingHF.h:334
 AliAnalysisVertexingHF.h:335
 AliAnalysisVertexingHF.h:336
 AliAnalysisVertexingHF.h:337
 AliAnalysisVertexingHF.h:338
 AliAnalysisVertexingHF.h:339
 AliAnalysisVertexingHF.h:340
 AliAnalysisVertexingHF.h:341
 AliAnalysisVertexingHF.h:342
 AliAnalysisVertexingHF.h:343
 AliAnalysisVertexingHF.h:344
 AliAnalysisVertexingHF.h:345
 AliAnalysisVertexingHF.h:346
 AliAnalysisVertexingHF.h:347
 AliAnalysisVertexingHF.h:348
 AliAnalysisVertexingHF.h:349
 AliAnalysisVertexingHF.h:350
 AliAnalysisVertexingHF.h:351
 AliAnalysisVertexingHF.h:352
 AliAnalysisVertexingHF.h:353
 AliAnalysisVertexingHF.h:354
 AliAnalysisVertexingHF.h:355
 AliAnalysisVertexingHF.h:356
 AliAnalysisVertexingHF.h:357
 AliAnalysisVertexingHF.h:358
 AliAnalysisVertexingHF.h:359
 AliAnalysisVertexingHF.h:360
 AliAnalysisVertexingHF.h:361
 AliAnalysisVertexingHF.h:362
 AliAnalysisVertexingHF.h:363
 AliAnalysisVertexingHF.h:364
 AliAnalysisVertexingHF.h:365
 AliAnalysisVertexingHF.h:366
 AliAnalysisVertexingHF.h:367
 AliAnalysisVertexingHF.h:368
 AliAnalysisVertexingHF.h:369
 AliAnalysisVertexingHF.h:370
 AliAnalysisVertexingHF.h:371
 AliAnalysisVertexingHF.h:372
 AliAnalysisVertexingHF.h:373
 AliAnalysisVertexingHF.h:374
 AliAnalysisVertexingHF.h:375
 AliAnalysisVertexingHF.h:376
 AliAnalysisVertexingHF.h:377
 AliAnalysisVertexingHF.h:378
 AliAnalysisVertexingHF.h:379
 AliAnalysisVertexingHF.h:380
 AliAnalysisVertexingHF.h:381
 AliAnalysisVertexingHF.h:382
 AliAnalysisVertexingHF.h:383
 AliAnalysisVertexingHF.h:384
 AliAnalysisVertexingHF.h:385
 AliAnalysisVertexingHF.h:386
 AliAnalysisVertexingHF.h:387
 AliAnalysisVertexingHF.h:388
 AliAnalysisVertexingHF.h:389
 AliAnalysisVertexingHF.h:390
 AliAnalysisVertexingHF.h:391
 AliAnalysisVertexingHF.h:392
 AliAnalysisVertexingHF.h:393
 AliAnalysisVertexingHF.h:394
 AliAnalysisVertexingHF.h:395
 AliAnalysisVertexingHF.h:396
 AliAnalysisVertexingHF.h:397
 AliAnalysisVertexingHF.h:398
 AliAnalysisVertexingHF.h:399
 AliAnalysisVertexingHF.h:400
 AliAnalysisVertexingHF.h:401