]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/MUON/OnlineAnalysis/AliHLTMUONFullTracker.h
Reversing commit 42768 because it apparently breaks the code. Waiting for proper...
[u/mrichter/AliRoot.git] / HLT / MUON / OnlineAnalysis / AliHLTMUONFullTracker.h
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 ///
8 ///  @file   AliHLTMUONFullTracker.h
9 ///  @author Indranil Das <indra.das@saha.ac.in> | <indra.ehep@gmail.com>
10 ///  @date   09 Feb 2010
11 ///  @brief  For full tracking in the dimuon HLT.
12 ///
13
14
15 /**********************************************************************
16  Created on : 08/12/2009
17  Purpose    : First version implementation of the Full tracker for dHLT.
18  Author     : Indranil Das, HEP Division, SINP
19  Email      : indra.das@saha.ac.in | indra.ehep@gmail.com
20 **********************************************************************/
21
22 #include <iostream>
23 #include <map>
24
25 #include "TMatrixD.h"
26 #include "TMath.h"
27 #include "AliHLTMUONDataBlockReader.h"
28 #include "AliHLTMUONDataBlockWriter.h"
29
30
31 class   AliHLTMUONUtils;
32 class   AliHLTMUONDataTypes;
33 class   AliHLTMUONConstants;
34 class   AliRunInfo;
35 class   AliLog;
36 class   AliCDBEntry;
37 class   AliMpDEIterator;
38 class   AliMpCDB;
39 class   AliMpSegmentation;
40 class   AliMpDDLStore;
41 class   AliMUONTrackParam;
42 class   AliMUONGeometryTransformer;
43 class   AliMUONConstants;
44 class   AliMUONTrackExtrap;
45 class   AliMUONTrackParam;
46 class   AliMUONTrackExtrap;
47 class   AliGRPObject;
48 class   AliGeomManager;
49 class   AliCDBStorage;
50 class   AliCDBManager;
51 class   AliMagF;
52 class   TGeoGlobalMagField;
53 class   TVector3;
54 class   TString;
55 class   TMap;
56
57 typedef std::map<Int_t, Int_t> DetElemList;
58
59 class AliHLTMUONFullTracker : public AliHLTLogging
60 {
61
62  public :
63   ///Constructor
64   AliHLTMUONFullTracker() ;
65   ///Destructor
66   virtual ~AliHLTMUONFullTracker();
67
68   ///Print message
69   void Print();
70   ///Set the input of trigrec blocks
71   Bool_t SetInput(AliHLTInt32_t ddl, const AliHLTMUONTriggerRecordStruct  *data, AliHLTInt32_t size);
72   ///Set the input of rechit blocks
73   Bool_t SetInput(AliHLTInt32_t ddl, const AliHLTMUONRecHitStruct  *data, AliHLTInt32_t size);
74   ///Main method to run and compute the tracks
75   Bool_t Run(AliHLTEventID_t iEvent,AliHLTMUONTrackStruct *data, AliHLTUInt32_t& size);
76   ///To be called once from DoInit method of component
77   Bool_t Init();
78   ///Max number of points per chamber
79   int MaxNofPointsPerCh(){return fgkMaxNofPointsPerCh;}
80   ///Set for fast tracking bypass  default Kalman tracking
81   void FastTracking(Bool_t isFast){fFastTracking = isFast;}
82   ///Getter for fast tracking
83   Bool_t FastTracking(){return fFastTracking;}
84
85  protected:
86
87   /// copy constructor
88   AliHLTMUONFullTracker(const AliHLTMUONFullTracker& rhs); 
89   /// assignment operator
90   AliHLTMUONFullTracker& operator=(const AliHLTMUONFullTracker& rhs); 
91
92  private :
93
94   /// intger pair needed for QuadTrackSeg method
95   struct IntPair{
96     Int_t fFirst,fSecond; /// First and second
97   };
98   ///Structure for internal track segments
99   struct TrackSeg{
100     Int_t fIndex[4]; /// index array for cluster address
101     AliHLTInt32_t fTrigRec; /// trigrec
102   };
103
104   ///Sructure for clusters
105   struct Cluster{
106     Float_t fX,fY,fZ;  /// position
107     Float_t fErrX2,fErrY2;  /// error in position
108   };
109  
110   ///Sructure for clusters
111   struct HalfTrack{
112     Float_t fPx,fPy,fPz; /// momentum
113     Int_t fCharge;       /// charge
114   }; 
115
116   static const Float_t fgkTrackDetCoordinate[3]; /// set the constant value for third station position and size
117   
118   static const Double_t fgkAbsoedge[4] ;     /// edge of the absorber
119   static const Double_t fgkRadLen[3] ;       /// radiation length of the main three matirials of the front absorber
120   static const Double_t fgkRho[3] ;          /// density of the main three matirials of the front absorber
121   static const Double_t fgkAtomicZ[3] ;      /// atomic number the main three matirials of the front absorber
122   static const Double_t fgkAtomicA[3] ;      /// atomic mass of the main three matirials of the front absorber
123
124   static const Int_t fgkMaxNofCellsPerCh ;      /// maximum number of cell are allowed to create
125   static const Int_t fgkMaxNofPointsPerCh ;     /// maximim number of points per chamber
126   static const Int_t fgkMaxNofCh ;              /// maximum number of chambrs
127   static const Int_t fgkMaxNofTracks;           /// maximum number of allowed tracks
128   static const Int_t fgkMaxNofConnectedTracks;  /// maximum number of back to front connected tracks
129   static const Int_t fgkMaxNofTriggers;         /// maximum number of triggers (condition comes from simulation prediction)
130   
131   AliMUONGeometryTransformer *fChamberGeometryTransformer; /// Pointer to AliMUONGeometryTransformer
132   AliHLTMUONRecHitStruct ***fChPoint; /// array of pointer to rechit data
133   AliHLTMUONTriggerRecordStruct **fChPoint11; ///array of pointer to trigrec data
134   TrackSeg *fBackTrackSeg; /// track segments at the rear part of the spectrometer
135   TrackSeg *fFrontTrackSeg; /// track segments close the part of interaction point  of ALICE
136   Float_t *fExtrapSt3X ; /// Extrapolated x position in third station
137   Float_t *fExtrapSt3Y ; /// Extrapolated y position in third station
138   Float_t *fInclinationBack; /// values of inclination angle of back track segments
139
140   Int_t *fNofConnectedfrontTrackSeg ; /// nof connected tracks in front direction for each back track segments
141   Int_t **fBackToFront; /// Pointer to back to front segment mapping
142   Int_t *fNofPoints ; /// Number of points for each stations
143   AliMUONTrackParam *fTrackParam ; /// track parameters;
144   HalfTrack *fHalfTrack; /// momentum parameters for the tracks which doesnot have tracksegment in quadrants
145
146   Int_t fTotNofPoints; /// Total number of points received from all rechit source
147   Int_t fTotTrackSeg; /// Total number of track segments
148   Int_t fNofCells[2]; // nof cell count per station
149   Int_t fNofbackTrackSeg; /// number of back track segments
150   Int_t fNoffrontTrackSeg; /// number of front track segments
151   Int_t fNofConnected ; /// number of connected track segments
152   AliHLTUInt32_t fNofTracks; /// number of connected track segments
153   DetElemList fDetElemList; ///Map for valid detelem
154   Bool_t fFastTracking ; ///flag for fast tracking avoiding kalman
155   Int_t   fNofInputs; /// Nof inputs
156   Int_t   fNofTriggerInputs; /// Nof inputs
157   Int_t   fNofTrackerInputs; /// Nof inputs
158   Bool_t  fIsMagfield ; /// checks the status of magfield
159
160   ///  Cross Check the inputs
161   Bool_t CheckInput(AliHLTEventID_t iEvent);
162   /// Slat Track segments 
163   Bool_t SlatTrackSeg();
164   /// Calculate preliminary momentum
165   Bool_t PrelimMomCalc();
166   /// Quad Track segments 
167   Bool_t QuadTrackSeg();
168   /// Kalman Chi2 test
169   Bool_t KalmanChi2Test();
170   /// track extrapolation through  dipole magnet to connect front and back track seg
171   Bool_t SelectFront();
172   /// Propagate tracks
173   void PropagateTracks(Double_t charge, Float_t& px, Float_t& py, Float_t& pz, 
174                        Float_t& xr, Float_t& yr, Float_t& zr, Float_t zprop);
175   /// extrapolate to origin
176   Bool_t ExtrapolateToOrigin();
177   /// Clean after each run
178   Bool_t Clear();
179
180
181   /// Angle calculate
182   inline Double_t Angle(const AliHLTMUONRecHitStruct *v1, const AliHLTMUONRecHitStruct *v2);
183   /// Subtracktion of two point
184   inline void Sub(const AliHLTMUONRecHitStruct *v1, const AliHLTMUONRecHitStruct *v2, AliHLTMUONRecHitStruct *v3) const;
185   /// Kalman Filter
186   inline Double_t KalmanFilter(AliMUONTrackParam &trackParamAtCluster, Cluster *cluster);
187   /// Try onecluster
188   inline Double_t TryOneCluster(const AliMUONTrackParam &trackParam, Cluster* cluster,
189                                 AliMUONTrackParam &trackParamAtCluster, Bool_t updatePropagator);
190   inline Bool_t TryOneClusterFast(const AliMUONTrackParam &trackParam, const Cluster* cluster);
191
192   /// MCS effect correction
193   inline void CorrectMCSEffectInAbsorber(AliMUONTrackParam* param,
194                                          Double_t xVtx, Double_t yVtx, Double_t zVtx,
195                                          Double_t absZBeg, 
196                                          Double_t f1, Double_t f2);
197   /// Covariant handling function
198   inline void Cov2CovP(const TMatrixD &param, TMatrixD &cov);
199   /// Covariant handling function
200   inline void CovP2Cov(const TMatrixD &param, TMatrixD &covP);
201   /// Energy loss coreection in front absorber
202   inline void CorrectELossEffectInAbsorber(AliMUONTrackParam* param, Double_t eLoss);
203   /// Linear Extrapolation to Z position
204   inline void LinearExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd);
205   /// Energy loss
206   inline Double_t EnergyLossFluctuation2(Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicA, Double_t atomicZ);
207   /// Bethe Bloch formula of enrgy loss
208   inline Double_t BetheBloch(Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicA, Double_t atomicZ);
209   
210   /// Runge Kutta method of track extrapolation through mag field
211   inline void OneStepRungekutta(Double_t charge, Double_t step, const Double_t* vect, Double_t* vout);
212   /// Helix3 method of track extrapolation through mag field
213   inline void OneStepHelix3(Double_t field, Double_t step, const Double_t *vect, Double_t *vout) const;                           
214   /// Fill the tracks to output pointer
215   Bool_t FillOutData(AliHLTMUONTrackStruct *data, AliHLTUInt32_t& size);
216   
217 };
218 #endif // ALIHLTMUONMANSOTRACKERFSM_H