PandaRoot
PndFtsHoughSpace.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
37 #ifndef PndFtsHoughSpace_H
38 #define PndFtsHoughSpace_H
39 
40 #include "PndFtsHoughTrackerTask.h"
41 #include "PndFtsHoughSpacePeak.h"
43 
44 #include "TH2.h"
45 #include <cmath>
46 #include <vector>
47 #include "Rtypes.h" // for Double_t, Int_t, etc
48 #include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
49 
50 #include "PndTrackCandHit.h"
51 #include "TClonesArray.h"
52 
53 // magnetic field
54 #include "FairField.h"
55 #include "TVector3.h"
56 
57 #include "PndFtsHoughTracklet.h"
58 #include "PndFtsHoughTrackCand.h"
59 
60 // For error throwing
61 #include "TString.h"
62 #include <stdexcept>
63 
64 // for saving the path through the Hough space
65 typedef std::vector<Int_t> IdxPath; // helper -- path for one hit
66 typedef std::map<Int_t, IdxPath> HitIdxPathMap; // that is the one I need -- map of hit indices to path for the corresponding hits
67 typedef std::pair<Int_t, IdxPath> HitIdxPathPair; // helper for inserting one pair of hit index and path in map
68 
69 // cout for the above types
70 std::ostream &operator<<(std::ostream &os, const IdxPath &outVector);
71 std::ostream &operator<<(std::ostream &os, const HitIdxPathMap &outMap);
72 
73 class PndFtsHoughSpace : public TH2S {
74  public:
75  // Constructors/Destructors ---------
76  PndFtsHoughSpace(const char *name = nullptr, const Int_t refIndex = -1,
77 
79 
80  Double_t zRefPos = 0., Double_t interceptZx = 0.,
81 
82  PndFtsHoughTrackCand *associatedTrackCand = nullptr,
83 
84  PndFtsHoughTrackerTask *trackerTask = nullptr);
86 
87  // TH2S ExportTH2S();
88 
89  // General info for peak finders
90  // The peaks contain the following information
91  // peakTheta = theta for peak
92  // peakSecond = Q/pzx for peak (parabola HT)
93  // peakSecond = intercept for peak (line HT) (in z-x- or z-y-plane)
94  //
95  // peakThetaHw = half width of peak in theta
96  // peakSecondHw = half width of peak in second value (see above for what it stands for)
97  //
98  // actualHeight height of the peak in the histogram (in counts)
99 
105  std::vector<PndFtsHoughTracklet> FindAllPeaksScanPathsMergeBins(const UInt_t minHeight);
106 
112  std::vector<PndFtsHoughTracklet> FindAllPeaksScanPathsMergeBinsCalculatingPaths(const UInt_t minHeight);
113 
119  std::vector<PndFtsHoughTracklet> FindAllPeaksWithTSpectrum2(const UInt_t minHeight);
120 
126  std::vector<PndFtsHoughTracklet> FindAllPeaksBinsWoMergingWithSearchWindow(const UInt_t minHeight, const Int_t vicinityLength = 0);
127 
133  std::vector<PndFtsHoughTracklet> FindAllPeaksBlanko(const UInt_t minHeight);
134 
151  void FillHoughSpace();
152 
153  // operators
154  // PndFtsHoughSpace are the same if they contain the same hits, that means if the PndTrackCand are the same, therefore no need to implement that operator here
155 
156  // Accessors -----------------------
157  inline void Print() const;
158 
159  void setVerbose(Int_t verbose) { fVerbose = verbose; };
160 
161  Double_t getInterceptZx() const { return fInterceptZx; };
162 
163  Double_t getZRefPos() const { return fZRefPos; };
164 
165  inline UInt_t GetNHits() const { return fHitId.size(); };
166 
167  private:
168  // for PandaRoot input/output
169  PndFtsHoughTrackerTask *fTrackerTask;
170 
172  void throwError(const TString s) const { throw std::runtime_error(s.Data()); };
173 
174  Bool_t setParametersForHsOption(); // set parameters according to the kind of Hough transform I want to do
175  void filterInputHits(); // copies input hits (based on z coordinate and skewed/non-skewed) from fFtsHitArray (all FTS hits) to fHitId (only the hits that qualify for the specific
176  // Hough transform)
177  inline void AddHitToHS(UInt_t hitId, Double_t rho);
178  inline void AddHitToHS(FairLink link, Double_t rho);
179  inline Bool_t IsHitFromTubeIdAlreadyAdded(const Int_t tubeIdToAdd);
180  inline const PndFtsHit *getHitFromHS(UInt_t index) const; // gets the FTS hit corresponding to index
181  inline Int_t getHitIdFromHS(UInt_t index) const { return fHitId.at(index).GetHitId(); }; // gets the FTS hit Id corresponding to index
182  inline const TVector3 CalculateHitPosFromIntersectionsWithZxTrackModel(const PndFtsHit *const myHit) const;
183  inline const TVector3 GetRawOrCalculatedHitPos(const PndFtsHit *const myHit) const;
184 
185  void AddHitsToTrackletByCalculating(PndFtsHoughTracklet *currentTracklet, Int_t locmax);
186 
192  HitIdxPathMap fHitThetaYIdxPath;
193 
194  // Private Data Members ------------
195  Int_t fVerbose;
196  const Int_t fRefIndex;
197 
199  Int_t fFtsBranchId;
200  std::vector<PndTrackCandHit> fHitId;
201 
203  PndFtsHoughTrackCand *fAssociatedTrackCand;
204 
205  // only hits with a z value (in the laboratory system, zreal) between onlyusehitsfromz and onlyusehitsuptoz will be used for building the houghspace
206  Double_t fOnlyUseHitsFromZ;
207  Double_t fOnlyUseHitsUpToZ;
208  Bool_t fUseNonSkewedStraws; // if kTRUE, then hits from non-skewed straws are used for Hough transform
209  Bool_t fUseSkewedStraws; // if kTRUE, then hits from skewed straws are used for Hough transform
210 
211  Bool_t fKeepBConstant; // set only to kFALSE for testing
212  // If kTRUE the y-component of the B-field is not used in the parabola hough transform
213  // if kFALSE the parabola's shape will be adjusted based on the magnetic field maps
214 
215  // at which z value the values are calculated
216  Double_t fZRefPos; // is used to redefine an origin for the coordinate system (so that the angle definition gives meaningful theta values)
217  // Double_t interceptZx; // cannot be constant, because might need to be reset if set incorrectly (has to be 0 for line HT)
218  // is used to shift the true x values of hits so that they hit the point (zOffset|0) in z-x-plane (value is determined by line fit on chambers1+2)
219  // (zreal=zOffset, xreal=interceptZx) = (zshifted = 0, xshifted = 0)
220  // zshifted = zreal - zRefPos
221  // xshifted = xreal - interceptZx
222  // zreal = zshifted + zRefPos
223  // xreal = xshifted + interceptZx
224  // For the z-x-plane parabola, a shift in x (hitshiftinx) needs to be set (which should be the result of the straight line hough transform before the dipole field)
225  // For the straight line (stations before dipole field) hitshiftinx HAS TO BE ZERO
226 
227  // is used only for parabola in zx plane to shift the true x values of hits so that they hit the point (zOffset|0) in z-x-plane (value is determined by line fit on chambers1+2)
228  Double_t fInterceptZx;
229 
230  // for B field access
231  FairField *fField;
232 
233  // for HoughTransform
235  // if kTRUE will correct the pzx prediction according to values which should be obtained from a line fit mc truth momentum VS. reco momentum with high statistics
236  static const Bool_t fCorrectPzx = kFALSE;
237 
252  inline Bool_t FillHoles(const Int_t lastBinX, const Int_t lastBinY, const Int_t currentBinY, IdxPath *ptrThetaYIdxPathVec);
253 
254  inline Double_t equationParabola(Double_t thetaRad, Double_t hitZShifted, Double_t hitXShifted, Double_t By)
255  {
256  // for parabola equation
257  const Double_t n = 1.;
258  const Double_t e = 1.;
259  const Double_t c = n * e / 2.;
260 
261  Double_t yVal = 1. / c / By * (-hitZShifted * sin(thetaRad) + hitXShifted * cos(thetaRad)) / pow((hitZShifted * cos(thetaRad) + hitXShifted * sin(thetaRad)), 2);
262 
263  // next line is with rotation as in paper (I believe it is incorrect)
264  // value = 1. / c / By * (hitZshifted * sin(realtheta) - hitXshifted * cos(realtheta))/ pow((hitZshifted * cos(realtheta) + hitXshifted * sin(realtheta)), 2);
265  return yVal;
266  };
267 
268  inline Double_t equationParabolaPz(Double_t thetaRad, Double_t hitZShifted, Double_t hitXShifted, Double_t By)
269  {
270  // for parabola equation
271  const Double_t n = 1.;
272  const Double_t e = 1.;
273  const Double_t c = n * e / 2.;
274 
275  // Use shifted x and shifted z for parabola
276  Double_t yVal = c * By * pow((hitZShifted * cos(thetaRad) + hitXShifted * sin(thetaRad)), 2) / (-hitZShifted * sin(thetaRad) + hitXShifted * cos(thetaRad));
277 
278  // next line is with rotation as in paper (I believe it is incorrect)
279  // value = c*By*pow((hitZshifted * cos(realtheta) + hitXshifted * sin(realtheta)), 2)/(hitZshifted * sin(realtheta) - hitXshifted * cos(realtheta));
280 
281  return yVal;
282  };
283 
284  inline Double_t equationLineZxOrZy(Double_t thetaRad, Double_t hitZShifted, Double_t hitXLabSys)
285  { // TODO: Check if that also works in zy plane
286  // calculate b which is the distance of point on line at z = zOffset from z axis
287  Double_t yVal = -tan(thetaRad) * hitZShifted + hitXLabSys;
288  return yVal;
289  }
290 
291  //------
292  // DEBUG
293  //-----
296  inline TString GetDebugOutName(TString title = "", Int_t param = -1) const;
299  void WriteHistoOfHoughSpace() const;
302  TH2S MakeEmptyHistoOfSameDimensions(TString specifier = "", Int_t index = -1) const;
305  void WriteHistoOfAllPeaks(const std::vector<PndFtsHoughSpacePeak> &peaksToPlot) const;
308  void WriteHistoOfAllPaths() const;
312  void WriteHistoOfAllPathsForEachMcTruthTrack() const;
313  inline void PrintFoundTracklets(const std::vector<PndFtsHoughTracklet> &tracklets) const;
314  inline Double_t getByFromBField(Double_t hitXLabSys, Double_t hitYLabSys, Double_t hitZLabSys);
315 
316  public:
318 };
319 
320 // inline functions
322 {
323  std::cout << "=========== PndFtsHoughSpace::Print() ==========" << '\n';
324  std::cout << "fZRefPos = " << fZRefPos << '\n';
325  std::cout << "fInterceptZx = " << fInterceptZx << 'n';
326  TH2S::Print();
327 }
328 
329 void PndFtsHoughSpace::PrintFoundTracklets(const std::vector<PndFtsHoughTracklet> &tracklets) const
330 {
331  std::cout << tracklets.size() << " peaks found for " << GetName() << '\n';
332  if (10 < fVerbose) {
333  for (UInt_t i = 0; i < tracklets.size(); ++i) {
334  tracklets[i].Print();
335  }
336  }
337 }
338 
339 const PndFtsHit *PndFtsHoughSpace::getHitFromHS(UInt_t index) const
340 {
341  if (GetNHits() < index)
342  throwError("index too high");
343 
344  Int_t hitIndex = fHitId.at(index).GetHitId();
345  const PndFtsHit *const myHit = (PndFtsHit *)fTrackerTask->GetFtsHit(hitIndex);
346  return myHit;
347 }
348 
349 const TVector3 PndFtsHoughSpace::GetRawOrCalculatedHitPos(const PndFtsHit *const myHit) const
350 {
351  // get hit position
352  TVector3 hitPos;
353  if (kFALSE == myHit->GetSkewed()) {
354  myHit->Position(hitPos);
355  } else {
356  hitPos = CalculateHitPosFromIntersectionsWithZxTrackModel(myHit);
357  }
358  return hitPos;
359 }
360 
361 const TVector3 PndFtsHoughSpace::CalculateHitPosFromIntersectionsWithZxTrackModel(const PndFtsHit *const myHit) const
362 {
363  if (0 == fAssociatedTrackCand)
364  throwError("No track cand associated with Hough space. Cannot calculate possible hit pos. from track model.");
365 
366  const PndFtsTube *ftsTube = fTrackerTask->GetFtsTube(myHit);
367  const TVector3 wireDirection = ftsTube->GetWireDirection();
368  const TVector3 wireCenter = ftsTube->GetPosition();
369 
370  const Double_t hitZLabSys = wireCenter.Z();
371 
372  // calculate xTM according to track model at hitZLabSys
373  const Double_t xTM = fAssociatedTrackCand->getXLabSys(hitZLabSys);
374  // calculate param which is needed for xStraw = xTM
375  // xTM = wireCenter.X() + param*wireDirection.X()
376  const Double_t param = (xTM - wireCenter.X()) / wireDirection.X();
377  // calculate corresponding yStraw
378  const Double_t yStraw = wireCenter.Y() + param * wireDirection.Y();
379 
380  TVector3 crossingPosition(xTM, yStraw, hitZLabSys);
381  return crossingPosition;
382 }
383 
384 TString PndFtsHoughSpace::GetDebugOutName(TString title, Int_t param) const
385 {
386  TString debugOut = "/home/plots/";
387  debugOut += fTrackerTask->GetEventNr();
388  debugOut += title;
389  if (-1 != fRefIndex) {
390  debugOut += "-rId";
391  debugOut += fRefIndex;
392  }
393  if (-1 != param) {
394  debugOut += param;
395  }
396  debugOut += ".rtg"; // root textual graphics ;) -- actually just a macro // png does not work in this way
397  return debugOut;
398 }
399 
400 Double_t PndFtsHoughSpace::getByFromBField(Double_t hitXLabSys, Double_t hitYLabSys, Double_t hitZLabSys)
401 {
402  Double_t By = 0.;
403  // set y component of B field constant
404  if (kTRUE == fKeepBConstant) {
405  // do not take B field into account
406  By = 1.;
407  } else {
408  // Use B field information
409  Double_t po[3], BB[3];
410  po[0] = hitXLabSys; // Use magnetic field at real (not shifted) x position
411  po[1] = hitYLabSys;
412  po[2] = hitZLabSys;
413  fField->GetFieldValue(po, BB); // return value in KG (G3)
414  By = BB[1] / 10.; // By is y-component of magnetic field in Tesla
415  }
416  return By;
417 }
418 
419 #endif
void setVerbose(Int_t verbose)
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:131
std::vector< PndFtsHoughTracklet > FindAllPeaksWithTSpectrum2(const UInt_t minHeight)
Finds all peaks that satisfy the minimum height requirement minHeight.
const PndFtsTube * GetFtsTube(const PndFtsHit *const myHit) const
Returns pointer to the FTS tube corresponding to input FTS hit.
std::vector< PndFtsHoughTracklet > FindAllPeaksBinsWoMergingWithSearchWindow(const UInt_t minHeight, const Int_t vicinityLength=0)
Finds all peaks that satisfy the minimum height requirement minHeight.
Interface between PandaRoot (data input and output) and PndFtsHoughTrackFinder (implementation of PR ...
Double_t getZRefPos() const
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:130
TVector3 GetPosition() const
const PndFtsHit * GetFtsHit(UInt_t hitId) const
Returns pointer to the hit with index hitId in the FTS hit array.
std::vector< PndFtsHoughTracklet > FindAllPeaksScanPathsMergeBins(const UInt_t minHeight)
Finds all peaks that satisfy the minimum height requirement minHeight.
std::vector< PndFtsHoughTracklet > FindAllPeaksScanPathsMergeBinsCalculatingPaths(const UInt_t minHeight)
Finds all peaks that satisfy the minimum height requirement minHeight.
Int_t GetSkewed() const
Definition: PndFtsHit.h:86
unsigned int i
Definition: P4_F32vec4.h:33
UInt_t GetEventNr() const
Returns the save debug flag.
Class for saving the result of one Hough transform for FTS PR.
std::ostream & operator<<(std::ostream &os, const IdxPath &outVector)
Double_t getXLabSys(const Double_t zLabSys) const
Class for Hough space based on TH2S (for the moment). Saves the hits which enter this Hough space and...
UInt_t GetNHits() const
void Print() const
std::map< Int_t, IdxPath > HitIdxPathMap
Helper class for Hough space containing binning. Created: 09.02.2015.
PndFtsHoughSpace(const char *name=nullptr, const Int_t refIndex=-1, PndFtsHoughSpaceBinning binning=PndFtsHoughSpaceBinning(), Double_t zRefPos=0., Double_t interceptZx=0., PndFtsHoughTrackCand *associatedTrackCand=nullptr, PndFtsHoughTrackerTask *trackerTask=nullptr)
TVector3 GetWireDirection() const
Class for saving a FTS track cand. for Hough transform based FTS PR.
std::vector< PndFtsHoughTracklet > FindAllPeaksBlanko(const UInt_t minHeight)
Finds all peaks that satisfy the minimum height requirement minHeight.
ClassDef(PndFtsHoughSpace, 1)
Double_t getInterceptZx() const
std::pair< Int_t, IdxPath > HitIdxPathPair
void FillHoughSpace()
Fills the Hough space using the equation which corresponds to the name of the Hough space...
std::vector< Int_t > IdxPath