ROOT logo
#ifndef ALIRELALIGNERKALMAN_H
#define ALIRELALIGNERKALMAN_H

/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice                               */

///////////////////////////////////////////////////////////////////////////////
//
//     Relative alignment of two tracking volumes (default ITS and TPC)
//     (see AliRelAlignerKalman.cxx for details)
//
//     Origin: Mikolaj Krzewicki, Nikhef, Mikolaj.Krzewicki@cern.ch
//
//////////////////////////////////////////////////////////////////////////////

#include <TMath.h>
#include <TVectorD.h>
#include <TMatrix.h>

class TObject;
class AliExternalTrackParam;
class AliESDEvent;
class AliESDtrack;
class TArrayI;
class TObjArray;

class AliRelAlignerKalman : public TObject {

public:
    AliRelAlignerKalman();
    virtual ~AliRelAlignerKalman();
    AliRelAlignerKalman& operator= (const AliRelAlignerKalman& a );
    AliRelAlignerKalman(const AliRelAlignerKalman& a);

    //User methods:
    Bool_t AddCosmicEvent( const AliESDEvent* pEvent );
    Bool_t AddTrackParams( const AliExternalTrackParam* p1, const AliExternalTrackParam* p2 );
    
    void Print(Option_t* option="") const;

    Double_t GetPsi() const {return (*fPX)(0);}
    Double_t GetTheta() const {return (*fPX)(1);}
    Double_t GetPhi() const {return (*fPX)(2);}
    Double_t GetX() const {return (*fPX)(3);}
    Double_t GetY() const {return (*fPX)(4);}
    Double_t GetZ() const {return (*fPX)(5);}
    Double_t GetTPCvdCorr() const {return (*fPX)(6);}
    Double_t GetTPCt0() const {return (*fPX)(7);}
    Double_t GetTPCvdY() const {if (fgkNSystemParams>8) return (*fPX)(8); else return 0.0;}
    Double_t GetPsiErr() const {return TMath::Sqrt((*fPXcov)(0,0));}
    Double_t GetThetaErr() const {return TMath::Sqrt((*fPXcov)(1,1));}
    Double_t GetPhiErr() const {return TMath::Sqrt((*fPXcov)(2,2));}
    Double_t GetXErr() const {return TMath::Sqrt((*fPXcov)(3,3));}
    Double_t GetYErr() const {return TMath::Sqrt((*fPXcov)(4,4));}
    Double_t GetZErr() const {return TMath::Sqrt((*fPXcov)(5,5));}
    Double_t GetTPCvdCorrErr() const {return TMath::Sqrt((*fPXcov)(6,6));}
    Double_t GetTPCt0Err() const {return TMath::Sqrt((*fPXcov)(7,7));}
    Double_t GetTPCvdYErr() const {if (fgkNSystemParams>8) return TMath::Sqrt((*fPXcov)(8,8)); else return 0.0;}
    void GetMeasurement( TVectorD& mes ) const { mes = *fPMeasurement; }
    TVectorD* GetMeasurement() { return fPMeasurement; }
    void GetMeasurementCov( TMatrixDSym& cov ) const { cov = *fPMeasurementCov; }
    void SetMeasurement( const TVectorD& mes ) {*fPMeasurement = mes;}
    void SetMeasurementCov( const TMatrixDSym& cov ) {*fPMeasurementCov = cov;}
    TMatrixDSym* GetMeasurementCov() const { return fPMeasurementCov; }
    void GetState( TVectorD& state ) const { state = *fPX; }
    TVectorD* GetState() const { return fPX; }
    void GetStateCov ( TMatrixDSym& cov ) const { cov = *fPXcov; }
    void SetState( const TVectorD& param ) {*fPX = param;}
    void SetStateCov (const TMatrixDSym& cov ) {*fPXcov = cov;}
    TMatrixDSym* GetStateCov() const { return fPXcov; }
    void GetSeed( TVectorD& seed, TMatrixDSym& seedCov ) const { seed = *fPX; seedCov = *fPXcov; }
    void SetSeed( const TVectorD& seed, const TMatrixDSym& seedCov ) {*fPX = seed; *fPXcov = seedCov; }
    Bool_t Merge( const AliRelAlignerKalman* al );
    Long64_t Merge( TCollection* list );

    //Expert methods:
    Bool_t AddESDevent( const AliESDEvent* pEvent );
    Bool_t AddESDtrack( const AliESDtrack* pTrack );
    void SetMagField( const Double_t f ) { fMagField=f; }
    Double_t GetMagField() const { return fMagField; }
    Bool_t FindCosmicTrackletNumbersInEvent( TArrayI& outITStracksTArr, TArrayI& outTPCtracksTArr, const AliESDEvent* pEvent );
    Int_t FindMatchingTracks(TObjArray& arrITS, TObjArray& arrTPC, AliESDEvent* pESD);
    Bool_t Update();
    void PrintCorrelationMatrix();
    //void PrintCovarianceCorrection();
    void PrintSystemMatrix();
  Int_t CheckCovariance(); // check covariance matrix
    void Reset();
    void ResetCovariance( const Double_t number=0. );
    void ResetTPCparamsCovariance( const Double_t number=0. );
    Double_t* GetStateArr() const { return fPX->GetMatrixArray(); }
    Double_t* GetStateCovArr() const { return fPXcov->GetMatrixArray(); }
    Double_t* GetMeasurementArr() const { return fPMeasurement->GetMatrixArray(); }
    Double_t* GetMeasurementCovArr() const { return fPMeasurementCov->GetMatrixArray(); }
    TMatrixD* GetH() const { return fPH; }
    TVectorD* GetMeasurementPrediction() const {return fPMeasurementPrediction;}
    const Double_t* GetDeltaArr() const {return fDelta;}
    void SetNumericalParanoia( const Bool_t mode=kFALSE ) { fNumericalParanoia=mode; }
    void SetCorrectionMode( const Bool_t mode=kTRUE ) { fCorrectionMode=mode; }
    void SetOutRejSigma( const Double_t a=2. ) { fOutRejSigmas = a; }
    void SetRejectOutliers( const Bool_t r=kTRUE ) {fRejectOutliers = r;}
    void SetRejectOutliersSigma2Median( const Bool_t b=kTRUE );
    void SetOutRejSigma2Median( const Double_t s ) {fOutRejSigma2Median = s;}
    Bool_t SetTrackParams( const AliExternalTrackParam* exparam1, const AliExternalTrackParam* exparam2 );
    const AliExternalTrackParam* GetTrackParams1() const {return fPTrackParam1;}
    const AliExternalTrackParam* GetTrackParams2() const {return fPTrackParam2;}
    void SetMinPointsVol1( const Int_t min ) {fMinPointsVol1=min;}
    void SetMinPointsVol2( const Int_t min ) {fMinPointsVol2=min;}
    void SetRequireMatchInTPC( const Bool_t s=kTRUE ) {fRequireMatchInTPC = s;}
    void SetQ( const Double_t Q = 1e-10 ) { fQ = Q; }
    void SetMaxMatchingDistance( const Double_t m ) {fMaxMatchingDistance=m;}
    void SetMaxMatchingAngle( const Double_t m ) {fMaxMatchingAngle=m;}
    void SetTPCvd( const Float_t v ) {fTPCvd=v;}
    void SetTPCZLengthA( const Double_t l ) {fTPCZLengthA=l;}
    void SetTPCZLengthC( const Double_t l ) {fTPCZLengthC=l;}
    Bool_t CorrectTrack( AliExternalTrackParam* tr, const TVectorD& misalignment ) const;
    Bool_t MisalignTrack( AliExternalTrackParam* tr, const TVectorD& misalignment ) const;
    static void Angles( TVectorD &angles, const TMatrixD &rotmat );
    static void RotMat( TMatrixD& R, const TVectorD& angles );
    static void TMatrixDSymFromTMatrixD( TMatrixDSym& matsym, const TMatrixD& mat );
    Bool_t IsPositiveDefinite( const TMatrixD& mat ) const;
    void SetTimeStamp( const UInt_t ts ) { fTimeStamp = ts; }
    UInt_t GetTimeStamp() const {return fTimeStamp;}
    void SetRunNumber( const Int_t rn ) { fRunNumber = rn; }
    Int_t GetRunNumber() const {return fRunNumber;}
    Int_t Compare(const TObject *obj) const;
    Bool_t IsSortable() const { return kTRUE; }
    Int_t GetNTracks() const {return fNTracks;}
    Int_t GetNUpdates() const {return fNUpdates;}
    Int_t GetNOutliers() const {return fNOutliers;}
    Int_t GetNOutliersSigma2Median() const {return fNOutliersSigma2Median;}
    Int_t GetNMerges() const {return fNMerges;}
    Int_t GetNMergesFailed() const {return fNMergesFailed;}
    void SetPoint2Track( Bool_t o );
    
protected:
    Bool_t UpdateEstimateKalman();
    Bool_t PrepareMeasurement();
    Bool_t PrepareSystemMatrix();
    Bool_t PreparePrediction();
    Bool_t PredictMeasurement( TVectorD& z, const TVectorD& x );
    Bool_t IsOutlier( const TVectorD& update, const TMatrixDSym& covmatrix );
    Bool_t IsOutlierSigma2Median( const AliExternalTrackParam* pITS, const AliExternalTrackParam* pTPC );

private:
    static const Int_t fgkNSystemParams=9;                //how many fit parameters
    static const Int_t fgkNtracksSigma2Median=500;        //how many sets for median and rms
    
    //Track parameters
    AliExternalTrackParam* fPTrackParam1;   //!local track parameters
    AliExternalTrackParam* fPTrackParam2;   //!local track parameters
    Double_t fMagField; //magnetic field

    //Kalman filter related stuff
    Int_t fNMeasurementParams;           //how many measurables
    TVectorD* fPX; //System (fit) parameters (phi, theta, psi, x, y, z, driftcorr, driftoffset )
    TMatrixDSym* fPXcov; //covariance matrix of system parameters
    TMatrixD* fPH;      //!System measurement matrix
    Double_t fQ;        //!measure for system noise
    TVectorD* fPMeasurement; //!the measurement vec for Kalman filter (theta,phi,x,z)
    TMatrixDSym* fPMeasurementCov; //!measurement vec cvariance
    TVectorD* fPMeasurementPrediction; //!prediction of the measurement
    Double_t fOutRejSigmas; //number of sigmas for outlier rejection
    Double_t fOutRejSigma2Median; //nsigmas to median of input residual distribution
    Double_t fDelta[fgkNSystemParams]; //array with differentials for calculating derivatives for every parameter(see PrepareSystemMatrix())
    Double_t* fResArrSigma2Median[4]; //!holds residuals for median based outlier removal

    //Control
    Bool_t fYZOnly;   //whether to consider only yz without directions.
    Bool_t fNumericalParanoia; //whether to perform additional checks for numerical stability
    Bool_t fRejectOutliers; //whether to do outlier rejection in the Kalman filter
    Bool_t fRejectOutliersSigma2Median; //whether to reject input based on distance to median
    Bool_t fRequireMatchInTPC;  //when looking for a cosmic in event, require that TPC has 2 matching segments
    Bool_t fCuts;    //track cuts?
    Int_t fMinPointsVol1;   //mininum number of points in volume 1
    Int_t fMinPointsVol2;   //mininum number of points in volume 2
    Double_t fMinPt;       //min momentum of track for track cuts
    Double_t fMaxPt;       //max momentum of track for track cuts
    Double_t fMaxMatchingAngle; //cuts
    Double_t fMaxMatchingDistance; //cuts
    Bool_t fCorrectionMode; //calculate corrective transform for TPC (or monitor actual TPC misal params)
    
    //Counters
    Int_t fNTracks; //number of processed tracks
    Int_t fNUpdates; //number of successful Kalman updates
    Int_t fNOutliers; //number of outliers
    Int_t fNOutliersSigma2Median; //number of rejected inputs
    Int_t fNMatchedCosmics; //number of cosmic events with matching tracklets (good cosmics)
    Int_t fNMatchedTPCtracklets;//number of cosmic events with 2 matching TPC tracklets
    Int_t fNProcessedEvents; //number of processed events
    UInt_t fTimeStamp;    //time stamp
    Int_t fRunNumber;    //run number
    Int_t fNMerges;      //how many succesful merges
    Int_t fNMergesFailed; //how many merges failed

    //TPC stuff
    Double_t fTPCvd; //TPC drift velocity
    Double_t fTPCZLengthA; //TPC length side A
    Double_t fTPCZLengthC; //TPC length side C
    
    ClassDef(AliRelAlignerKalman,4)     //AliRelAlignerKalman class
};

#endif


 AliRelAlignerKalman.h:1
 AliRelAlignerKalman.h:2
 AliRelAlignerKalman.h:3
 AliRelAlignerKalman.h:4
 AliRelAlignerKalman.h:5
 AliRelAlignerKalman.h:6
 AliRelAlignerKalman.h:7
 AliRelAlignerKalman.h:8
 AliRelAlignerKalman.h:9
 AliRelAlignerKalman.h:10
 AliRelAlignerKalman.h:11
 AliRelAlignerKalman.h:12
 AliRelAlignerKalman.h:13
 AliRelAlignerKalman.h:14
 AliRelAlignerKalman.h:15
 AliRelAlignerKalman.h:16
 AliRelAlignerKalman.h:17
 AliRelAlignerKalman.h:18
 AliRelAlignerKalman.h:19
 AliRelAlignerKalman.h:20
 AliRelAlignerKalman.h:21
 AliRelAlignerKalman.h:22
 AliRelAlignerKalman.h:23
 AliRelAlignerKalman.h:24
 AliRelAlignerKalman.h:25
 AliRelAlignerKalman.h:26
 AliRelAlignerKalman.h:27
 AliRelAlignerKalman.h:28
 AliRelAlignerKalman.h:29
 AliRelAlignerKalman.h:30
 AliRelAlignerKalman.h:31
 AliRelAlignerKalman.h:32
 AliRelAlignerKalman.h:33
 AliRelAlignerKalman.h:34
 AliRelAlignerKalman.h:35
 AliRelAlignerKalman.h:36
 AliRelAlignerKalman.h:37
 AliRelAlignerKalman.h:38
 AliRelAlignerKalman.h:39
 AliRelAlignerKalman.h:40
 AliRelAlignerKalman.h:41
 AliRelAlignerKalman.h:42
 AliRelAlignerKalman.h:43
 AliRelAlignerKalman.h:44
 AliRelAlignerKalman.h:45
 AliRelAlignerKalman.h:46
 AliRelAlignerKalman.h:47
 AliRelAlignerKalman.h:48
 AliRelAlignerKalman.h:49
 AliRelAlignerKalman.h:50
 AliRelAlignerKalman.h:51
 AliRelAlignerKalman.h:52
 AliRelAlignerKalman.h:53
 AliRelAlignerKalman.h:54
 AliRelAlignerKalman.h:55
 AliRelAlignerKalman.h:56
 AliRelAlignerKalman.h:57
 AliRelAlignerKalman.h:58
 AliRelAlignerKalman.h:59
 AliRelAlignerKalman.h:60
 AliRelAlignerKalman.h:61
 AliRelAlignerKalman.h:62
 AliRelAlignerKalman.h:63
 AliRelAlignerKalman.h:64
 AliRelAlignerKalman.h:65
 AliRelAlignerKalman.h:66
 AliRelAlignerKalman.h:67
 AliRelAlignerKalman.h:68
 AliRelAlignerKalman.h:69
 AliRelAlignerKalman.h:70
 AliRelAlignerKalman.h:71
 AliRelAlignerKalman.h:72
 AliRelAlignerKalman.h:73
 AliRelAlignerKalman.h:74
 AliRelAlignerKalman.h:75
 AliRelAlignerKalman.h:76
 AliRelAlignerKalman.h:77
 AliRelAlignerKalman.h:78
 AliRelAlignerKalman.h:79
 AliRelAlignerKalman.h:80
 AliRelAlignerKalman.h:81
 AliRelAlignerKalman.h:82
 AliRelAlignerKalman.h:83
 AliRelAlignerKalman.h:84
 AliRelAlignerKalman.h:85
 AliRelAlignerKalman.h:86
 AliRelAlignerKalman.h:87
 AliRelAlignerKalman.h:88
 AliRelAlignerKalman.h:89
 AliRelAlignerKalman.h:90
 AliRelAlignerKalman.h:91
 AliRelAlignerKalman.h:92
 AliRelAlignerKalman.h:93
 AliRelAlignerKalman.h:94
 AliRelAlignerKalman.h:95
 AliRelAlignerKalman.h:96
 AliRelAlignerKalman.h:97
 AliRelAlignerKalman.h:98
 AliRelAlignerKalman.h:99
 AliRelAlignerKalman.h:100
 AliRelAlignerKalman.h:101
 AliRelAlignerKalman.h:102
 AliRelAlignerKalman.h:103
 AliRelAlignerKalman.h:104
 AliRelAlignerKalman.h:105
 AliRelAlignerKalman.h:106
 AliRelAlignerKalman.h:107
 AliRelAlignerKalman.h:108
 AliRelAlignerKalman.h:109
 AliRelAlignerKalman.h:110
 AliRelAlignerKalman.h:111
 AliRelAlignerKalman.h:112
 AliRelAlignerKalman.h:113
 AliRelAlignerKalman.h:114
 AliRelAlignerKalman.h:115
 AliRelAlignerKalman.h:116
 AliRelAlignerKalman.h:117
 AliRelAlignerKalman.h:118
 AliRelAlignerKalman.h:119
 AliRelAlignerKalman.h:120
 AliRelAlignerKalman.h:121
 AliRelAlignerKalman.h:122
 AliRelAlignerKalman.h:123
 AliRelAlignerKalman.h:124
 AliRelAlignerKalman.h:125
 AliRelAlignerKalman.h:126
 AliRelAlignerKalman.h:127
 AliRelAlignerKalman.h:128
 AliRelAlignerKalman.h:129
 AliRelAlignerKalman.h:130
 AliRelAlignerKalman.h:131
 AliRelAlignerKalman.h:132
 AliRelAlignerKalman.h:133
 AliRelAlignerKalman.h:134
 AliRelAlignerKalman.h:135
 AliRelAlignerKalman.h:136
 AliRelAlignerKalman.h:137
 AliRelAlignerKalman.h:138
 AliRelAlignerKalman.h:139
 AliRelAlignerKalman.h:140
 AliRelAlignerKalman.h:141
 AliRelAlignerKalman.h:142
 AliRelAlignerKalman.h:143
 AliRelAlignerKalman.h:144
 AliRelAlignerKalman.h:145
 AliRelAlignerKalman.h:146
 AliRelAlignerKalman.h:147
 AliRelAlignerKalman.h:148
 AliRelAlignerKalman.h:149
 AliRelAlignerKalman.h:150
 AliRelAlignerKalman.h:151
 AliRelAlignerKalman.h:152
 AliRelAlignerKalman.h:153
 AliRelAlignerKalman.h:154
 AliRelAlignerKalman.h:155
 AliRelAlignerKalman.h:156
 AliRelAlignerKalman.h:157
 AliRelAlignerKalman.h:158
 AliRelAlignerKalman.h:159
 AliRelAlignerKalman.h:160
 AliRelAlignerKalman.h:161
 AliRelAlignerKalman.h:162
 AliRelAlignerKalman.h:163
 AliRelAlignerKalman.h:164
 AliRelAlignerKalman.h:165
 AliRelAlignerKalman.h:166
 AliRelAlignerKalman.h:167
 AliRelAlignerKalman.h:168
 AliRelAlignerKalman.h:169
 AliRelAlignerKalman.h:170
 AliRelAlignerKalman.h:171
 AliRelAlignerKalman.h:172
 AliRelAlignerKalman.h:173
 AliRelAlignerKalman.h:174
 AliRelAlignerKalman.h:175
 AliRelAlignerKalman.h:176
 AliRelAlignerKalman.h:177
 AliRelAlignerKalman.h:178
 AliRelAlignerKalman.h:179
 AliRelAlignerKalman.h:180
 AliRelAlignerKalman.h:181
 AliRelAlignerKalman.h:182
 AliRelAlignerKalman.h:183
 AliRelAlignerKalman.h:184
 AliRelAlignerKalman.h:185
 AliRelAlignerKalman.h:186
 AliRelAlignerKalman.h:187
 AliRelAlignerKalman.h:188
 AliRelAlignerKalman.h:189
 AliRelAlignerKalman.h:190
 AliRelAlignerKalman.h:191
 AliRelAlignerKalman.h:192
 AliRelAlignerKalman.h:193
 AliRelAlignerKalman.h:194
 AliRelAlignerKalman.h:195
 AliRelAlignerKalman.h:196
 AliRelAlignerKalman.h:197
 AliRelAlignerKalman.h:198
 AliRelAlignerKalman.h:199
 AliRelAlignerKalman.h:200
 AliRelAlignerKalman.h:201
 AliRelAlignerKalman.h:202
 AliRelAlignerKalman.h:203
 AliRelAlignerKalman.h:204
 AliRelAlignerKalman.h:205
 AliRelAlignerKalman.h:206