PandaRoot
PndHoughData.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 
14 // PndHoughData
15 // Data structure for fast access to STT data
17 
28 /*
29 This class is similar to the data structure of the CellTrackFinder.
30 */
31 
32 #ifndef PndHoughData_H_
33 #define PndHoughData_H_
34 
35 #include <vector>
36 #include <map>
37 
38 #include "PndSttGeometryMap.h"
39 #include "PndSttStrawMap.h"
40 #include "FairHit.h"
41 #include "PndHoughSpace.h"
42 
43 class PndHoughData {
44 
45  public:
46  PndHoughData(TClonesArray *fTubeArray);
47 
48  virtual ~PndHoughData()
49  {
50  delete fGeometryMap;
51  delete fPndHoughSpace;
52  for (size_t i = 0; i < fHits.size(); ++i)
53  delete fHits.at(i);
54  }
55 
57  void AddHits(TClonesArray *hits, TString branchName);
59  void AddHit(FairLink link);
65  void Init(int NBins1 = 450, int NBins2 = 450)
66  {
67  fNbins1 = NBins1;
68  fNbins2 = NBins2;
69  fPndHoughSpace->SetNBins1(fNbins1);
70  fPndHoughSpace->SetNBins2(fNbins2);
71  fPndHoughSpace->Init();
72  }
74  void clear()
75  {
76  fHits.clear();
77  fLinks.clear();
78  fSTTLinks.clear();
79  fGEMLinks.clear();
80  fhittedTubes.clear();
81  fMVDHits.clear();
82 
83  fMapFairLinktoFairHit.clear();
84  fMapFairLinktoIsochrone.clear();
85  fMapFairLinktoIsochroneError.clear();
86  fMapFairLinktoTubeId.clear();
87  fMapTubetoHit.clear();
88  fMapTubetoLink.clear();
89 
90  fGEMNeighbors.clear();
91  fMapNumberhittedSTTNeighbors.clear();
92 
93  fPndHoughSpace->clear();
94  }
96  std::vector<FairHit *> GetHits() const { return fHits; }
98  std::vector<FairLink> GetLinks() const { return fLinks; }
100  PndSttGeometryMap *GetGeometryMap() const { return fGeometryMap; }
102  PndSttStrawMap *GetStrawMap() const { return fStrawMap; }
104  std::map<FairLink, FairHit *> GetMapFairLinktoFairHit() const { return fMapFairLinktoFairHit; }
106  std::map<FairLink, Double_t> GetMapFairLinktoIsochrone() const { return fMapFairLinktoIsochrone; }
108  std::map<FairLink, Double_t> GetMapFairLinktoIsochroneError() const { return fMapFairLinktoIsochroneError; }
110  std::map<FairLink, Int_t> GetMapFairLinktoTubeId() const { return fMapFairLinktoTubeId; }
112  std::map<int, FairHit *> GetMapTubetoHit() const { return fMapTubetoHit; }
114  std::map<int, FairLink> GetMapTubetoLink() const { return fMapTubetoLink; }
116  std::vector<int> GetHittedTubes() { return fhittedTubes; };
118  std::vector<FairLink> GetGEMNeighbors(FairLink link) { return fGEMNeighbors[link]; };
120  int GetNumOfGEMNeighbors(FairLink link) { return fGEMNeighbors[link].size(); };
122  double GetNumOfSTTNeighbors(FairLink link) { return fMapNumberhittedSTTNeighbors[link]; };
124  PndHoughSpace *GetHoughSpace() { return fPndHoughSpace; };
126  std::vector<FairLink> GetMVDHits() { return fMVDHits; };
128  std::vector<FairLink> GetGEMHits() { return fGEMLinks; };
130  std::vector<FairLink> GetSTTHits() { return fSTTLinks; };
131 
132  private:
133  PndSttStrawMap *fStrawMap = nullptr; // for getting more information about the tubes
134  PndSttGeometryMap *fGeometryMap = nullptr; // for initializing the neighbors of each tube
135  PndHoughSpace *fPndHoughSpace = nullptr;
136 
137  std::map<FairLink, FairHit *> fMapFairLinktoFairHit;
138  std::map<FairLink, Double_t> fMapFairLinktoIsochrone;
139  std::map<FairLink, Double_t> fMapFairLinktoIsochroneError;
140  std::map<FairLink, Int_t> fMapFairLinktoTubeId;
141  std::map<FairLink, std::vector<FairLink>> fGEMNeighbors;
142  std::map<FairLink, double> fMapNumberhittedSTTNeighbors;
143  std::map<int, FairHit *> fMapTubetoHit; // maps hitted tube to FairHit
144  std::map<int, FairLink> fMapTubetoLink; // maps hitted tube to FairLink
145 
146  std::vector<FairHit *> fHits;
147  std::vector<FairLink> fGEMLinks;
148  std::vector<FairLink> fSTTLinks;
149  std::vector<FairLink> fLinks; // vector with selected hits of an event
150  std::vector<int> fhittedTubes;
151  std::vector<FairLink> fMVDHits;
152 
153  int fNbins1;
154  int fNbins2;
155 
156  ClassDef(PndHoughData, 1);
157 };
158 
159 #endif /* PndHoughData_H_ */
void AddHit(FairLink link)
Adds a sigle hit to the used data structure for faster access. Data are stored in different maps...
void AddHits(TClonesArray *hits, TString branchName)
Adds hits to the used data structure for faster access. Data are stored in different maps...
std::map< FairLink, Int_t > GetMapFairLinktoTubeId() const
Returns the map linking FairLinks to corresponding tubeId.
Definition: PndHoughData.h:110
std::map< FairLink, Double_t > GetMapFairLinktoIsochroneError() const
Returns the map linking FairLinks to isochrone radius errors.
Definition: PndHoughData.h:108
virtual ~PndHoughData()
Definition: PndHoughData.h:48
void CreateSTTNeighborhoodData()
For all STT hits all other STT neighbors are counted and stored in a map.
std::vector< FairLink > GetSTTHits()
Returns all STT Hits.
Definition: PndHoughData.h:130
std::vector< FairHit * > GetHits() const
Returns a vector of all hits in the event.
Definition: PndHoughData.h:96
void Init()
Initializes the Hough space.
std::map< FairLink, FairHit * > GetMapFairLinktoFairHit() const
Returns the map linking FairLinks to FairHits.
Definition: PndHoughData.h:104
unsigned int i
Definition: P4_F32vec4.h:33
PndSttStrawMap * GetStrawMap() const
Returns the straw map for a fast access.
Definition: PndHoughData.h:102
PndHoughData(TClonesArray *fTubeArray)
int GetNumOfGEMNeighbors(FairLink link)
Returns the number of GEM hits in a distance d < 1.5 cm.
Definition: PndHoughData.h:120
virtual void SetNBins1(double n)
Sets the number of bins in x direction of the Hough space.
Definition: PndHoughSpace.h:47
PndSttGeometryMap * GetGeometryMap() const
Returns the geometry map for a fast access.
Definition: PndHoughData.h:100
std::map< int, FairHit * > GetMapTubetoHit() const
Returns the map linking a tube to a FairHit.
Definition: PndHoughData.h:112
std::map< FairLink, Double_t > GetMapFairLinktoIsochrone() const
Returns the map linking FairLinks to isochrone radii.
Definition: PndHoughData.h:106
PndHoughSpace * GetHoughSpace()
Returns the Hough space.
Definition: PndHoughData.h:124
void Init(int NBins1=450, int NBins2=450)
Initializes the Hough space.
Definition: PndHoughData.h:65
void clear()
Clears all data maps.
Definition: PndHoughData.h:74
std::vector< int > GetHittedTubes()
Returns a vector of all hitted tubes.
Definition: PndHoughData.h:116
double GetNumOfSTTNeighbors(FairLink link)
Returns the number of neighbored STT hits. The return value is a double (not a int) to avoid a necess...
Definition: PndHoughData.h:122
std::vector< FairLink > GetMVDHits()
Returns all MVD Hits.
Definition: PndHoughData.h:126
virtual void SetNBins2(double n)
Sets the number of bins in y direction of the Hough space.
Definition: PndHoughSpace.h:49
std::vector< FairLink > GetGEMNeighbors(FairLink link)
Returns a vector of all GEM hits in a distance d < 1.5 cm from a certain GEM hit. ...
Definition: PndHoughData.h:118
std::vector< FairLink > GetLinks() const
Returns a vector of all FairLinks in the event.
Definition: PndHoughData.h:98
void CreateGEMNeighborhoodData()
For all GEM hits all other GEM hits in a certain distance (hier 1.5 cm) are counted and stored as nei...
std::map< int, FairLink > GetMapTubetoLink() const
Returns the map linking a tube to a FairLink.
Definition: PndHoughData.h:114
std::vector< FairLink > GetGEMHits()
Returns all GEM Hits.
Definition: PndHoughData.h:128
void clear()
Clears the Hough space.
Definition: PndHoughSpace.h:61