]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliTwoTrackRes.h
V0 rescaling in MC and multiple bins in correct.C
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / AliTwoTrackRes.h
1 ////////////////////////////////////////////////////////////////////////////////
2 //345678901234567890123456789012345678901234567890123456789012345678901234567890
3 //       1         2         3         4         5         6         7         8
4 //
5 // Class AliTwoTrackRes 
6 // J. Mercado. Last modified: 28.12.2009
7 //
8 ////////////////////////////////////////////////////////////////////////////////
9
10 #ifndef ALITWOTRACKRES_H
11 #define ALITWOTRACKRES_H
12
13 #include <TNtuple.h>
14 #include <TObject.h>
15 #include <TFile.h>
16 #include <TTree.h>
17 #include "TMath.h"
18 #include "TBits.h"
19 #include "TVector3.h"
20 #include "TLorentzVector.h"
21 #include "AliAnalysisTask.h"
22 #include "AliESD.h"
23 #include "AliESDEvent.h"
24 #include "AliESDtrack.h"
25 #include "AliESDtrackCuts.h"
26 #include "AliExternalTrackParam.h"
27
28 //______________________________________________________________________________
29 class AliTwoTrackRes : public AliAnalysisTask {
30
31 public:
32   
33   AliTwoTrackRes() : AliAnalysisTask("",""), fChain(0), fESDEvent(0), fOutContainer(0), 
34     fTrackCuts(0), fNTuple1(0), fNTuple2(0), fP1(), fP2(), fPb1(), fPb2(), fP(), 
35     fQ(), fTpcEnt1(), fTpcEnt2(), fTpcDist() {}
36   AliTwoTrackRes(const char *name);
37   AliTwoTrackRes(const AliTwoTrackRes& aTwoTrackRes);
38   virtual ~AliTwoTrackRes();
39
40   AliTwoTrackRes& operator=(const AliTwoTrackRes& aTwoTrackRes);
41
42   virtual void ConnectInputData(Option_t *);
43   virtual void CreateOutputObjects();
44   virtual void Exec(Option_t *opt = "");
45   virtual void Terminate(Option_t *opt = "");
46
47   void   SetTr1(double pt1, double eta1, double phi1, double m);
48   void   SetTr2(double pt2, double eta2, double phi2, double m);
49   void   SetTpcEnt1(double x1, double y1, double z1);
50   void   SetTpcEnt2(double x2, double y2, double z2);
51   double Qinv2()         {fQ = fP2 - fP1; return -1.*fQ.M2();}
52   double Qinv()          {return TMath::Sqrt(TMath::Abs(Qinv2()));}
53   double MinDist(AliExternalTrackParam *trk1, AliExternalTrackParam *trk2);
54   int    GetNSha(TBits cl, TBits sh);
55   double Corr(TBits cl1, TBits cl2, TBits sh1, TBits sh2);
56   double Qfac(TBits cl1, TBits cl2, TBits sh1, TBits sh2);
57   double RotTr2Phi();
58   double Dist()    {fTpcDist = fTpcEnt2 - fTpcEnt1; return fTpcDist.Mag();}
59   double DEta()    const {return TMath::Abs(fP2.Eta()-fP1.Eta());}
60   double DTheta()  const {return fP2.Theta()-fP1.Theta();}
61   double DPhi()    const {return TVector2::Phi_mpi_pi(fP2.Phi()-fP1.Phi());}
62   void   NoSwap()        {fPb1 = fP1; fPb2 = fP2;}
63   void   Swap()          {fPb1 = fP2; fPb2 = fP1;}
64   void   FillNTuple1(double minsep, double sep, double corr, double qf, 
65                      int ns1, int ns2);
66   void   FillNTuple2(double minsep, double sep, double corr, double qf, 
67                      int ns1, int ns2);
68
69 private:
70
71   TTree           *fChain;        // Chain of ESD trees
72   AliESDEvent     *fESDEvent;     // Leaves
73   TObjArray       *fOutContainer; // Output data container
74   AliESDtrackCuts *fTrackCuts;    // Track cuts
75   TNtuple         *fNTuple1;      // True pairs
76   TNtuple         *fNTuple2;      // Mixed pairs
77   TLorentzVector  fP1;            // Four-momentum of track 1 (lab)
78   TLorentzVector  fP2;            // Four-momentum of track 2 (lab)
79   TLorentzVector  fPb1;           // Buffer single-track 1 for swapping
80   TLorentzVector  fPb2;           // Buffer single-track 2 for swapping
81   TLorentzVector  fP;             // Total four-momentum (lab)
82   TLorentzVector  fQ;             // Four-momentum difference (lab)
83   TVector3        fTpcEnt1;       // Nominal TPC entrance point track 1
84   TVector3        fTpcEnt2;       // Nominal TPC entrance point track 2
85   TVector3        fTpcDist;       // Nominal TPC entrance separation 
86
87   ClassDef(AliTwoTrackRes, 0);
88 };
89 #endif
90
91 //______________________________________________________________________________
92 // EOF
93