]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibCosmic.h
Setting compilaton directory to current one
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibCosmic.h
1 #ifndef ALITPCCALIBCOSMIC_H
2 #define ALITPCCALIBCOSMIC_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 #include "AliTPCcalibBase.h"
8 class TH2F;
9 class TH1F;
10 class TList;
11 class AliESDEvent;
12 class AliESDtrack;
13 class THnSparse;
14
15 class AliTPCcalibCosmic:public AliTPCcalibBase {
16 public:
17   AliTPCcalibCosmic(); 
18   AliTPCcalibCosmic(const Text_t *name, const Text_t *title);
19   virtual ~AliTPCcalibCosmic();
20   
21   virtual void      Process(AliESDEvent *event);
22   virtual Long64_t  Merge(TCollection *const li);
23   virtual void      Analyze();
24   void              Add(const AliTPCcalibCosmic* cosmic);
25   //
26   //
27   void              Init();
28   void              FindPairs(AliESDEvent *event);
29   void              FindCosmicPairs(AliESDEvent * event);
30
31   Bool_t            IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1);
32   static void       CalculateBetheParams(TH2F *hist, Double_t * initialParam);
33   static Double_t   CalculateMIPvalue(TH1F * hist);
34   AliExternalTrackParam *MakeTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1);
35   AliExternalTrackParam *MakeCombinedTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1);
36
37   void UpdateTrack(AliExternalTrackParam &track0, const AliExternalTrackParam &track1);
38   //
39   void FillHistoPerformance(const AliExternalTrackParam *par0, const AliExternalTrackParam *par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam *inner1, AliTPCseed *seed0,  AliTPCseed *seed1, const AliExternalTrackParam *param0Combined, Int_t cross);
40   void MaterialBudgetDump(AliExternalTrackParam *const par0, AliExternalTrackParam *const par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam *inner1, AliTPCseed *const seed0,  AliTPCseed *const seed1, AliExternalTrackParam *const param0Combined, AliExternalTrackParam *const param1Combined);
41   static void MakeFitTree(TTree * treeInput, TTreeSRedirector *pcstream, const TObjArray * corrArray, Int_t step, Int_t run);
42   TTree * GetCosmicTree() const {return fCosmicTree;}
43   //
44   TH1F   *          GetHistNTracks() const {return fHistNTracks;};
45   TH1F   *          GetHistClusters() const {return fClusters;};
46   TH2F   *          GetHistAcorde()const {return fModules;};
47   TH1F   *          GetHistPt() const {return fHistPt;};
48   TH2F   *          GetHistDeDx() const {return fDeDx;};
49   TH1F   *          GetHistMIP() const {return fDeDxMIP;};
50   //
51   Double_t          GetMIPvalue()const {return fMIPvalue;};
52   //
53   static void       BinLogX(TH1 *const h);   // method for correct histogram binning
54   static void       BinLogX(THnSparse *const h, Int_t axisDim);   // method for correct histogram binning
55
56   void     Process(AliESDtrack *const track, Int_t runNo=-1) {AliTPCcalibBase::Process(track,runNo);};
57   void     Process(AliTPCseed *const track)  {return AliTPCcalibBase::Process(track);}
58   virtual void  Terminate();
59   static Double_t GetDeltaTime(Double_t rmin0, Double_t rmax0, Double_t rmin1, Double_t rmax1, Double_t tmin0, Double_t tmax0, Double_t tmin1, Double_t tmax1, Double_t dcaR, TVectorD& vectorDT);
60 public:  
61   //
62   // Performance histograms
63   //
64   THnSparse   *fHistoDelta[6];  // histograms of tracking performance delta
65   THnSparse   *fHistoPull[6];   // histograms of tracking performance pull
66   THnSparse   *fHistodEdxMax[4];   // histograms of dEdx perfomance - max charge
67   THnSparse   *fHistodEdxTot[4];   // histograms of dEdx perfomance - tot charge
68   static void AddTree(TTree* treeOutput, TTree * treeInput);
69 private:
70   
71   void              FillAcordeHist(AliESDtrack *upperTrack);
72
73   
74
75   TH1F  *fHistNTracks;            //  histogram showing number of ESD tracks per event
76   TH1F  *fClusters;               //  histogram showing the number of clusters per track
77   TH2F  *fModules;                //  2d histogram of tracks which are propagated to the ACORDE scintillator array
78   TH1F  *fHistPt;                 //  Pt histogram of reconstructed tracks
79   TH2F  *fDeDx;                   //  dEdx spectrum showing the different particle types
80   TH1F  *fDeDxMIP;                //  TPC signal close to the MIP region of muons 0.4 < p < 0.45 GeV
81
82   Double_t fMIPvalue;             //  MIP value calculated via a fit to fDeDxMIP
83   //
84   // cuts
85   //
86   Float_t fCutMaxD;     // maximal distance in rfi ditection
87   Float_t fCutMaxDz;     // maximal distance in z ditection
88   Float_t fCutTheta;    // maximal distance in theta ditection
89   Float_t fCutMinDir;   // direction vector products
90
91   TTree  *fCosmicTree;  // tree with the cosmic tracks
92   AliTPCcalibCosmic(const AliTPCcalibCosmic&); 
93   AliTPCcalibCosmic& operator=(const AliTPCcalibCosmic&); 
94
95   ClassDef(AliTPCcalibCosmic, 3); 
96 };
97
98 #endif
99