PandaRoot
PndHoughTrackFinder.h
Go to the documentation of this file.
1 // PndHoughTrackFinder
3 // Finds Track for one event
5 
16 #ifndef PndHoughTrackFinder_H_
17 #define PndHoughTrackFinder_H_
18 
19 #include <vector>
20 
21 #include "PndSttCA.h"
22 
23 #include "PndHoughData.h"
24 #include "PndHoughUtilities.h"
25 #include "PndHoughTransformation.h"
26 #include "PndHoughSegment.h"
28 #include "PndHoughPreprocessing.h"
29 #include "PndHoughMerge.h"
30 
32  public:
34 
35  PndHoughTrackFinder(TClonesArray *tubeArray) : fNbins1(450), fNbins2(450), fNBinsSeg(90.), fBz(2.), fCutMergeByHoughSpace(5.), fWithGPU(kFALSE), fWithGhostReduction(kTRUE)
36  {
37  // Generate TrackFinderData-Object
38  fData = new PndHoughData(tubeArray);
39  fCATrackFinder = new PndSttCA(tubeArray);
40  ioman = FairRootManager::Instance();
41  };
42 
45  {
46  delete fData;
47  delete fCATrackFinder;
48 
49  delete fPndHoughTransformation;
50  delete fPndHoughPreprocessing;
51  delete fPndHoughMerge;
52  delete fPndHoughSegment;
53  delete fPndHoughUtilities;
54  delete fPndHoughTrackCorrection;
55  }
56 
58  virtual void SetBinningX(double bin) { fNbins1 = bin; };
60  virtual void SetBinningY(double bin) { fNbins2 = bin; };
62  virtual void SetNBinsSeg(Int_t n) { fNBinsSeg = n; };
64  virtual void SetBz(double B) { fBz = B; };
66  virtual void SetCutMergeByHoughSpace(float cut){ fCutMergeByHoughSpace = cut; }
68  virtual void SetWithGPU(bool cuda) { fWithGPU = cuda; }
70  virtual void SetWithGhostReduction(bool ghostred) { fWithGhostReduction = ghostred; }
71  virtual void SetPreselectedTrackCands(std::vector<PndTrackCand> cands) { fPreselectedTrackCands = cands; };
72  virtual void SetCombinedSkewed(TClonesArray *skewed) { fCombinedSkewed = skewed; };
74  void Init()
75  {
76  fData->Init(fNbins1, fNbins2);
77  fCATrackFinder->SetUseGPU(fWithGPU);
78  }
79 
81  void Reset()
82  {
83  fData->clear();
84 
85  fCATrackFinder->SetUseGPU(fWithGPU);
86  fCATrackFinder->Reset();
87 
88  fPreselectedTrackCands.clear();
89  fPreselectedTrackCircles.clear();
90  }
91 
93  void AddHits(TClonesArray *hits, TString branchName);
94 
97  {
100  }
102  void FindTracks();
103 
105  void Preselection();
106 
108  void FindHoughMaxima();
109  void AddSkewedHits(int i);
111  PndTrackCand GetPreselectedTrackCand(int i) { return fPreselectedTrackCands[i]; };
112 
114  int GetNumPreselectedTrackCands() { return fPreselectedTrackCands.size(); };
115 
117  TVector3 GetPreselectedTrackCircles(int i) { return fPreselectedTrackCircles[i]; };
118 
120  int GetNumPreselectedTrackCircles() { return fPreselectedTrackCircles.size(); };
121 
123  PndTrack GetMergedTrack(int i) { return fPndHoughMerge->GetMergedTrack(i); };
124 
126  int GetNumMergedTracks() { return fPndHoughMerge->GetNumMergedTracks(); };
127 
129  // PndTrackCand GetMergedTrackCand(int i) { return fPndHoughMerge->GetMergedTrackCand(i); };
130 
132  // int GetNumMergedTrackCands() { return fPndHoughMerge->GetNumMergedTrackCands(); };
133 
135  std::map<FairLink, FairHit *> GetMapFairLinktoFairHit() { return fData->GetMapFairLinktoFairHit(); };
136 
138  PndHoughData *GetData() { return fData; };
139 
140  private:
141  FairRootManager *ioman = nullptr;
142 
143  PndSttCA *fCATrackFinder = nullptr;
144 
145  PndHoughData *fData = nullptr;
146  PndHoughTransformation *fPndHoughTransformation = nullptr;
147  PndHoughPreprocessing *fPndHoughPreprocessing = nullptr;
148  PndHoughMerge *fPndHoughMerge = nullptr;
149  PndHoughSegment *fPndHoughSegment = nullptr;
150  PndHoughUtilities *fPndHoughUtilities = nullptr;
151  PndHoughTrackCorrection *fPndHoughTrackCorrection = nullptr;
152  TClonesArray *fCombinedSkewed = nullptr;
153  std::vector<std::vector<TVector3>> fHoughSpacesToTracks;
154  std::vector<PndTrackCand> fHoughSpacesToTrackCands;
155  std::vector<PndTrackCand> fHoughSpacesToTrackCandTots;
156 
157  std::vector<FairLink> fFoundHits;
158 
159  std::vector<PndTrackCand> fPreselectedTrackCands;
160  std::vector<TVector3> fPreselectedTrackCircles;
161  std::vector<PndTrack> fApolloniusMergedTracks;
162  std::vector<PndTrackCand> fApolloniusMergedTrackCands;
163 
164  Double_t fBz;
165  float fCutMergeByHoughSpace;
166  Int_t fNbins1;
167  Int_t fNbins2;
168  Int_t fNBinsSeg;
169  bool fWithGPU;
170  bool fWithGhostReduction;
171  double fDistanceThresholdSTTCombinedSkewed = 1.;
172  std::map<FairLink, FairHit *> fMapFairLinktoFairHit;
173  std::map<FairLink, double> fMapFairLinktoIsochrone;
174 
175  ClassDef(PndHoughTrackFinder, 1);
176 };
177 
178 #endif /*PndHoughTrackFinder_H_*/
virtual void SetWithGhostReduction(bool ghostred)
Sets a bool for deciding to use a ghost reduction.
PndTrack GetMergedTrack(int i)
Returns a specific found track.
int GetNumMergedTracks()
Returns the number of found tracks.
void CreateSTTNeighborhoodData()
For all STT hits all other STT neighbors are counted and stored in a map.
int GetNumPreselectedTrackCircles()
Returns the number of preselected circles for track candidates.
int GetNumMergedTracks()
Returns a specific merged track candidate.
Definition: PndHoughMerge.h:70
void Reset()
Resets the data for a new event.
std::map< FairLink, FairHit * > GetMapFairLinktoFairHit() const
Returns the map linking FairLinks to FairHits.
Definition: PndHoughData.h:90
void Preselection()
Here all data are preselected in smaller tracklets.
void Init()
Initializes the data structure and the Hough space of the HoughtrackFinder.
virtual void SetCutMergeByHoughSpace(float cut)
Sets the cut value for the distance of two maxima in the Hough space which have to be merged...
virtual void SetPreselectedTrackCands(std::vector< PndTrackCand > cands)
unsigned int i
Definition: P4_F32vec4.h:21
virtual void SetNBinsSeg(Int_t n)
Sets the number of bins used for the segmentation preselection algorithm.
TVector3 GetPreselectedTrackCircles(int i)
Returns a the corresponding circle for a specific preselected track candidate.
virtual void SetWithGPU(bool cuda)
Sets a bool for deciding to use cuda.
virtual void SetBinningY(double bin)
Sets the number of bins in y direction for the Hough space.
void CreateNeighborhood()
Creates the neighborhood parameters for all investigated hits.
void Init(int NBins1=450, int NBins2=450)
Initializes the Hough space.
Definition: PndHoughData.h:53
int GetNumPreselectedTrackCands()
Returns the number of preselected track candidates.
virtual void SetCombinedSkewed(TClonesArray *skewed)
void clear()
Clears all data maps.
Definition: PndHoughData.h:62
void AddSkewedHits(int i)
std::map< FairLink, FairHit * > GetMapFairLinktoFairHit()
Returns a specific found track candidate.
void FindHoughMaxima()
Performs a Hough transformation for all preselected tracklets. The found track parameters are stored ...
PndTrack GetMergedTrack(int i)
Returns a specific merged track.
Definition: PndHoughMerge.h:66
PndHoughData * GetData()
Returns the data class of the HoughTrackFinder.
void Reset()
Definition: PndSttCA.h:76
PndTrackCand GetPreselectedTrackCand(int i)
Returns a specific preselected track candidate.
PndHoughTrackFinder(TClonesArray *tubeArray)
virtual void SetBz(double B)
Sets the z component of the magnetic field.
void AddHits(TClonesArray *hits, TString branchName)
Adds hits to the data structure of the HoughTrackFinder.
void FindTracks()
Main function of the HoughTrackFinder, which finds the tracks.
virtual void SetBinningX(double bin)
Sets the number of bins in x direction for the Hough space.
void CreateGEMNeighborhoodData()
For all GEM hits all other GEM hits in a certain distance (hier 1.5 cm) are counted and stored as nei...
void SetUseGPU(Bool_t val)
Definition: PndSttCA.h:49