FairRoot/PandaRoot
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | List of all members
PndFtsHoughTrackFinder Class Reference

Implementation of the Hough transform based FTS PR. Creates Hough spaces, finds peaks (=tracklets) and combines them to track candidates. More...

#include <PndFtsHoughTrackFinder.h>

Inheritance diagram for PndFtsHoughTrackFinder:
PndFtsHoughTrackFinderQA

Public Member Functions

 PndFtsHoughTrackFinder (PndFtsHoughTrackerTask *trackerTask)
 Set pointer to tracker task (super important as it provides an I/O interface to PandaRoot) More...
 
virtual ~PndFtsHoughTrackFinder ()
 Destructor. More...
 
void OverwriteTrackFinderParams (PndFtsHoughTrackFinderParams newParams)
 
virtual void FindTracks ()
 Performs the track finding. More...
 
Int_t NTracks () const
 
PndTrack GetPndTrack (int i)
 Returns the number of found tracks. More...
 
PndTrackCand GetPndTrackCand (int i)
 Returns the track cand. with index i. More...
 
PndFtsHoughTrackCand GetHoughTrack (int i) const
 Returns the track cand. with index i. More...
 
Int_t getNLinesBeforeDipoleFound () const
 
Int_t getNLinesBehindDipoleFound () const
 
Int_t getNParabolasFound () const
 
Int_t getNTracksFound () const
 

Protected Member Functions

void throwError (const TString s) const
 For error reporting. More...
 
Bool_t FilterTrackletsBasedOnSharedHits (UInt_t maxAcceptableSharedHits, std::vector< PndFtsHoughTracklet > &tracklets)
 Filters a vector of tracklets based on the number of shared hits. More...
 
std::vector< PndFtsHoughTrackletFindLinesBehindDipoleZx ()
 
std::vector< PndFtsHoughTrackletFindLinesBeforeDipoleZx ()
 
void FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole (const std::vector< PndFtsHoughTracklet > &trackletsLineBeforeDipole, const std::vector< PndFtsHoughTracklet > &trackletsLineBehindDipole)
 
Bool_t LineBehindDipoleMatchesToLinePlusParabola (const PndFtsHoughTrackCand &lineParabola, const PndFtsHoughTracklet &lineBehindDipole) const
 
void FindZyLineMatchingToLineParabolaLineInZx ()
 

Protected Attributes

PndFtsHoughTrackerTaskfTrackerTask
 Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the constructor. More...
 
PndFtsHoughTrackFinderParams fParams
 
std::vector< PndFtsHoughTrackCandfHoughTrackCandsComplete
 For internal storing of complete track cands. More...
 
std::vector< PndFtsHoughTrackCandfHoughTrackCandsZxPlaneOnly
 For internal storing of track cands. (zx plane track model only) More...
 
const UInt_t fMinPeakHeightZxLineBeforeDipole
 Minimum required height for peaks in Hough spaces. More...
 
const UInt_t fMinPeakHeightZxParabola
 Minimum required height for peaks in Hough spaces. More...
 
const UInt_t fMinPeakHeightZxLineBehindDipole
 Minimum required height for peaks in Hough spaces. More...
 
const UInt_t fMinPeakHeightZyLine
 zy line More...
 
Int_t fNLinesBeforeDipoleFound
 
Int_t fNLinesBehindDipoleFound
 
Int_t fNParabolasFound
 
Int_t fNTracksFound
 

Static Protected Attributes

static const Double_t fThetaRadLineBehindDipoleMatchesToParabolaIfBelow = 5*TMath::DegToRad()
 Minimum required height for peaks in Hough spaces. More...
 

Private Member Functions

 ClassDef (PndFtsHoughTrackFinder, 1)
 

Detailed Description

Implementation of the Hough transform based FTS PR. Creates Hough spaces, finds peaks (=tracklets) and combines them to track candidates.

Author
Martin J. Galuska <martin [dot] j [dot] galuska (at) physik [dot] uni (minus) giessen [dot] de>

This is a class version of the HoughTest.C macro PR test implementation minus all the plotting stuff. Take a look at the notes of the macro version for further details.

Recent Changes Major code cleanup and deletion of test code / unneeded code use PndFtsHoughTrackCand to store information about track candidates and Hough transforms Find all peaks with a minimum height, analysing peak shapes. Moved all Hough space related code to PndFtsHoughSpace -> Major code simplification, better maintainability Fill PndTrackCands and PndTrack for output

This class is loosely modeled after the mvd/MvdTracking/PndRiemannTrackFinder sttmvdtracking/PndMvdSttGemRiemannTrackFinder classes

TODO Match straight line for stations 5+6 to parabola Add skewed hits Add drift circles Adaptive Hough

Created: 18.06.2013

Definition at line 69 of file PndFtsHoughTrackFinder.h.

Constructor & Destructor Documentation

PndFtsHoughTrackFinder::PndFtsHoughTrackFinder ( PndFtsHoughTrackerTask trackerTask)

Set pointer to tracker task (super important as it provides an I/O interface to PandaRoot)

Definition at line 10 of file PndFtsHoughTrackFinder.cxx.

References fTrackerTask, and PndFtsHoughTrackerTask::GetVerbose().

11 : fTrackerTask(trackerTask),
12  fParams(),
13 
14  // min peak heights
19 
20  // event statistics
24  fNTracksFound(0)
25 {
26  if (0==fTrackerTask){
27  std::cout << "PndFtsHoughTrackFinder FATAL ERROR Tracker task not set.\n";
28  } else {
29  if(3<fTrackerTask->GetVerbose()) std::cout << "PndFtsHoughTrackFinder called with tracker ptr " << fTrackerTask << '\n';
30  }
31 }
const UInt_t fMinPeakHeightZxLineBehindDipole
Minimum required height for peaks in Hough spaces.
PndFtsHoughTrackFinderParams fParams
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
const UInt_t fMinPeakHeightZyLine
zy line
const UInt_t fMinPeakHeightZxParabola
Minimum required height for peaks in Hough spaces.
const UInt_t fMinPeakHeightZxLineBeforeDipole
Minimum required height for peaks in Hough spaces.
PndFtsHoughTrackFinder::~PndFtsHoughTrackFinder ( )
virtual

Destructor.

Definition at line 34 of file PndFtsHoughTrackFinder.cxx.

References PndFtsHoughTrackerTask::fLogger, fTrackerTask, and PndFtsHoughTrackerTask::GetVerbose().

35 {
36  if(3<fTrackerTask->GetVerbose()) fTrackerTask->fLogger->Info(MESSAGE_ORIGIN,"Destructor of PndFtsHoughTrackFinder");
37 }
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
FairLogger * fLogger
Returns pointer to the B field.

Member Function Documentation

PndFtsHoughTrackFinder::ClassDef ( PndFtsHoughTrackFinder  ,
 
)
private
Bool_t PndFtsHoughTrackFinder::FilterTrackletsBasedOnSharedHits ( UInt_t  maxAcceptableSharedHits,
std::vector< PndFtsHoughTracklet > &  tracklets 
)
protected

Filters a vector of tracklets based on the number of shared hits.

It will only keep the heighest peaks if 2 or more peaks share more than maxAcceptableSharedHits hits. If two peaks have the same height, both are kept.

Parameters
maxAcceptableSharedHitsDefines how many hits two tracklets / peaks are allowed to share.
[in,out]trackletsVector containing tracklets which are supposed to be filtered. The vector will be modified.
Returns

Definition at line 372 of file PndFtsHoughTrackFinder.cxx.

References Double_t, fTrackerTask, PndFtsHoughTracklet::getNSharedHits(), PndFtsHoughTracklet::getPeakHeightFromPeakFinder(), and PndFtsHoughTrackerTask::GetVerbose().

376 {
377  // to store indices that I want to "delete"
378  std::set<UInt_t> indicesToDelete;
379  std::set<UInt_t>::iterator findIndex;
380 
381  // for all tracklets that share more than maxAcceptableSharedHits hits keep only the heighest peak
382  // compare all tracklets with each other, keep track of which tracklets should be "deleted"
383  for (UInt_t iTrackletLeft = 0; iTrackletLeft < tracklets.size(); ++iTrackletLeft)
384  {
385  for (UInt_t iTrackletRight = iTrackletLeft+1; iTrackletRight < tracklets.size(); ++iTrackletRight)
386  {
387  // check if two tracks share more than maxAcceptableSharedHits hits
388  PndFtsHoughTracklet trackletLeft = tracklets[iTrackletLeft];
389  PndFtsHoughTracklet trackletRight = tracklets[iTrackletRight];
390  const UInt_t nSharedHits = trackletLeft.getNSharedHits(trackletRight);
391  if (nSharedHits>maxAcceptableSharedHits){
392  // mark the peak with the lower "height" for deletion
393  const Double_t heightLeft = trackletLeft.getPeakHeightFromPeakFinder();
394  const Double_t heightRight = trackletRight.getPeakHeightFromPeakFinder();
395  if (heightLeft>heightRight){
396  indicesToDelete.insert(iTrackletRight);
397  } else if (heightLeft<heightRight){
398  indicesToDelete.insert(iTrackletLeft);
399  } else if (heightLeft==heightRight){
400  std::cout << "WARNING: Found two peaks of the same height that share " << nSharedHits << " hits.\n";
401  std::cout << "The max. number of shared hits was set to " << maxAcceptableSharedHits << " hits.\n";
402  std::cout << "Both peaks / tracklets will be kept. \n";
403  }
404  }
405 
406  }
407  } // end for loop: compare all tracklets with each other
408 
409  // if we have no tracklets to "delete", we don't need to do anything
410  if (0==indicesToDelete.size()){
411  if (1<fTrackerTask->GetVerbose()) { std::cout << "FilterTrackletsBasedOnSharedHits: No tracklets are marked for deletion.\n"; }
412  return kTRUE;
413  }
414 
415  // create new tracklets vector and copy all entries from input vector to it which are not marked for deletion
416  std::vector<PndFtsHoughTracklet> filteredTracklets;
417  for (UInt_t iTracklet = 0; iTracklet < tracklets.size(); ++iTracklet)
418  {
419  // check if the index is marked for "deletion"
420  findIndex = indicesToDelete.find(iTracklet);
421  if (findIndex == indicesToDelete.end()){
422  // index was not found -> entry is not marked for deletion -> I should copy it
423  filteredTracklets.push_back(tracklets[iTracklet]);
424  }
425  }
426  // replace input tracklets with filtered ones
427  tracklets = filteredTracklets;
428  return kTRUE;
429 }
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
Class for saving the result of one Hough transform for FTS PR.
Double_t
UInt_t getNSharedHits(const PndFtsHoughTracklet &rhs)
Double_t getPeakHeightFromPeakFinder() const
std::vector< PndFtsHoughTracklet > PndFtsHoughTrackFinder::FindLinesBeforeDipoleZx ( )
protected

Definition at line 82 of file PndFtsHoughTrackFinder.cxx.

References Double_t, PndFtsHoughSpace::FillHoughSpace(), fMinPeakHeightZxLineBeforeDipole, fParams, fTrackerTask, and PndFtsHoughTrackFinderParams::getZLineParabola().

Referenced by FindTracks().

82  {
83  //----------------------------------
84  // zx plane: Straight line Hough transform before dipole
85 
86  // parameters: for hits with MC truth coordinates with 6 min peak height
87  // stepsPerThetaDeg = 16
88  // binsInSecVal = 1000;
89 
90 
91  // parameters: for hits with center of straws coordinates with 6 min peak height
92  // stepsPerThetaDeg = 16
93  // binsInSecVal = 1000;
94 
95 
96  // the following parameters are valid at z position fZLineParabola
97 
98  // min. and max. theta angles (in zx plane) for lines
99  static const Int_t thetaDegLowHigh = 25; // - for low, + for high (in degree)
100 
101  // min. and max. distance from z axis in x direction
102  static const Double_t xLowHigh = 80; // - for low, + for high (in cm)
103 
104  static const Int_t stepsPerThetaDeg = 8; // greater number means finer scanning in theta
105  static const Int_t nBinsInSecVal = 160; // greater number means more bins for second value of Hough space
106 
107 
108  PndFtsHoughSpaceBinning binningZxLineBeforeDipole(
109  stepsPerThetaDeg,
110  thetaDegLowHigh,
111  nBinsInSecVal,
112  xLowHigh
113  );
114  PndFtsHoughSpace houghSpaceZxLineBeforeDipole("lineBeforeDipole", -1,
115  binningZxLineBeforeDipole,
117 
118  // Do straight line Hough transform on non-skewed hits from stations 1+2
119  try { houghSpaceZxLineBeforeDipole.FillHoughSpace(); }
120  catch (std::runtime_error& e) {
121  std::cerr
122  << "Hough Space could not be created! \n";
123  std::cerr << "runtime_error: " << e.what() << '\n';
124  }
125  // find peaks for line Hough space and store in vector
126  // std::vector<PndFtsHoughTracklet> trackletsLineBeforeDipole = houghSpaceZxLineBeforeDipole.FindAllPeaksBinsWoMergingWithSearchWindow(fMinPeakHeightZxLineBeforeDipole);
127  std::vector<PndFtsHoughTracklet> trackletsLineBeforeDipole = houghSpaceZxLineBeforeDipole.FindAllPeaksScanPathsMergeBins(fMinPeakHeightZxLineBeforeDipole);
128  return trackletsLineBeforeDipole;
129 }
PndFtsHoughTrackFinderParams fParams
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
Double_t
Class for Hough space based on TH2S (for the moment). Saves the hits which enter this Hough space and...
const UInt_t fMinPeakHeightZxLineBeforeDipole
Minimum required height for peaks in Hough spaces.
Helper class for Hough space containing binning. Created: 09.02.2015.
void FillHoughSpace()
Fills the Hough space using the equation which corresponds to the name of the Hough space...
std::vector< PndFtsHoughTracklet > PndFtsHoughTrackFinder::FindLinesBehindDipoleZx ( )
protected

Definition at line 40 of file PndFtsHoughTrackFinder.cxx.

References Double_t, PndFtsHoughSpace::FillHoughSpace(), fMinPeakHeightZxLineBehindDipole, fParams, fTrackerTask, and PndFtsHoughTrackFinderParams::getZParabolaLine().

Referenced by FindTracks().

40  {
41  // zx plane: Straight line Hough transform behind dipole
42 
43 
44 
45 
46 
47  // the following parameters are valid at z position fZParabolaLine
48 
49  // min. and max. theta angles (in zx plane) for lines
50  static const Int_t thetaDegLowHigh = 80; // - for low, + for high (in degree)
51 
52  // min. and max. distance from z axis in x direction
53  static const Double_t xLowHigh = 200; // - for low, + for high (in cm)
54 
55  static const Int_t stepsPerThetaDeg = 8; // greater number means finer scanning in theta
56  static const Int_t nBinsInSecVal = 400; // greater number means more bins for second value of Hough space
57 
58 
59  PndFtsHoughSpaceBinning binningZxLineBehindDipole(
60  stepsPerThetaDeg,
61  thetaDegLowHigh,
62  nBinsInSecVal,
63  xLowHigh
64  );
65  PndFtsHoughSpace houghSpaceZxLineBehindDipole("lineBehindDipole", -1,
66  binningZxLineBehindDipole,
68 
69  // Do straight line Hough transform on non-skewed hits from stations 1+2
70  try { houghSpaceZxLineBehindDipole.FillHoughSpace(); }
71  catch (std::runtime_error& e) {
72  std::cerr
73  << "Hough Space for zx line behind dipole could not be created! \n";
74  std::cerr << "runtime_error: " << e.what() << '\n';
75  }
76  // find peaks for line hough space and store in vector
77  std::vector<PndFtsHoughTracklet> trackletsLineBehindDipole = houghSpaceZxLineBehindDipole.FindAllPeaksScanPathsMergeBins(
79  return trackletsLineBehindDipole;
80 }
const UInt_t fMinPeakHeightZxLineBehindDipole
Minimum required height for peaks in Hough spaces.
PndFtsHoughTrackFinderParams fParams
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
Double_t
Class for Hough space based on TH2S (for the moment). Saves the hits which enter this Hough space and...
Helper class for Hough space containing binning. Created: 09.02.2015.
void FillHoughSpace()
Fills the Hough space using the equation which corresponds to the name of the Hough space...
void PndFtsHoughTrackFinder::FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole ( const std::vector< PndFtsHoughTracklet > &  trackletsLineBeforeDipole,
const std::vector< PndFtsHoughTracklet > &  trackletsLineBehindDipole 
)
protected

Definition at line 131 of file PndFtsHoughTrackFinder.cxx.

References Bool_t, Double_t, fHoughTrackCandsZxPlaneOnly, PndFtsHoughSpace::FillHoughSpace(), fMinPeakHeightZxParabola, fParams, fTrackerTask, PndFtsHoughTrackerTask::GetEventNr(), PndFtsHoughTrackFinderParams::getParabolaHwScan(), PndFtsHoughTrackFinderParams::getParabolaNBinsPzxInv(), PndFtsHoughTrackFinderParams::getParabolaQDivPzxArgMax(), PndFtsHoughTrackFinderParams::getParabolaStepsPerThetaDeg(), PndFtsHoughTrackerTask::GetVerbose(), PndFtsHoughTrackFinderParams::getZLineParabola(), LineBehindDipoleMatchesToLinePlusParabola(), PndFtsHoughTrackCand::SetZxLineBeforeDipole(), PndFtsHoughTrackCand::SetZxLineBehindDipole(), and PndFtsHoughTrackCand::SetZxParabola().

Referenced by FindTracks().

134  {
135  // loop over all line tracklets which were found by line HT before dipole field and find a matching parabola
136  // After that try to match a line after dipole to that track candidate.
137  // If that is not possible, accept short track candidate if particle is likely to have left detector acceptance, otherwise discard.
138  // LB4D = line before dipole
139 
140  // loop over lines before dipole
141  for (UInt_t iLB4D = 0; iLB4D < linesBeforeDipole.size(); ++iLB4D) {
142  const Double_t peakThetaRadLB4D = linesBeforeDipole[iLB4D].getThetaRadVal();
143  const Double_t peakInterceptLB4D = linesBeforeDipole[iLB4D].getSecondVal();
144  const Double_t peakThetaRadHwLB4D = linesBeforeDipole[iLB4D].getThetaRadHw();
145 
146  // determine where to look for parabola
147  const Double_t parabolaStepsPerThetaDeg = fParams.getParabolaStepsPerThetaDeg();
148  const Double_t parabolaHwScan = fParams.getParabolaHwScan();
149 
150  const Double_t parabolaThetaRadLow = peakThetaRadLB4D
151  - parabolaHwScan*peakThetaRadHwLB4D; // in rad
152  const Double_t parabolaThetaRadHigh = peakThetaRadLB4D
153  + parabolaHwScan*peakThetaRadHwLB4D; // in rad
154  if (parabolaThetaRadHigh == parabolaThetaRadLow) std::cout << "ERROR: low and high are the same for parabola!\n";
155  const Double_t parabolaThetaDegLow = parabolaThetaRadLow*TMath::RadToDeg(); // in deg
156  const Double_t parabolaThetaDegHigh = parabolaThetaRadHigh*TMath::RadToDeg(); // in deg
157 
158  UInt_t thetaBins = ceil(parabolaStepsPerThetaDeg * (parabolaThetaDegHigh - parabolaThetaDegLow));
159 
160  if (1<fTrackerTask->GetVerbose()) std::cout << "event: " << fTrackerTask->GetEventNr()
161  << " Line " << iLB4D
162  << "\n thetaRadLowParabola=" << parabolaThetaRadLow
163  << " thetaRadHighParabola=" << parabolaThetaRadHigh
164  << " thetaBins=" << thetaBins
165  << " peakThetaRadLB4D="
166  << peakThetaRadLB4D
167  << " peakThetaRadHwLB4D="
168  << peakThetaRadHwLB4D
169  << " peakInterceptLB4D=" << peakInterceptLB4D
170  << " peakThetaRadHwLB4D=" << peakThetaRadHwLB4D
171  << '\n';
172 
173  const Int_t parabolaNBinsPzxInv = fParams.getParabolaNBinsPzxInv();
174  const Double_t parabolaQDivPzxArgMax = fParams.getParabolaQDivPzxArgMax();
175 
176 
177  PndFtsHoughSpaceBinning binningZxParabola(
178  thetaBins,
179  parabolaThetaRadLow,
180  parabolaThetaRadHigh,
181  parabolaNBinsPzxInv,
182  -1*parabolaQDivPzxArgMax,
183  parabolaQDivPzxArgMax
184  );
185  PndFtsHoughSpace houghspaceZxParabola("parabola", iLB4D,
186  binningZxParabola,
187  fParams.getZLineParabola(), peakInterceptLB4D, 0, fTrackerTask);
188 
189  // Do parabola Hough transform for current line before dipole (shifts FTS hits by hitshiftinx) for non-skewed hits in stations 3+4+5
190  try { houghspaceZxParabola.FillHoughSpace(); }
191  catch (std::runtime_error& e) {
192  std::cerr << "Hough Space for zx parabola could not be created! \n";
193  std::cerr << "runtime_error: " << e.what() << '\n';
194  }
195 
196  // find peaks for zx parabola Hough space (for current line before dipole)
197  // std::vector<PndFtsHoughTracklet> zxParabolaTracklets = houghspaceZxParabola.FindAllPeaksScanPathsMergeBins(fMinPeakHeightZxParabola);
198  std::vector<PndFtsHoughTracklet> zxParabolaTracklets = houghspaceZxParabola.FindAllPeaksBinsWoMergingWithSearchWindow(fMinPeakHeightZxParabola);
199 
200 
201 
202 
203  // construct newHoughTrackCand from line+parabola+line (in zx)
204  if (1 < fTrackerTask->GetVerbose()) {
205  std::cout
206  << "Create track candidates from line in zx before dipole "
207  << iLB4D << " and all parabolas within dipole.\n"
208  << "And all lines behind dipole which match line+parabola (if exist)\n";
209  }
210  // loop over all parabolas which were found for one line before the dipole
211  for (UInt_t iParabola = 0; iParabola < zxParabolaTracklets.size(); ++iParabola) {
212  // Try to find a matching line behind the dipole
213 
214  PndFtsHoughTrackCand lineParabola(fTrackerTask);
215  lineParabola.SetZxLineBeforeDipole(linesBeforeDipole[iLB4D]);
216  lineParabola.SetZxParabola(zxParabolaTracklets[iParabola]);
217 
218  // loop over all lines behind dipole
219  // add track cands for each matching line behind dipole (LBhD)
220  Bool_t foundLineBehindDipole = kFALSE;
221  for (UInt_t iLBhD = 0; iLBhD < linesBehindDipole.size(); ++iLBhD) {
222  if (kTRUE == LineBehindDipoleMatchesToLinePlusParabola( lineParabola, linesBehindDipole[iLBhD]) ){
223  PndFtsHoughTrackCand lineParabolaLine(lineParabola); // make a copy
224  lineParabolaLine.SetZxLineBehindDipole(linesBehindDipole[iLBhD]); // add line behind dipole
225  fHoughTrackCandsZxPlaneOnly.push_back(lineParabolaLine); // add as track candidate
226  foundLineBehindDipole = kTRUE;
227  }
228  }
229 
230  // iif no lines behind dipole matched, add track cand without any line behind dipole
231  if (kFALSE == foundLineBehindDipole) fHoughTrackCandsZxPlaneOnly.push_back(lineParabola);
232  }
233  } // loop over all line tracklets which were found by line HT before dipole field
234 }
PndFtsHoughTrackFinderParams fParams
Bool_t LineBehindDipoleMatchesToLinePlusParabola(const PndFtsHoughTrackCand &lineParabola, const PndFtsHoughTracklet &lineBehindDipole) const
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
const UInt_t fMinPeakHeightZxParabola
Minimum required height for peaks in Hough spaces.
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsZxPlaneOnly
For internal storing of track cands. (zx plane track model only)
Double_t
Class for Hough space based on TH2S (for the moment). Saves the hits which enter this Hough space and...
Helper class for Hough space containing binning. Created: 09.02.2015.
Class for saving a FTS track cand. for Hough transform based FTS PR.
UInt_t GetEventNr() const
Returns the save debug flag.
void FillHoughSpace()
Fills the Hough space using the equation which corresponds to the name of the Hough space...
void PndFtsHoughTrackFinder::FindTracks ( )
virtual

Performs the track finding.

Reimplemented in PndFtsHoughTrackFinderQA.

Definition at line 295 of file PndFtsHoughTrackFinder.cxx.

References fHoughTrackCandsComplete, fHoughTrackCandsZxPlaneOnly, FindLinesBeforeDipoleZx(), FindLinesBehindDipoleZx(), FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole(), FindZyLineMatchingToLineParabolaLineInZx(), fMinPeakHeightZxLineBeforeDipole, fMinPeakHeightZxParabola, fNLinesBeforeDipoleFound, fNLinesBehindDipoleFound, fNParabolasFound, fNTracksFound, fTrackerTask, PndFtsHoughTrackerTask::GetEventNr(), PndFtsHoughTrackerTask::GetNFtsHits(), and PndFtsHoughTrackerTask::GetVerbose().

Referenced by PndFtsHoughTrackerTaskQA::Exec(), and PndFtsHoughTrackFinderQA::FindTracks().

295  {
296 
297  if (1<fTrackerTask->GetVerbose()) std::cout << "PndFtsHoughTrackFinder::FindTracks()\n";
298 
299  // reset
301  fHoughTrackCandsComplete.clear();
302 
303 
304  //--------------------------
305  // Do we have enough hits in the FTS?
307  if( minFtsHits > fTrackerTask->GetNFtsHits() ) {
308  if(0<fTrackerTask->GetVerbose()) std::cout << "Skip event " << fTrackerTask->GetEventNr() << " since we have too few hits in FTS. (" << fTrackerTask->GetNFtsHits() << " present, " << minFtsHits << " required)\n";
309  return;
310  }
311 
312 
313 
314 
315  // zx plane: Straight line Hough transform behind dipole
316  if(1<fTrackerTask->GetVerbose()) std::cout << "Lines behind dipole:\n";
317  std::vector<PndFtsHoughTracklet> linesBehindDipole = FindLinesBehindDipoleZx();
318  fNLinesBehindDipoleFound = linesBehindDipole.size();
319  if(1<fTrackerTask->GetVerbose()) std::cout << "event " << fTrackerTask->GetEventNr() << ": " << fNLinesBehindDipoleFound << " lines behind dipole found.\n";
320 
321  // zx plane: Straight line Hough transform before dipole
322  if(1<fTrackerTask->GetVerbose()) std::cout << "Lines before dipole:\n";
323  std::vector<PndFtsHoughTracklet> linesBeforeDipole = FindLinesBeforeDipoleZx();
324  fNLinesBeforeDipoleFound = linesBeforeDipole.size();
325  if(1<fTrackerTask->GetVerbose()) std::cout << "event " << fTrackerTask->GetEventNr() << ": " << fNLinesBeforeDipoleFound << " lines before dipole found.\n";
326 
327  // loop over all line tracklets which were found by line HT before dipole field and find a matching parabola
328  if(1<fTrackerTask->GetVerbose()) std::cout << "Matching parabolas to lines before dipole:\n";
329  FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole(linesBeforeDipole, linesBehindDipole);
331  if(0<fTrackerTask->GetVerbose()) std::cout << "event " << fTrackerTask->GetEventNr() << ": " << fNParabolasFound << " parabolas matched!\n";
332 
333  // loop over line+parabola+line (LPL) from zx plane and do straight line Hough transform in zy plane
334  if (1<fTrackerTask->GetVerbose()) std::cout << "Lines in zy plane:\n";
337  if(1<fTrackerTask->GetVerbose()) std::cout << "event " << fTrackerTask->GetEventNr() << ": " << fNTracksFound << " tracks found.\n";
338 
339 
340 
341  // HIER GEHT ES WEITER!
342 
343 
344 
345  // TODO: Convert momenta to physical units
346  // pzinvpeak = 1./pzinvpeak*0.00299792458;
347  // pzinvpeakWithBField = pzinvpeak/BMeanForParabola;
348 
349 
350  // if(1<fTrackerTask->GetVerbose()) {
351  //
352  // // TODO Add missing output
353  //
354  //
355  //
356  //
357  // // if (kTRUE == correctpz)
358  // // {
359  // // std::cout << "Watch out! p_z value in plot is already corrected!!!" << std::endl;
360  // // }
361  // }
362 }
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
std::vector< PndFtsHoughTracklet > FindLinesBehindDipoleZx()
Int_t GetNFtsHits() const
Returns the event number.
const UInt_t fMinPeakHeightZxParabola
Minimum required height for peaks in Hough spaces.
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsZxPlaneOnly
For internal storing of track cands. (zx plane track model only)
const UInt_t fMinPeakHeightZxLineBeforeDipole
Minimum required height for peaks in Hough spaces.
UInt_t GetEventNr() const
Returns the save debug flag.
std::vector< PndFtsHoughTracklet > FindLinesBeforeDipoleZx()
void FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole(const std::vector< PndFtsHoughTracklet > &trackletsLineBeforeDipole, const std::vector< PndFtsHoughTracklet > &trackletsLineBehindDipole)
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsComplete
For internal storing of complete track cands.
void PndFtsHoughTrackFinder::FindZyLineMatchingToLineParabolaLineInZx ( )
protected

Definition at line 236 of file PndFtsHoughTrackFinder.cxx.

References fHoughTrackCandsComplete, fHoughTrackCandsZxPlaneOnly, PndFtsHoughSpace::FillHoughSpace(), PndFtsHoughSpace::FindAllPeaksScanPathsMergeBins(), fMinPeakHeightZyLine, fParams, fTrackerTask, PndFtsHoughTrackerTask::GetVerbose(), PndFtsHoughTrackFinderParams::getZLineParabola(), and PndFtsHoughTrackCand::SetZyLine().

Referenced by FindTracks().

236  {
237  // loop over line+parabola+line track candidates in zx plane, pass track cand. to Hough space (one Hough transform needed for each line+parabola+line track cand.)
238  // Check all skewed hits whether they might belong to any of the track candidates
239  // Calculate (x,z) coordinate hypotheses for skewed hits which might belong to track candidates and run line Hough transform on these hypotheses
240  //----------------------------------
241  // zy plane: Straight line Hough transform
242  // loop over line+parabola+line (LPL) from zx plane
243  for (UInt_t iLPL = 0; iLPL < fHoughTrackCandsZxPlaneOnly.size(); ++iLPL) {
244  // determine where to look for line in zy plane
245  static const Int_t stepsPerThetaDegZyLine = 4; // greater number means finer scanning in theta
246  static const Int_t thetaDegLowZyLine = -18; // in degree
247  static const Int_t thetaDegHighZyLine = 18; // in degree
248  // create Hough space
249  PndFtsHoughSpaceBinning binningZyLine(
250  stepsPerThetaDegZyLine * (thetaDegHighZyLine - thetaDegLowZyLine),
251  thetaDegLowZyLine * TMath::DegToRad(), // in rad
252  thetaDegHighZyLine * TMath::DegToRad(), // in rad
253  stepsPerThetaDegZyLine * 16, // TODO: Check values
254  -80., // in cm // TODO: Check values
255  80. // in cm
256  );
257  PndFtsHoughSpace houghspaceZyLine("lineZy", iLPL,
258  binningZyLine,
260  fTrackerTask);
261  // Do line Hough transform for current line+parabola+line for skewed hits in all stations
262  try {
263  houghspaceZyLine.FillHoughSpace();
264  } catch (std::runtime_error& e) {
265  std::cerr
266  << "Hough Space for zy line before dipole could not be created! \n";
267  std::cerr << "runtime_error: " << e.what() << '\n';
268  }
269  // find peaks for line Hough space and store in vector
270  std::vector<PndFtsHoughTracklet> trackletsZyLine =
271  houghspaceZyLine.FindAllPeaksScanPathsMergeBins(
273  // construct newHoughTrackCand by adding zy line to line+parabola+line (in zx)
274  if (1 < fTrackerTask->GetVerbose()) {
275  std::cout
276  << "Create track candidates by adding line in zy to line+parabola+line in zx "
277  << iLPL << "\n"
278  << "If no zy line can be found, delete track candidate.\n";
279  }
280  // loop over all zy lines which were found for line+parabola+line
281  for (UInt_t iZyLine = 0; iZyLine < trackletsZyLine.size(); ++iZyLine) {
282  PndFtsHoughTrackCand fullTrackCand(
283  fHoughTrackCandsZxPlaneOnly[iLPL]); // make copy
284  fullTrackCand.SetZyLine(trackletsZyLine[iZyLine]); // add zy line
285  fHoughTrackCandsComplete.push_back(fullTrackCand); // add as track candidate
286  }
287  } // loop over all line+parabola+line from zx plane
288 }
PndFtsHoughTrackFinderParams fParams
PndFtsHoughTrackerTask * fTrackerTask
Task which handles PandaRoot input/output and provides settings for FTS PR. Has to be set using the c...
const UInt_t fMinPeakHeightZyLine
zy line
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsZxPlaneOnly
For internal storing of track cands. (zx plane track model only)
Class for Hough space based on TH2S (for the moment). Saves the hits which enter this Hough space and...
Helper class for Hough space containing binning. Created: 09.02.2015.
Class for saving a FTS track cand. for Hough transform based FTS PR.
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsComplete
For internal storing of complete track cands.
PndFtsHoughTrackCand PndFtsHoughTrackFinder::GetHoughTrack ( int  i) const
inline

Returns the track cand. with index i.

Note: For debugging only.

Parameters
iIndex of requested track cand.
Returns
Track cand. to index i.

Definition at line 103 of file PndFtsHoughTrackFinder.h.

References fHoughTrackCandsComplete, and i.

Referenced by PndFtsHoughTrackerTask::Exec().

103 { return fHoughTrackCandsComplete[i]; };
Int_t i
Definition: run_full.C:25
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsComplete
For internal storing of complete track cands.
Int_t PndFtsHoughTrackFinder::getNLinesBeforeDipoleFound ( ) const
inline

Definition at line 106 of file PndFtsHoughTrackFinder.h.

References fNLinesBeforeDipoleFound.

Int_t PndFtsHoughTrackFinder::getNLinesBehindDipoleFound ( ) const
inline

Definition at line 107 of file PndFtsHoughTrackFinder.h.

References fNLinesBehindDipoleFound.

Int_t PndFtsHoughTrackFinder::getNParabolasFound ( ) const
inline

Definition at line 108 of file PndFtsHoughTrackFinder.h.

References fNParabolasFound.

Referenced by PndFtsHoughTrackerTaskQA::Exec().

Int_t PndFtsHoughTrackFinder::getNTracksFound ( ) const
inline

Definition at line 109 of file PndFtsHoughTrackFinder.h.

References fNTracksFound.

109 { return fNTracksFound; };
PndTrack PndFtsHoughTrackFinder::GetPndTrack ( int  i)
inline

Returns the number of found tracks.

Returns the track cand. with index i.

Note: Method calculates first and last parameters of the PndTrack object, but uses an empty PndTrackCand which has to be set lateron using SetTrackCandRef!

Parameters
iIndex of requested track cand.
Returns
Track cand. to index i.

Definition at line 89 of file PndFtsHoughTrackFinder.h.

References fHoughTrackCandsComplete, and i.

Referenced by PndFtsHoughTrackerTask::Exec().

89 { return fHoughTrackCandsComplete[i].getPndTrack(); };
Int_t i
Definition: run_full.C:25
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsComplete
For internal storing of complete track cands.
PndTrackCand PndFtsHoughTrackFinder::GetPndTrackCand ( int  i)
inline

Returns the track cand. with index i.

Note: Use this to add a PndTrackCand to the corresponding PndTrack object with SetTrackCandRef!

Parameters
iIndex of requested track cand.
Returns
Track cand. to index i.

Definition at line 96 of file PndFtsHoughTrackFinder.h.

References fHoughTrackCandsComplete, and i.

Referenced by PndFtsHoughTrackerTask::Exec().

96 { return fHoughTrackCandsComplete[i].getPndTrackCand(); };
Int_t i
Definition: run_full.C:25
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsComplete
For internal storing of complete track cands.
Bool_t PndFtsHoughTrackFinder::LineBehindDipoleMatchesToLinePlusParabola ( const PndFtsHoughTrackCand lineParabola,
const PndFtsHoughTracklet lineBehindDipole 
) const
inlineprotected

Definition at line 189 of file PndFtsHoughTrackFinder.h.

References Double_t, fabs(), fThetaRadLineBehindDipoleMatchesToParabolaIfBelow, PndFtsHoughTracklet::getThetaRadVal(), PndFtsHoughTrackCand::getThetaZyRad(), and PndFtsHoughTracklet::getZRefLabSys().

Referenced by FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole().

192  {
193 
194  // the angles should be compared at the z coordinate where I switch from parabola to line behind dipole
195  Double_t zParabolaLine = lineBehindDipole.getZRefLabSys();
196  return fabs( lineParabola.getThetaZyRad(zParabolaLine) - lineBehindDipole.getThetaRadVal() ) < fThetaRadLineBehindDipoleMatchesToParabolaIfBelow;
197 }
Double_t getThetaRadVal() const
Double_t getZRefLabSys() const
Double_t
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
static const Double_t fThetaRadLineBehindDipoleMatchesToParabolaIfBelow
Minimum required height for peaks in Hough spaces.
Double_t getThetaZyRad(const Double_t) const
Int_t PndFtsHoughTrackFinder::NTracks ( ) const
inline

Definition at line 82 of file PndFtsHoughTrackFinder.h.

References fHoughTrackCandsComplete.

Referenced by PndFtsHoughTrackerTask::Exec().

82 { return fHoughTrackCandsComplete.size(); };
std::vector< PndFtsHoughTrackCand > fHoughTrackCandsComplete
For internal storing of complete track cands.
void PndFtsHoughTrackFinder::OverwriteTrackFinderParams ( PndFtsHoughTrackFinderParams  newParams)
inline

Definition at line 76 of file PndFtsHoughTrackFinder.h.

References fParams.

Referenced by PndFtsHoughTrackerTaskQA::Exec().

76 { fParams = newParams; };
PndFtsHoughTrackFinderParams fParams
void PndFtsHoughTrackFinder::throwError ( const TString  s) const
inlineprotected

For error reporting.

Definition at line 117 of file PndFtsHoughTrackFinder.h.

117 { throw std::runtime_error(s.Data()); };
TLorentzVector s
Definition: Pnd2DStar.C:50

Member Data Documentation

std::vector<PndFtsHoughTrackCand> PndFtsHoughTrackFinder::fHoughTrackCandsComplete
protected

For internal storing of complete track cands.

Definition at line 129 of file PndFtsHoughTrackFinder.h.

Referenced by FindTracks(), FindZyLineMatchingToLineParabolaLineInZx(), GetHoughTrack(), GetPndTrack(), GetPndTrackCand(), and NTracks().

std::vector<PndFtsHoughTrackCand> PndFtsHoughTrackFinder::fHoughTrackCandsZxPlaneOnly
protected

For internal storing of track cands. (zx plane track model only)

Definition at line 130 of file PndFtsHoughTrackFinder.h.

Referenced by FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole(), FindTracks(), and FindZyLineMatchingToLineParabolaLineInZx().

const UInt_t PndFtsHoughTrackFinder::fMinPeakHeightZxLineBeforeDipole
protected

Minimum required height for peaks in Hough spaces.

zx line before dipole field

Definition at line 137 of file PndFtsHoughTrackFinder.h.

Referenced by FindLinesBeforeDipoleZx(), and FindTracks().

const UInt_t PndFtsHoughTrackFinder::fMinPeakHeightZxLineBehindDipole
protected

Minimum required height for peaks in Hough spaces.

zx line after dipole field

Definition at line 141 of file PndFtsHoughTrackFinder.h.

Referenced by FindLinesBehindDipoleZx().

const UInt_t PndFtsHoughTrackFinder::fMinPeakHeightZxParabola
protected

Minimum required height for peaks in Hough spaces.

zx parabola within dipole field

Definition at line 139 of file PndFtsHoughTrackFinder.h.

Referenced by FindMatchingParabolaToLineBeforeDipoleZxAndAddLineBehindDipole(), and FindTracks().

const UInt_t PndFtsHoughTrackFinder::fMinPeakHeightZyLine
protected

zy line

Definition at line 143 of file PndFtsHoughTrackFinder.h.

Referenced by FindZyLineMatchingToLineParabolaLineInZx().

Int_t PndFtsHoughTrackFinder::fNLinesBeforeDipoleFound
protected

Definition at line 148 of file PndFtsHoughTrackFinder.h.

Referenced by FindTracks(), and getNLinesBeforeDipoleFound().

Int_t PndFtsHoughTrackFinder::fNLinesBehindDipoleFound
protected

Definition at line 149 of file PndFtsHoughTrackFinder.h.

Referenced by FindTracks(), and getNLinesBehindDipoleFound().

Int_t PndFtsHoughTrackFinder::fNParabolasFound
protected

Definition at line 150 of file PndFtsHoughTrackFinder.h.

Referenced by FindTracks(), and getNParabolasFound().

Int_t PndFtsHoughTrackFinder::fNTracksFound
protected

Definition at line 151 of file PndFtsHoughTrackFinder.h.

Referenced by FindTracks(), and getNTracksFound().

PndFtsHoughTrackFinderParams PndFtsHoughTrackFinder::fParams
protected
const Double_t PndFtsHoughTrackFinder::fThetaRadLineBehindDipoleMatchesToParabolaIfBelow = 5*TMath::DegToRad()
staticprotected

Minimum required height for peaks in Hough spaces.

Definition at line 134 of file PndFtsHoughTrackFinder.h.

Referenced by LineBehindDipoleMatchesToLinePlusParabola().

PndFtsHoughTrackerTask* PndFtsHoughTrackFinder::fTrackerTask
protected

The documentation for this class was generated from the following files: