]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibAlign.h
Adding the processing from tree (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibAlign.h
1 #ifndef ALITPCCALIBALIGN_H
2 #define ALITPCCALIBALIGN_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 ////
8 ////
9 ////
10
11 #include "TObject.h"
12 #include "TObjArray.h"
13 #include "TLinearFitter.h"
14 #include "AliTPCcalibBase.h"
15 #include "TH1.h"
16
17 class AliExternalTrackParam;
18 class AliTPCseed;
19 class TGraphErrors;
20 class TTree;
21
22
23 class AliTPCcalibAlign:public AliTPCcalibBase {
24 public:
25   enum HistoType {kY=0, kZ =1, kPhi=2, kTheta=3, 
26                   kYPhi=4, kZTheta=5, 
27                   kYz=6,kZz=7,kPhiZ=8,kThetaZ=9};
28   enum FitType{ k6=0, k9=1, k12=2};
29   AliTPCcalibAlign();
30   AliTPCcalibAlign(const Text_t *name, const Text_t *title);
31   AliTPCcalibAlign(const AliTPCcalibAlign &align);
32   //
33   virtual ~AliTPCcalibAlign();
34   virtual void Process(AliTPCseed *track);
35   virtual void Analyze();
36   virtual void Terminate();  
37   virtual Long64_t Merge(TCollection* list);
38   //
39   virtual void EvalFitters();
40   TH1 * GetHisto(HistoType type, Int_t s1, Int_t s2, Bool_t force=kFALSE);
41   void  MakeTree(const char *fname="alignTree.root");
42   TGraphErrors * MakeGraph(Int_t sec0, Int_t sec1, Int_t dsec, 
43                            Int_t i0, Int_t i1, FitType type); 
44   void ProcessTracklets(const AliExternalTrackParam &t1,
45                         const AliExternalTrackParam &t2,
46                         const AliTPCseed * seed,
47                         Int_t s1,Int_t s2);
48   inline Int_t GetIndex(Int_t s1,Int_t s2){return 72*s1+s2;}
49   //
50   inline const TMatrixD     * GetTransformation(Int_t s1,Int_t s2, Int_t fitType);
51   //
52   inline TLinearFitter* GetFitter12(Int_t s1,Int_t s2);
53   inline TLinearFitter* GetFitter9(Int_t s1,Int_t s2);
54   inline TLinearFitter* GetFitter6(Int_t s1,Int_t s2);
55   //
56   Bool_t GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a);
57   Bool_t GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a);
58   Bool_t GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a);
59   Bool_t AcceptTracklet(const AliExternalTrackParam &tp1,
60                         const AliExternalTrackParam &tp2);
61
62   void ProcessDiff(const AliExternalTrackParam &t1,
63                    const AliExternalTrackParam &t2,
64                    const AliTPCseed *seed,
65                    Int_t s1,Int_t s2);
66   void ProcessAlign(Double_t * t1, Double_t * t2, Int_t s1,Int_t s2);
67
68 //   Bool_t GetTransformationCovar12(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
69 //   Bool_t GetTransformationCovar9(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
70 //   Bool_t GetTransformationCovar6(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
71   void Add(AliTPCcalibAlign * align);
72   Int_t *GetPoints() {return fPoints;}
73   void     Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
74   void     Process(AliESDEvent *event){AliTPCcalibBase::Process(event);}
75   TLinearFitter* GetOrMakeFitter12(Int_t s1,Int_t s2);
76   TLinearFitter* GetOrMakeFitter9(Int_t s1,Int_t s2);
77   TLinearFitter* GetOrMakeFitter6(Int_t s1,Int_t s2);
78   void Process12(const Double_t *t1, const Double_t *t2,
79                  TLinearFitter *fitter);
80   void Process9(Double_t *t1, Double_t *t2, TLinearFitter *fitter);
81   void Process6(Double_t *t1, Double_t *t2, TLinearFitter *fitter);
82   void ProcessTree(TTree * tree);
83   //
84   // For visualization and test purposes
85   //
86   Double_t Correct(Int_t type, Int_t value, Int_t s1, Int_t s2, Double_t x, Double_t y, Double_t z, Double_t phi,Double_t theta); 
87   static Double_t SCorrect(Int_t type, Int_t value, Int_t s1, Int_t s2, Double_t x, Double_t y, Double_t z, Double_t phi,Double_t theta){return Instance()->Correct(type,value,s1,s2,x,y,z,phi,theta);}
88   static AliTPCcalibAlign* Instance();
89   void SetInstance(AliTPCcalibAlign*param){fgInstance = param;}
90 private:
91   
92   void FillHisto(const AliExternalTrackParam &t1,
93                  const AliExternalTrackParam &t2,
94                  Int_t s1,Int_t s2);
95
96   TObjArray fDphiHistArray;    // array of residual histograms  phi      -kPhi
97   TObjArray fDthetaHistArray;  // array of residual histograms  theta    -kTheta
98   TObjArray fDyHistArray;      // array of residual histograms  y        -kY
99   TObjArray fDzHistArray;      // array of residual histograms  z        -kZ
100   //
101   TObjArray fDyPhiHistArray;      // array of residual histograms  y     -kYPhi
102   TObjArray fDzThetaHistArray;    // array of residual histograms  z-z   -kZTheta
103   //
104   TObjArray fDphiZHistArray;      // array of residual histograms  phi   -kPhiz
105   TObjArray fDthetaZHistArray;    // array of residual histograms  theta -kThetaz
106   TObjArray fDyZHistArray;        // array of residual histograms  y     -kYz
107   TObjArray fDzZHistArray;        // array of residual histograms  z     -kZz
108   //
109   //
110   TObjArray fFitterArray12;    // array of fitters
111   TObjArray fFitterArray9;     // array of fitters
112   TObjArray fFitterArray6;     // array of fitters
113   //
114   TObjArray fMatrixArray12;    // array of transnformtation matrix
115   TObjArray fMatrixArray9;     // array of transnformtation matrix
116   TObjArray fMatrixArray6;     // array of transnformtation matrix  
117   //
118   Int_t fPoints[72*72];        // number of points in the fitter 
119   static AliTPCcalibAlign*   fgInstance; //! Instance of this class (singleton implementation)
120   ClassDef(AliTPCcalibAlign,1)
121 };
122
123
124 TLinearFitter* AliTPCcalibAlign::GetFitter12(Int_t s1,Int_t s2) {
125   return static_cast<TLinearFitter*>(fFitterArray12[GetIndex(s1,s2)]);
126 }
127 TLinearFitter* AliTPCcalibAlign::GetFitter9(Int_t s1,Int_t s2) {
128   return static_cast<TLinearFitter*>(fFitterArray9[GetIndex(s1,s2)]);
129 }
130 TLinearFitter* AliTPCcalibAlign::GetFitter6(Int_t s1,Int_t s2) {
131   return static_cast<TLinearFitter*>(fFitterArray6[GetIndex(s1,s2)]);
132 }
133
134 const TMatrixD * AliTPCcalibAlign::GetTransformation(Int_t s1,Int_t s2, Int_t fitType){
135   if (fitType==0) return static_cast<TMatrixD*>(fMatrixArray6[GetIndex(s1,s2)]);
136   if (fitType==1) return static_cast<TMatrixD*>(fMatrixArray9[GetIndex(s1,s2)]);
137   if (fitType==2) return static_cast<TMatrixD*>(fMatrixArray12[GetIndex(s1,s2)]);
138 }
139
140
141
142 #endif