1 #ifndef ALIHLTMUONFULLTRACKER_H
2 #define ALIHLTMUONFULLTRACKER_H
3 /* This file is property of and copyright by the ALICE HLT Project *
4 * ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
6 /**********************************************************************
7 // Created on : 08/12/2009
8 // Purpose : First version implementation of the Full tracker for dHLT.
9 // Author : Indranil Das, HEP Division, SINP
10 // Email : indra.das@saha.ac.in | indra.ehep@gmail.com
11 **********************************************************************/
13 #include "AliHLTLogging.h"
16 class AliHLTMUONConstants;
17 class AliHLTMUONMansoTrackStruct;
18 class AliHLTMUONTriggerRecordStruct;
19 class AliHLTMUONRecHitStruct;
20 class AliMUONGeometryTransformer;
21 class AliMUONTrackParam;
23 class AliHLTMUONFullTracker : public AliHLTLogging
28 AliHLTMUONFullTracker() ;
30 ~AliHLTMUONFullTracker();
34 ///Set the input of trigrec blocks
35 Bool_t SetInput(AliHLTInt32_t ddl, const AliHLTMUONTriggerRecordStruct *data, AliHLTInt32_t size);
36 ///Set the input of rechit blocks
37 Bool_t SetInput(AliHLTInt32_t ddl, const AliHLTMUONRecHitStruct *data, AliHLTInt32_t size);
38 ///Main method to run and compute the tracks
39 Bool_t Run(int iEvent,AliHLTMUONMansoTrackStruct *data, AliHLTUInt32_t& size);
40 ///To be called once from DoInit method of component
46 AliHLTMUONFullTracker(const AliHLTMUONFullTracker& rhs);
47 /// assignment operator
48 AliHLTMUONFullTracker& operator=(const AliHLTMUONFullTracker& rhs);
52 /// intger pair needed for QuadTrackSeg method
54 Int_t fFirst,fSecond; /// First and second element of the pair.
56 ///Structure for internal track segments
58 Int_t fIndex[4]; /// Index of the chamber rec hit point in fChPoint.
59 AliHLTInt32_t fTrigRec; /// The corresponding trigger record ID.
62 ///Sructure for clusters
64 Float_t fX,fY,fZ; /// Cluster point X, Y, Z.
65 Float_t fErrX2,fErrY2; /// The square of the uncertainty (variance) on X and Y.
69 static const Float_t fgkTrackDetCoordinate[3]; /// set the constant value for third station position and size
71 static const Double_t fgkAbsoedge[4] ; /// edge of the absorber
72 static const Double_t fgkRadLen[3] ; /// radiation length of the main three matirials of the front absorber
73 static const Double_t fgkRho[3] ; /// density of the main three matirials of the front absorber
74 static const Double_t fgkAtomicZ[3] ; /// atomic number the main three matirials of the front absorber
75 static const Double_t fgkAtomicA[3] ; /// atomic mass of the main three matirials of the front absorber
77 static const Int_t fgkMaxNofCellsPerCh ; /// maximum number of cell are allowed to create
78 static const Int_t fgkMaxNofPointsPerCh ; /// maximim number of points per chamber
79 static const Int_t fgkMaxNofCh ; /// maximum number of chambrs
80 static const Int_t fgkMaxNofTracks; /// maximum number of allowed tracks
81 static const Int_t fgkMaxNofConnectedTracks; /// maximum number of back to front connected tracks
83 AliMUONGeometryTransformer *fChamberGeometryTransformer; /// Pointer to AliMUONGeometryTransformer
85 AliHLTMUONRecHitStruct ***fChPoint; /// array of pointer to rechit data
86 AliHLTMUONTriggerRecordStruct **fChPoint11; ///array of pointer to trigrec data
87 TrackSeg *fBackTrackSeg; /// track segments at the rear part of the spectrometer
88 TrackSeg *fFrontTrackSeg; /// track segments close the part of interaction point of ALICE
90 Float_t *fExtrapSt3X ; /// Extrapolated x position in third station
91 Float_t *fExtrapSt3Y ; /// Extrapolated y position in third station
92 Float_t *fInclinationBack; /// values of inclination angle of back track segments
94 Int_t *fNofConnectedfrontTrackSeg ; /// nof connected tracks in front direction for each back track segments
95 Int_t **fBackToFront; /// Pointer to back to front segment mapping
96 Float_t *fCharge; /// Charge of the tracks
97 Int_t *fNofPoints ; /// Number of points for each stations
98 AliMUONTrackParam *fTrackParam ; /// track parameters;
100 Int_t fTotNofPoints; /// Total number of points received from all rechit source
101 Int_t fTotTrackSeg; /// Total number of track segments
102 Int_t fNofCells[2]; /// Number of cells per station in QuadSeg
103 Bool_t fOverflowed; /// Check if overflowed
104 Int_t fNofbackTrackSeg; /// number of back track segments
105 Int_t fNoffrontTrackSeg; /// number of front track segments
106 Int_t fNofConnected ; /// number of connected track segments
108 /// Slat Track segments
109 Bool_t SlatTrackSeg();
110 /// Quad Track segments
111 Bool_t QuadTrackSeg();
113 Bool_t KalmanChi2Test();
114 /// track extrapolation through dipole magnet to connect front and back track seg
115 Bool_t SelectFront();
117 void PropagateTracks(Double_t charge, Float_t& px, Float_t& py, Float_t& pz,
118 Float_t& xr, Float_t& yr, Float_t& zr, Float_t zprop) const;
119 /// extrapolate to origin
120 Bool_t ExtrapolateToOrigin(Bool_t extrap);
121 /// Clean after each run
126 Double_t Angle(const AliHLTMUONRecHitStruct *v1, const AliHLTMUONRecHitStruct *v2) const;
127 /// Subtracktion of two point
128 void Sub(const AliHLTMUONRecHitStruct *v1, const AliHLTMUONRecHitStruct *v2, AliHLTMUONRecHitStruct *v3) const;
130 Double_t KalmanFilter(AliMUONTrackParam &trackParamAtCluster, Cluster *cluster);
132 Double_t TryOneCluster(const AliMUONTrackParam &trackParam, Cluster* cluster,
133 AliMUONTrackParam &trackParamAtCluster, Bool_t updatePropagator);
134 Bool_t TryOneClusterFast(const AliMUONTrackParam &trackParam, const Cluster* cluster) const;
136 /// MCS effect correction
137 void CorrectMCSEffectInAbsorber(AliMUONTrackParam* param,
138 Double_t xVtx, Double_t yVtx, Double_t zVtx,
140 Double_t f1, Double_t f2);
141 /// Covariant handling function
142 void Cov2CovP(const TMatrixD ¶m, TMatrixD &cov);
143 /// Covariant handling function
144 void CovP2Cov(const TMatrixD ¶m, TMatrixD &covP);
145 /// Energy loss coreection in front absorber
146 void CorrectELossEffectInAbsorber(AliMUONTrackParam* param, Double_t eLoss);
147 /// Linear Extrapolation to Z position
148 void LinearExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd) const;
150 Double_t EnergyLossFluctuation2(Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicA, Double_t atomicZ);
151 /// Bethe Bloch formula of enrgy loss
152 Double_t BetheBloch(Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicA, Double_t atomicZ);
154 /// Runge Kutta method of track extrapolation through mag field
155 void OneStepRungekutta(Double_t charge, Double_t step, const Double_t* vect, Double_t* vout) const;
156 /// Helix3 method of track extrapolation through mag field
157 void OneStepHelix3(Double_t field, Double_t step, const Double_t *vect, Double_t *vout) const;
158 /// Initialise GRP when running without reconstruction chain
160 /// Fill the tracks to output pointer
161 Bool_t FillOutData(AliHLTMUONMansoTrackStruct *data, AliHLTUInt32_t& size) const;
164 #endif // ALIHLTMUONMANSOTRACKERFSM_H