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