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