]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibAlign.h
Add OCDB FEE macro by Frederick
[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
21 class AliTPCcalibAlign:public AliTPCcalibBase {
22 public:
23   enum HistoType {kY=0, kZ =1, kPhi=2, kTheta=3, 
24                   kYPhi=4, kZTheta=5, 
25                   kYz=6,kZz=7,kPhiZ=8,kThetaZ=9};
26   enum FitType{ k6=0, k9=1, k12=2};
27   AliTPCcalibAlign();
28   AliTPCcalibAlign(const Text_t *name, const Text_t *title);
29   AliTPCcalibAlign(const AliTPCcalibAlign &align);
30   //
31   virtual ~AliTPCcalibAlign();
32   virtual void Process(AliTPCseed *track);
33   virtual void Analyze();
34   virtual void Terminate();  
35   virtual Long64_t Merge(TCollection* list);
36   //
37   virtual void EvalFitters();
38   TH1 * GetHisto(HistoType type, Int_t s1, Int_t s2, Bool_t force=kFALSE);
39   void  MakeTree(const char *fname="alignTree.root");
40   TGraphErrors * MakeGraph(Int_t sec0, Int_t sec1, Int_t dsec, 
41                            Int_t i0, Int_t i1, FitType type); 
42   void ProcessTracklets(const AliExternalTrackParam &t1,
43                         const AliExternalTrackParam &t2,
44                         const AliTPCseed * seed,
45                         Int_t s1,Int_t s2);
46   inline Int_t GetIndex(Int_t s1,Int_t s2){return 72*s1+s2;}
47   //
48   inline TLinearFitter* GetFitter12(Int_t s1,Int_t s2);
49   inline TLinearFitter* GetFitter9(Int_t s1,Int_t s2);
50   inline TLinearFitter* GetFitter6(Int_t s1,Int_t s2);
51   //
52   Bool_t GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a);
53   Bool_t GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a);
54   Bool_t GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a);
55   
56   void ProcessDiff(const AliExternalTrackParam &t1,
57                    const AliExternalTrackParam &t2,
58                    const AliTPCseed *seed,
59                    Int_t s1,Int_t s2);
60
61 //   Bool_t GetTransformationCovar12(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
62 //   Bool_t GetTransformationCovar9(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
63 //   Bool_t GetTransformationCovar6(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
64   void Add(AliTPCcalibAlign * align);
65   Int_t *GetPoints() {return fPoints;}
66   void     Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
67   void     Process(AliESDEvent *event){AliTPCcalibBase::Process(event);}
68 private:
69   void FillHisto(const AliExternalTrackParam &t1,
70                         const AliExternalTrackParam &t2,
71                         Int_t s1,Int_t s2);
72
73   void Process12(const Double_t *t1, const Double_t *t2,
74                  TLinearFitter *fitter);
75   void Process9(Double_t *t1, Double_t *t2, TLinearFitter *fitter);
76   void Process6(Double_t *t1, Double_t *t2, TLinearFitter *fitter);
77   TLinearFitter* GetOrMakeFitter12(Int_t s1,Int_t s2);
78   TLinearFitter* GetOrMakeFitter9(Int_t s1,Int_t s2);
79   TLinearFitter* GetOrMakeFitter6(Int_t s1,Int_t s2);
80   TObjArray fDphiHistArray;    // array of residual histograms  phi      -kPhi
81   TObjArray fDthetaHistArray;  // array of residual histograms  theta    -kTheta
82   TObjArray fDyHistArray;      // array of residual histograms  y        -kY
83   TObjArray fDzHistArray;      // array of residual histograms  z        -kZ
84   //
85   TObjArray fDyPhiHistArray;      // array of residual histograms  y     -kYPhi
86   TObjArray fDzThetaHistArray;    // array of residual histograms  z-z   -kZTheta
87   //
88   TObjArray fDphiZHistArray;      // array of residual histograms  phi   -kPhiz
89   TObjArray fDthetaZHistArray;    // array of residual histograms  theta -kThetaz
90   TObjArray fDyZHistArray;        // array of residual histograms  y     -kYz
91   TObjArray fDzZHistArray;        // array of residual histograms  z     -kZz
92   //
93   //
94   TObjArray fFitterArray12;    // array of fitters
95   TObjArray fFitterArray9;     // array of fitters
96   TObjArray fFitterArray6;     // array of fitters
97   Int_t fPoints[72*72];        // number of points in the fitter
98   ClassDef(AliTPCcalibAlign,1)
99 };
100
101
102 TLinearFitter* AliTPCcalibAlign::GetFitter12(Int_t s1,Int_t s2) {
103   return static_cast<TLinearFitter*>(fFitterArray12[GetIndex(s1,s2)]);
104 }
105 TLinearFitter* AliTPCcalibAlign::GetFitter9(Int_t s1,Int_t s2) {
106   return static_cast<TLinearFitter*>(fFitterArray9[GetIndex(s1,s2)]);
107 }
108 TLinearFitter* AliTPCcalibAlign::GetFitter6(Int_t s1,Int_t s2) {
109   return static_cast<TLinearFitter*>(fFitterArray6[GetIndex(s1,s2)]);
110 }
111
112
113
114
115
116 #endif