]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRelAlignerKalman.h
fNEvents initialized in the member initialization list
[u/mrichter/AliRoot.git] / STEER / AliRelAlignerKalman.h
1 #ifndef ALIRELALIGNERKALMAN_H
2 #define ALIRELALIGNERKALMAN_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 //     Relative alignment of two tracking volumes (default ITS and TPC)
10 //     (see AliRelAlignerKalman.cxx for details)
11 //
12 //     Origin: Mikolaj Krzewicki, Nikhef, Mikolaj.Krzewicki@cern.ch
13 //
14 //////////////////////////////////////////////////////////////////////////////
15
16 #include <TMath.h>
17 #include <TObject.h>
18 #include <TVectorD.h>
19 #include <TMatrix.h>
20
21 class AliExternalTrackParam;
22 class AliESDEvent;
23 class AliESDtrack;
24
25 class AliRelAlignerKalman : public TObject {
26
27 public:
28     AliRelAlignerKalman();
29     virtual ~AliRelAlignerKalman();
30     AliRelAlignerKalman& operator= (const AliRelAlignerKalman& a );
31     AliRelAlignerKalman(const AliRelAlignerKalman& a);
32
33     //User methods:
34     Bool_t AddESDTrack(  const AliESDtrack* pTrack );
35     Bool_t AddCosmicEvent( const AliESDEvent* pEvent );
36     
37     void Print(Option_t* option="") const;
38
39     Double_t GetPsi() const {return (*fPX)(0);}
40     Double_t GetTheta() const {return (*fPX)(1);}
41     Double_t GetPhi() const {return (*fPX)(2);}
42     Double_t GetX() const {return (*fPX)(3);}
43     Double_t GetY() const {return (*fPX)(4);}
44     Double_t GetZ() const {return (*fPX)(5);}
45     Double_t GetTPCvdCorr() const {return (*fPX)(6);}
46     Double_t GetTPCvdY() const {return (*fPX)(7);}
47     Double_t GetTPCt0() const {return (*fPX)(8);}
48     Double_t GetPsiErr() const {return TMath::Sqrt((*fPXcov)(0,0));}
49     Double_t GetThetaErr() const {return TMath::Sqrt((*fPXcov)(1,1));}
50     Double_t GetPhiErr() const {return TMath::Sqrt((*fPXcov)(2,2));}
51     Double_t GetXErr() const {return TMath::Sqrt((*fPXcov)(3,3));}
52     Double_t GetYErr() const {return TMath::Sqrt((*fPXcov)(4,4));}
53     Double_t GetZErr() const {return TMath::Sqrt((*fPXcov)(5,5));}
54     Double_t GetTPCvdCorrErr() const {return TMath::Sqrt((*fPXcov)(6,6));}
55     Double_t GetTPCvdYErr() const {return TMath::Sqrt((*fPXcov)(7,7));}
56     Double_t GetTPCt0Err() const {return TMath::Sqrt((*fPXcov)(8,8));}
57     void GetMeasurement( TVectorD& mes ) const { mes = *fPMeasurement; }
58     void GetMeasurementCov( TMatrixDSym& cov ) const { cov = *fPMeasurementCov; }
59     void GetState( TVectorD& state ) const { state = *fPX; }
60     void GetStateCov ( TMatrixDSym& cov ) const { cov = *fPXcov; }
61     void GetSeed( TVectorD& seed, TMatrixDSym& seedCov ) const { seed = *fPX; seedCov = *fPXcov; }
62     void SetMeasurement( const TVectorD& mes ) {*fPMeasurement = mes;}
63     void SetMeasurementCov( const TMatrixDSym& cov ) {*fPMeasurementCov = cov;}
64     void SetState( const TVectorD& param ) {*fPX = param;}
65     void SetStateCov (const TMatrixDSym& cov ) {*fPXcov = cov;}
66     void SetSeed( const TVectorD& seed, const TMatrixDSym& seedCov ) {*fPX = seed; *fPXcov = seedCov; }
67
68     //Expert methods:
69     Bool_t FindCosmicTrackletNumbersInEvent( TArrayI& outITStracksTArr, TArrayI& outTPCtracksTArr, const AliESDEvent* pEvent );
70     Bool_t PrepareUpdate();
71     Bool_t Update();
72     void SetRefSurface( const Double_t x, const Double_t alpha );
73     void PrintDebugInfo();
74     void PrintCorrelationMatrix();
75     void PrintCovarianceCorrection();
76     void PrintSystemMatrix();
77     void Reset();
78     void ResetCovariance( const Double_t number=0. );
79     void ResetTPCparamsCovariance( const Double_t number=0. );
80     Double_t* GetStateArr() const { return fPX->GetMatrixArray(); }
81     Double_t* GetStateCovArr() const { return fPXcov->GetMatrixArray(); }
82     Double_t* GetMeasurementArr() const { return fPMeasurement->GetMatrixArray(); }
83     Double_t* GetMeasurementCovArr() const { return fPMeasurementCov->GetMatrixArray(); }
84     const Double_t* GetDeltaArr() const {return fDelta;}
85     TH1D* GetMes0Hist() const {return fPMes0Hist;}
86     TH1D* GetMes1Hist() const {return fPMes1Hist;} 
87     TH1D* GetMes2Hist() const {return fPMes2Hist;}
88     TH1D* GetMes3Hist() const {return fPMes3Hist;}
89     TH1D* GetMesErr0Hist() const {return fPMesErr0Hist;}
90     TH1D* GetMesErr1Hist() const {return fPMesErr1Hist;}
91     TH1D* GetMesErr2Hist() const {return fPMesErr2Hist;}
92     TH1D* GetMesErr3Hist() const {return fPMesErr3Hist;}
93     Bool_t SetCalibrationMode( const Bool_t cp=kTRUE );
94     void SetCorrectionMode( const Bool_t mode=kTRUE ){ fCorrectionMode=mode; }
95     void SetApplyCovarianceCorrection( const Bool_t s=kTRUE ) {fApplyCovarianceCorrection = s;}
96     void SetOutRejSigma( const Double_t a=2. ) { fOutRejSigmas = a; }
97     void SetRejectOutliers( const Bool_t r=kTRUE ) {fRejectOutliers = r;}
98     void SetCovarianceCorrection( const TMatrixDSym& c ) {*fPMeasurementCovCorr = c;}
99     void GetCovarianceCorrection( TMatrixDSym& c ) {c=*fPMeasurementCovCorr;}
100     void SetTrackParams1( const AliExternalTrackParam* exparam );
101     void SetTrackParams2( const AliExternalTrackParam* exparam );
102     void SetMinPointsVol1( const Int_t min ) {fMinPointsVol1=min;}
103     void SetMinPointsVol2( const Int_t min ) {fMinPointsVol2=min;}
104     void SetRequireMatchInTPC( const Bool_t s=kTRUE ) {fRequireMatchInTPC = s;}
105     void SetNHistogramBins( const Int_t n ) {fNHistogramBins = n;}
106     void SetQ( const Double_t Q = 1e-10 ) { fQ = Q; }
107     void SetMaxMatchingDistance( const Double_t m ) {fMaxMatchingDistance=m;}
108     void SetMaxMatchingAngle( const Double_t m ) {fMaxMatchingAngle=m;}
109     void SetTPCvd( const Float_t v ) {fTPCvd=v;}
110     void SetTPCZLengthA( const Double_t l ) {fTPCZLengthA=l;}
111     void SetTPCZLengthC( const Double_t l ) {fTPCZLengthC=l;}
112     Bool_t CorrectTrack( AliExternalTrackParam* tr, const TVectorD& misalignment );
113     Bool_t MisalignTrack( AliExternalTrackParam* tr, const TVectorD& misalignment );
114     static void Angles( TVectorD &angles, const TMatrixD &rotmat );
115     static void RotMat( TMatrixD& R, const TVectorD& angles );
116     static void TMatrixDSymFromTMatrixD( TMatrixDSym& matsym, const TMatrixD& mat );
117     
118 protected:
119     Bool_t UpdateCalibration();
120     Bool_t UpdateEstimateKalman();
121     Bool_t PrepareMeasurement();
122     Bool_t PrepareSystemMatrix();
123     Bool_t PredictMeasurement( TVectorD& z, const TVectorD& x );
124     Bool_t IsOutlier( const TVectorD& update, const TMatrixDSym& covmatrix );
125     Bool_t CalculateCovarianceCorrection();
126
127 private:
128     static const Int_t fgkNMeasurementParams = 4;           //how many measurables
129     static const Int_t fgkNSystemParams = 9;                //how many fit parameters
130     
131     //Track parameters
132     Double_t fAlpha;       //rotation angle between the local and global coordinate system like in AliExternalTrackParam
133     Double_t fLocalX;      //local x coordinate of reference plane = r(global)
134     const AliExternalTrackParam* fkPTrackParam1;   //local track parameters (theta,phi,y,z)
135     const AliExternalTrackParam* fkPTrackParam2;   //local track parameters
136
137     //Kalman filter related stuff
138     TVectorD* fPX; //System (fit) parameters (phi, theta, psi, x, y, z, driftcorr, driftoffset )
139     TMatrixDSym* fPXcov; //covariance matrix of system parameters
140     TMatrixD* fPH;      //System measurement matrix
141     Double_t fQ;        //measure for system noise
142     TVectorD* fPMeasurement; //the measurement vec for Kalman filter (theta,phi,x,z)
143     TMatrixDSym* fPMeasurementCov; //measurement vec cvariance
144     Double_t fOutRejSigmas; //number of sigmas for outlier rejection
145     Double_t fDelta[fgkNSystemParams]; //array with differentials for calculating derivatives for every parameter(see PrepareSystemMatrix())
146
147     //Control
148     Bool_t fRejectOutliers; //whether to do outlier rejection in the Kalman filter
149     Bool_t fCalibrationMode;            //are we running in calibration mode?
150     Bool_t fFillHistograms;     //whether to fill the histograms with residuals
151     Bool_t fRequireMatchInTPC;  //when looking for a cosmic in event, require that TPC has 2 matching segments
152     Bool_t fApplyCovarianceCorrection;         //add the correction to the covariance of measurement
153     Bool_t fCuts;    //track cuts?
154     Int_t fMinPointsVol1;   //mininum number of points in volume 1
155     Int_t fMinPointsVol2;   //mininum number of points in volume 2
156     Double_t fMinMom;       //min momentum of track for track cuts
157     Double_t fMaxMom;       //max momentum of track for track cuts
158     Double_t fMaxMatchingAngle; //cuts
159     Double_t fMaxMatchingDistance; //cuts
160     Bool_t fCorrectionMode; //calculate corrective transform for TPC (or monitor actual TPC misal params)
161     
162     //Counters
163     Int_t fNTracks; //number of processed tracks
164     Int_t fNUpdates; //number of successful Kalman updates
165     Int_t fNOutliers; //number of outliers
166     Int_t fNMatchedCosmics; //number of cosmic events with matching tracklets (good cosmics)
167     Int_t fNMatchedTPCtracklets;//number of cosmic events with 2 matching TPC tracklets
168     Int_t fNProcessedEvents; //number of processed events
169
170     //Calibration histograms
171     Int_t fNHistogramBins;   //how many bins in control histograms
172     TH1D* fPMes0Hist;  //histo of x measurement
173     TH1D* fPMes1Hist;  //histo of y measurement
174     TH1D* fPMes2Hist; //histo of phi measurement
175     TH1D* fPMes3Hist; //histo of theta measurement
176     TH1D* fPMesErr0Hist;  //histogram of the covariance of a fit parameter, used in calibration
177     TH1D* fPMesErr1Hist;  //histogram of the covariance of a fit parameter, used in calibration
178     TH1D* fPMesErr2Hist;  //histogram of the covariance of a fit parameter, used in calibration
179     TH1D* fPMesErr3Hist;  //histogram of the covariance of a fit parameter, used in calibration
180     TMatrixDSym* fPMeasurementCovCorr; //correction to be added to the measurement covariance
181
182     //TPC stuff
183     Double_t fTPCvd; //TPC drift velocity
184     Double_t fTPCZLengthA; //TPC length side A
185     Double_t fTPCZLengthC; //TPC length side C
186     
187     ClassDef(AliRelAlignerKalman,2)     //AliRelAlignerKalman class
188 };
189
190 #endif
191