]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliRelAlignerKalman.h
added function returning number in truncated gaussian (user passes mean, sigma and...
[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     void GetMeasurement( TVectorD& mes ) { mes = *fPMeasurement; }
52     void GetMeasurementCov( TMatrixDSym& cov ) { cov = *fPMeasurementCov; }
53     void GetState( TVectorD& state ) { state = *fPX; }
54     void GetStateCov ( TMatrixDSym& cov ) { cov = *fPXcov; }
55     void GetSeed( TVectorD& seed, TMatrixDSym& seedCov ) { seed = *fPX; seedCov = *fPXcov; }
56     void SetMeasurement( const TVectorD& mes ) {*fPMeasurement = mes;}
57     void SetMeasurementCov( const TMatrixDSym& cov ) {*fPMeasurementCov = cov;}
58     void SetState( const TVectorD& param ) {*fPX = param;}
59     void SetStateCov (const TMatrixDSym& cov ) {*fPXcov = cov;}
60     void SetSeed( const TVectorD& seed, const TMatrixDSym& seedCov ) {*fPX = seed; *fPXcov = seedCov; }
61
62     //Expert methods:
63     Bool_t FindCosmicTrackletNumbersInEvent( Int_t& ITSgood1, Int_t& TPCgood1, Int_t& ITSgood2, Int_t& TPCgood2, const AliESDEvent* pEvent );
64     Bool_t PrepareUpdate();
65     Bool_t Update();
66     void SetRefSurface( const Double_t x, const Double_t alpha );
67     void PrintDebugInfo();
68     void PrintCorrelationMatrix();
69     void PrintCovarianceCorrection();
70     void PrintSystemMatrix();
71     void ResetCovariance();
72     Double_t* GetStateArr() { return fPX->GetMatrixArray(); }
73     Double_t* GetStateCovArr() { return fPXcov->GetMatrixArray(); }
74     Double_t* GetMeasurementArr() { return fPMeasurement->GetMatrixArray(); }
75     Double_t* GetMeasurementCovArr() { return fPMeasurementCov->GetMatrixArray(); }
76     Double_t* GetDeltaArr() {return fDelta;}
77     TH1D* GetMes0Hist() {return fPMes0Hist;}
78     TH1D* GetMes1Hist() {return fPMes1Hist;} 
79     TH1D* GetMes2Hist() {return fPMes2Hist;}
80     TH1D* GetMes3Hist() {return fPMes3Hist;}
81     TH1D* GetMesErr0Hist() {return fPMesErr0Hist;}
82     TH1D* GetMesErr1Hist() {return fPMesErr1Hist;}
83     TH1D* GetMesErr2Hist() {return fPMesErr2Hist;}
84     TH1D* GetMesErr3Hist() {return fPMesErr3Hist;}
85     Bool_t SetCalibrationMode( Bool_t cp=kTRUE );
86     void SetApplyCovarianceCorrection( const Bool_t s=kTRUE ) {fApplyCovarianceCorrection = s;}
87     void SetOutRejSigma( const Double_t a=2. ) { fOutRejSigmas = a; }
88     void SetRejectOutliers( const Bool_t r=kTRUE ) {fRejectOutliers = r;}
89     void SetCovarianceCorrection( const TMatrixDSym& c ) {*fPMeasurementCovCorr = c;}
90     void GetCovarianceCorrection( TMatrixDSym& c ) {c=*fPMeasurementCovCorr;}
91     void SetTrackParams1( const AliExternalTrackParam* exparam );
92     void SetTrackParams2( const AliExternalTrackParam* exparam );
93     void SetQ( const Double_t Q = 1e-10 ) { fQ = Q; }
94     static Bool_t MisalignTrack( AliExternalTrackParam* tr, const TVectorD& misalignment );
95     static void Angles( TVectorD &angles, const TMatrixD &rotmat );
96     static void RotMat( TMatrixD& R, const TVectorD& angles );
97     static void TMatrixDSym_from_TMatrixD( TMatrixDSym& matsym, const TMatrixD& mat );
98     
99 protected:
100     Bool_t UpdateCalibration();
101     Bool_t UpdateEstimateKalman();
102     Bool_t PrepareMeasurement();
103     Bool_t PrepareSystemMatrix();
104     Bool_t PredictMeasurement( TVectorD& z, const TVectorD& x );
105     Bool_t IsOutlier( const TVectorD& update, const TMatrixDSym& covmatrix );
106     Bool_t CalculateCovarianceCorrection();
107
108 private:
109     static const Int_t fgkNMeasurementParams = 4;           //how many measurables
110     static const Int_t fgkNSystemParams = 8;                //how many fit parameters
111     
112     //Track parameters
113     Double_t fAlpha;       //rotation angle between the local and global coordinate system like in AliExternalTrackParam
114     Double_t fLocalX;      //local x coordinate of reference plane = r(global)
115     const AliExternalTrackParam* fPTrackParam1;   //local track parameters (theta,phi,y,z)
116     const AliExternalTrackParam* fPTrackParam2;   //local track parameters
117
118     //Kalman filter related stuff
119     TVectorD* fPX; //System (fit) parameters (phi, theta, psi, x, y, z, driftcorr, driftoffset )
120     TMatrixDSym* fPXcov; //covariance matrix of system parameters
121     TMatrixD* fPH;      //System measurement matrix
122     Double_t fQ;        //measure for system noise
123     TVectorD* fPMeasurement; //the measurement vec for Kalman filter (theta,phi,x,z)
124     TMatrixDSym* fPMeasurementCov; //measurement vec cvariance
125     Double_t fOutRejSigmas; //number of sigmas for outlier rejection
126     Double_t fDelta[fgkNSystemParams]; //array with differentials for calculating derivatives for every parameter(see PrepareSystemMatrix())
127
128     //Control
129     Bool_t fRejectOutliers; //whether to do outlier rejection in the Kalman filter
130     Bool_t fCalibrationMode;            //are we running in calibration mode?
131     Bool_t fFillHistograms;     //whether to fill the histograms with residuals
132     Bool_t fRequireDoubleTPCtrack;
133     Bool_t fApplyCovarianceCorrection;         //add the correction to the covariance of measurement
134     Bool_t fCuts;    //track cuts?
135     Int_t fMinPointsVol1;   //mininum number of points in volume 1
136     Int_t fMinPointsVol2;   //mininum number of points in volume 2
137     Double_t fMinMom;       //min momentum of track for track cuts
138     Double_t fMaxMom;       //max momentum of track for track cuts
139     Double_t fMinAbsSinPhi; //more cuts
140     Double_t fMaxAbsSinPhi; //more cuts
141     Double_t fMinSinTheta;  //cuts
142     Double_t fMaxSinTheta;  //cuts
143     Double_t fMaxMatchingAngle; //cuts
144     Double_t fMaxMatchingDistance; //cuts
145     
146     //Counters
147     Int_t fNTracks; //number of processed tracks
148     Int_t fNUpdates; //number of successful Kalman updates
149     Int_t fNOutliers; //number of outliers
150     Int_t fNMatchedCosmics; //number of cosmic events with matching tracklets
151     Int_t fNProcessedEvents; //number of processed events
152
153     //Calibration histograms
154     TH1D* fPMes0Hist;  //histo of x measurement
155     TH1D* fPMes1Hist;  //histo of y measurement
156     TH1D* fPMes2Hist; //histo of phi measurement
157     TH1D* fPMes3Hist; //histo of theta measurement
158     TH1D* fPMesErr0Hist;  //histogram of the covariance of a fit parameter, used in calibration
159     TH1D* fPMesErr1Hist;  //histogram of the covariance of a fit parameter, used in calibration
160     TH1D* fPMesErr2Hist;  //histogram of the covariance of a fit parameter, used in calibration
161     TH1D* fPMesErr3Hist;  //histogram of the covariance of a fit parameter, used in calibration
162     TMatrixDSym* fPMeasurementCovCorr; //correction to be added to the measurement covariance
163     
164     ClassDef(AliRelAlignerKalman,1)     //AliRelAlignerKalman class
165 };
166
167 #endif
168