Importing the code for relative ITS-TPC alignment (Mikolaj)
[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 "TMath.h"
18 #include "TMatrixD.h"
19 #include "TMatrixDSym.h"
20 #include "TVectorD.h"
21 #include "TVectorD.h"
22 #include "TDecompLU.h"
23 #include "TVector3.h"
24 #include "TArrayI.h"
25
26 #include "AliESDtrack.h"
27 #include "AliTrackPointArray.h"
28 #include "AliGeomManager.h"
29 #include "AliTrackFitterKalman.h"
30 #include "AliTrackFitterRieman.h"
31 #include "AliESDfriendTrack.h"
32 #include "AliESDEvent.h"
33 #include "AliESDVertex.h"
34 #include "TH1D.h"
35 #include "TF1.h"
36
37 class AliRelAlignerKalman {
38
39 public:
40     AliRelAlignerKalman();
41
42     //User methods:
43     Bool_t AddESDTrack( AliESDtrack* pTrack );
44     Bool_t AddTrackPointArray( AliTrackPointArray* pTrackPointArr );
45     Bool_t AddCosmicEventSeparateTracking( AliESDEvent* pEvent );
46     
47     void Print();
48     void PrintCorrelationMatrix();
49
50     void GetMeasurementCov( TMatrixDSym& cov ) { cov = *fPMeasurementCov; }
51     void GetMeasurement( TVectorD& mes ) { mes = *fPMeasurement; }
52     void GetState( TVectorD& state ) { state = *fPX; }
53     void GetStateCov ( TMatrixDSym& cov ) { cov = *fPXcov; }
54     void GetSeed( TVectorD& seed, TMatrixDSym& seedCov ) { seed = *fPX; seedCov = *fPXcov; }
55     Double_t* GetStateArr() { return fPX->GetMatrixArray(); }
56     Double_t* GetStateCovArr() { return fPXcov->GetMatrixArray(); }
57     Double_t* GetMeasurementArr() { return fPMeasurement->GetMatrixArray(); }
58     Double_t* GetMeasurementCovArr() { return fPMeasurementCov->GetMatrixArray(); }
59     
60     Bool_t SetCalibrationPass( Bool_t cp=kTRUE );
61
62     //Experts only:
63     void SetOutRejSigma( const Double_t a=2. ) { fOutRejSigmas = a; }
64     void SetMeasurementCov( const TMatrixDSym& cov ) {*fPMeasurementCov = cov;}
65     void SetMeasurement( const TVectorD& mes ) {*fPMeasurement = mes;}
66     void SetState( const TVectorD& param ) {*fPX = param;}
67     void SetStateCov (const TMatrixDSym& cov ) {*fPXcov = cov;}
68     void SetSeed( const TVectorD& seed, const TMatrixDSym& seedCov ) {*fPX = seed; *fPXcov = seedCov; }
69     void SetTrackParams1( const TVector3& point, const TMatrixDSym& pointcov, const TVector3& dir, const TMatrixDSym& dircov );
70     void SetTrackParams2( const TVector3& point, const TMatrixDSym& pointcov, const TVector3& dir, const TMatrixDSym& dircov );
71     void SetTrackParams2b( const TVector3& point, const TMatrixDSym& pointcov, const TVector3& dir, const TMatrixDSym& dircov );
72     void SetRefSurface( const TVector3& point, const TVector3& dir );
73     void SetRefSurface( const AliTrackPoint& point, const Bool_t horizontal=kTRUE );
74     void SetRefSurface( const Double_t y );
75     Bool_t PrepareUpdate();
76     Bool_t Update();
77     Bool_t IsReady() {return fFitted;}
78     void SetQ( const Double_t Q = 1e-10 ) { fQ = Q; }
79     void Set3TracksMode( Bool_t mode=kTRUE );
80     void SetCosmicEvent( Bool_t cosmic=kTRUE ) {Set3TracksMode(cosmic);}
81     void PrintDebugInfo();
82     
83 protected:
84     AliRelAlignerKalman& operator= (const AliRelAlignerKalman& aligner );
85     AliRelAlignerKalman(const AliRelAlignerKalman&);
86
87     Bool_t UpdateEstimateKalman();
88     Bool_t FillMeasurement();
89     Bool_t FillMeasurementMatrix();
90     Bool_t PredictMeasurement( TVectorD& z, const TVectorD& x );
91     void GetThetaPhiCov( TMatrixDSym& cov, const TVector3& vec, const TMatrixDSym& vecCov );
92     void Angles( TVectorD &angles, const TMatrixD &rotmat );
93     void RotMat( TMatrixD& R, const TVectorD& angles );
94     Bool_t IntersectionLinePlane( TVector3& intersection, const TVector3& base, const TVector3& dir, const TVector3& p, const TVector3& n );
95     void TMatrixDSym_from_TMatrixD( TMatrixDSym& matsym, const TMatrixD& mat );
96     void SanitizeExtendedCovMatrix( TMatrixDSym& covout, const TMatrixDSym& pointcov, const TVector3& dir );
97     Bool_t ProcessTrackPointArrays();
98
99     Bool_t ExtractSubTrackPointArray( AliTrackPointArray* pTrackOut, TArrayI* pVolIDArrayOut, const AliTrackPointArray* pTrackIn, const Int_t firstlayer, const Int_t lastlayer );
100     Bool_t FitTrackPointArrayRieman( TVector3* pPoint, TMatrixDSym* pPointCov, TVector3* pDir, TMatrixDSym* pDirCov, AliTrackPointArray* pTrackPointArray, TArrayI* pVolIDs, AliTrackPoint* pRefPoint );
101     Bool_t FitTrackPointArrayKalman( TVector3* pPoint, TMatrixDSym* pPointCov, TVector3* pDir, TMatrixDSym* pDirCov, AliTrackPointArray* pTrackPointArray, TArrayI* pVolIDs, AliTrackPoint* pRefPoint );
102     Bool_t FindCosmicInEvent( AliESDEvent* pEvent );
103     void SortTrackPointArrayWithY( AliTrackPointArray* pout, const AliTrackPointArray* pin, const Bool_t descending=kFALSE );
104     void JoinTrackArrays( AliTrackPointArray* pout, const AliTrackPointArray* pin1, const AliTrackPointArray* pin2 );
105     Bool_t UpdateCalibration();
106
107 private:
108     AliTrackPointArray* fPTrackPointArray1;  //track point array in first volume (ITS)
109     AliTrackPointArray* fPTrackPointArray2;  //track point array in second volume (TPC)
110     AliTrackPointArray* fPTrackPointArray2b; //optionally second trackpoint array in second volume (TPC)
111     TArrayI* fPVolIDArray1;                 //array with VolID's of the points in the trackpointarray 
112     TArrayI* fPVolIDArray2;                 //array with VolID's of the points in the trackpointarray
113     TArrayI* fPVolIDArray2b;                //array with VolID's of the points in the trackpointarray
114     TVector3* fPDir1;  //track 1 direction
115     TVector3* fPDir2;  //track 2 direction
116     TVector3* fPDir2b; //track 2b direction (for cosmics)
117     TMatrixDSym* fPDir1Cov;     //Covariance matrix of dir vec track 1
118     TMatrixDSym* fPDir2Cov;     //Covariance matrix of dir vec track 2
119     TMatrixDSym* fPDir2bCov;    //Covariance matrix of dir vec track 2b
120     TMatrixDSym* fPDir1ThetaPhiCov; //Covariance matrix of the dir vector in spherical coord (angles only)
121     TMatrixDSym* fPDir2ThetaPhiCov; //Covariance matrix of the dir vector in spherical coord (angles only)
122     TMatrixDSym* fPDir2bThetaPhiCov; //Covariance matrix of the dir vector in spherical coord (angles only)
123     TVector3* fPPoint1;  //track 1 extrapolated point at reference surface
124     TVector3* fPPoint2;  //track 2 extrapolated point at reference surface
125     TVector3* fPPoint2b; //track 2b extrapol. (for cosmics)
126     TMatrixDSym* fPPoint1Cov; //covariance matrix of point 1
127     TMatrixDSym* fPPoint2Cov; //covariance matrix of point 2
128     TMatrixDSym* fPPoint2bCov; //covariance matrix of point 2b
129     TVector3* fPRefPoint;    //a point on the reference surface
130     TVector3* fPRefPointDir; //reference surfaces' direction vec
131     AliTrackPoint* fPRefAliTrackPoint; //reference point as an AliTrackPoint (orientation stored in cov matrix)
132     TMatrixD* fPRotMat; // system rotation matrix for estimated alignment angles
133     
134     //
135     //Kalman filter related stuff
136     //
137     TVectorD* fPX; //System (fit) parameters (phi, theta, psi, x, y, z, driftcorr, driftoffset )
138     TMatrixDSym* fPXcov; //covariance matrix of system parameters
139
140     TMatrixD* fPH;      //System measurement matrix
141     Double_t fQ;        //measure for system noise
142     
143     TVectorD* fPMeasurement; //the measurement vec for Kalman filter (theta,phi,x,z)
144     TMatrixDSym* fPMeasurementCov; //measurement vec cvariance
145
146     Int_t fNMeasurementParams; //how many numbers do we measure?
147     Int_t fNSystemParams;      //how many parameters do we fit?
148
149     Double_t fOutRejSigmas; //number of sigmas for outlier rejection
150     //
151     //
152     //
153
154     //Control and calibration histograms
155     TH1D* fPXMesHist;  //histo of x measurement
156     TH1D* fPZMesHist;  //histo of y measurement
157     TH1D* fPPhiMesHist; //histo of phi measurement
158     TH1D* fPThetaMesHist; //histo of theta measurement
159     TH1D* fPXMesHist2;  //histo of x measurement (3tracks mode)
160     TH1D* fPZMesHist2;  //histo of y measurement (3tracks mode)
161     TH1D* fPPhiMesHist2; //histo of phi measurement (3tracks mode)
162     TH1D* fPThetaMesHist2; //histo of theta measurement (3tracks mode)
163     TH1D* fPMesCov11Hist;  //histogram of the covariance of a fit parameter, used in calibration
164     TH1D* fPMesCov22Hist;  //histogram of the covariance of a fit parameter, used in calibration
165     TH1D* fPMesCov33Hist;  //histogram of the covariance of a fit parameter, used in calibration
166     TH1D* fPMesCov44Hist;  //histogram of the covariance of a fit parameter, used in calibration
167     TH1D* fPMesCov55Hist;  //histogram of the covariance of a fit parameter, used in calibration
168     TH1D* fPMesCov66Hist;  //histogram of the covariance of a fit parameter, used in calibration
169     TH1D* fPMesCov77Hist;  //histogram of the covariance of a fit parameter, used in calibration
170     TH1D* fPMesCov88Hist;  //histogram of the covariance of a fit parameter, used in calibration
171     TMatrixDSym* fPMeasurementCovCorr; //correction to be applied to the measurement covariance
172     
173     Bool_t f3TracksMode; //are we using 3 tracklets?
174     Bool_t fSortTrackPointsWithY; //whether to sort the points after processing
175     Bool_t fFixedMeasurementCovariance; //don't fiddle with measurement cov - supply it externally
176     Bool_t fCalibrationPass;
177     Bool_t fCalibrationUpdated;
178     
179     Bool_t fFitted; //did fitting finish without error?
180     Bool_t fCuts;    //track cuts?
181
182     Int_t fNTracks; //number of processed tracks
183     Int_t fNUpdates; //number of successful Kalman update
184     Int_t fNOutliers; //number of outliers
185     Int_t fNUnfitted; //number of trackpointarrays that were not fitted for some reason
186     Int_t fNMatchedCosmics; //number of matched tracklets for cosmics
187
188     Int_t fMinPointsVol1;   //mininum number of points in volume 1
189     Int_t fMinPointsVol2;   //mininum number of points in volume 2
190     Double_t fMinMom;       //min momentum of track for track cuts
191     Double_t fMaxMom;       //max momentum of track for track cuts
192     Double_t fMinAbsSinPhi; //more cuts
193     Double_t fMaxAbsSinPhi; //more cuts
194     Double_t fMinSinTheta;  //cuts
195     Double_t fMaxSinTheta;  //cuts
196     Double_t fMaxMatchingAngle; //cuts
197     Double_t fMaxMatchingDistance; //cuts
198     Int_t fFirstLayerVol1; //first layer of volume 1
199     Int_t fLastLayerVol1;  //last layer of volume 1
200     Int_t fFirstLayerVol2; //first layer of volume 2
201     Int_t fLastLayerVol2;  //last layer of volume 2
202
203     TVector3 * fPVec010; //vector pointing up
204
205     static const Int_t fgkNMeasurementParams2TrackMode = 4; //how many measurables in 2 track mode
206     static const Int_t fgkNMeasurementParams3TrackMode = 8; //how many measurables in 3 track mode
207     static const Int_t fgkNSystemParams = 8;                //how many fit parameters
208     
209     ClassDef(AliRelAlignerKalman,1)     //AliRelAlignerKalman class
210 };
211
212 #endif