ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
AtPRA.h
Go to the documentation of this file.
1 #ifndef AtPRA_H
2 #define AtPRA_H
3 
4 #include "AtTrack.h" // for AtTrack
5 #include "AtTrackTransformer.h"
6 
7 #include <Rtypes.h> // for Double_t, Float_t, Int_t, THashConsistencyHolder
8 #include <TObject.h> // for TObject
9 
10 #include <algorithm> // for max
11 #include <memory>
12 #include <type_traits> // for false_type, is_signed, true_type
13 #include <vector> // for vector
14 
15 class AtDigiPar;
16 class AtEvent;
17 class AtHit;
18 class AtPatternEvent;
19 class TBuffer;
20 class TClass;
21 class TMemberInspector;
22 
23 namespace AtPATTERN {
24 
32 class AtPRA : public TObject {
33 protected:
34  std::vector<AtTrack> fTrackCand; //< Candidate tracks
35 
37 
38  Int_t fMaxHits{5000};
39  Int_t fMinHits{0};
40  Float_t fMeanDistance{1e9};
41  Int_t fKNN{5}; //<! Number of nearest neighbors kNN
42  Double_t fStdDevMulkNN{0}; //<! Std dev multiplier for kNN
43  Double_t fkNNDist{10}; //<! Distance threshold for outlier rejection in kNN
44 
45  Bool_t kSetPrunning{false}; //<<! Enable prunning of tracks
46 
47  std::unique_ptr<AtTools::AtTrackTransformer> fTrackTransformer{std::make_unique<AtTools::AtTrackTransformer>()};
48  Double_t fClusterRadius{0}; //<! Radius of hit clusters
49  Double_t fClusterDistance{0}; //<! Distance between hit clusters
50 
51 public:
52  virtual ~AtPRA() = default;
53 
54  // Getters
55  virtual std::vector<AtTrack> GetTrackCand() const { return fTrackCand; }
56 
57  // Setters
58  void SetMaxHits(Int_t maxHits) { fMaxHits = maxHits; }
59  void SetMinHits(Int_t minHits) { fMinHits = minHits; }
60  void SetMeanDistance(Float_t meanDistance) { fMeanDistance = meanDistance; }
61  void SetkNN(Double_t knn) { fKNN = knn; }
62  void SetStdDevMulkNN(Double_t stdDevMul) { fStdDevMulkNN = stdDevMul; }
63  void SetkNNDist(Double_t dist) { fkNNDist = dist; }
64  void SetPrunning() { kSetPrunning = kTRUE; }
65  void SetClusterRadius(Double_t clusterRadius) { fClusterRadius = clusterRadius; }
66  void SetClusterDistance(Double_t clusterDistance) { fClusterDistance = clusterDistance; }
67 
68  virtual std::unique_ptr<AtPatternEvent> FindTracks(AtEvent &event) = 0;
69 
70  void PruneTrack(AtTrack &track);
71  bool kNN(const std::vector<std::unique_ptr<AtHit>> &hits, AtHit &hit, int k);
72 
73 protected:
74  // Functions that need to be moved to another class. They assume a curved track
75  /*
76  * Takes track and sets fGeo... parameters using SampleConsensus. I think these are then used
77  * as initial guesses for GenFit. Probably, this should be moved into the fitting classes instead of
78  * here. At the very least, it needs to be in a subclass that only deals with curved tracks, it
79  * would make no sense to apply this function to straight tracks.
80  */
82 
83  template <typename T>
84  inline constexpr int GetSign(T num, std::true_type is_signed)
85  {
86  return (T(0) < num) - (num < T(0));
87  }
88  template <typename T>
89  inline constexpr int GetSign(T num, std::false_type is_signed)
90  {
91  return (T(0) < num);
92  }
93  template <typename T>
94  inline constexpr int GetSign(T num)
95  {
96  return GetSign(num, std::is_signed<T>());
97  }
98 
99  ClassDef(AtPRA, 1)
100 };
101 
102 } // namespace AtPATTERN
103 
104 Double_t fitf(Double_t *x, Double_t *par);
105 
106 #endif
AtPATTERN::AtPRA::kNN
bool kNN(const std::vector< std::unique_ptr< AtHit >> &hits, AtHit &hit, int k)
Definition: AtPRA.cxx:228
AtPatternEvent
Definition: AtPatternEvent.h:19
AtPATTERN::AtPRA::fPar
AtDigiPar * fPar
parameter container
Definition: AtPRA.h:36
AtPATTERN::AtPRA::fkNNDist
Double_t fkNNDist
Definition: AtPRA.h:43
AtPATTERN::AtPRA::SetMinHits
void SetMinHits(Int_t minHits)
Definition: AtPRA.h:59
fitf
Double_t fitf(Double_t *x, Double_t *par)
Definition: AtPRA.cxx:197
AtPATTERN::AtPRA::GetSign
constexpr int GetSign(T num, std::false_type is_signed)
Definition: AtPRA.h:89
AtPATTERN
Definition: AtEventDrawTaskS800.h:41
AtPATTERN::AtPRA::kSetPrunning
Bool_t kSetPrunning
Definition: AtPRA.h:45
AtPATTERN::AtPRA::SetkNN
void SetkNN(Double_t knn)
Definition: AtPRA.h:61
AtPATTERN::AtPRA::SetPrunning
void SetPrunning()
Definition: AtPRA.h:64
AtPATTERN::AtPRA::SetMaxHits
void SetMaxHits(Int_t maxHits)
Definition: AtPRA.h:58
AtPATTERN::AtPRA::SetkNNDist
void SetkNNDist(Double_t dist)
Definition: AtPRA.h:63
AtPATTERN::AtPRA::fStdDevMulkNN
Double_t fStdDevMulkNN
Definition: AtPRA.h:42
AtPATTERN::AtPRA::fClusterDistance
Double_t fClusterDistance
Definition: AtPRA.h:49
AtPATTERN::AtPRA
Find patterns in hit clouds.
Definition: AtPRA.h:32
AtEvent
Definition: AtEvent.h:22
AtPATTERN::AtPRA::SetMeanDistance
void SetMeanDistance(Float_t meanDistance)
Definition: AtPRA.h:60
AtPATTERN::AtPRA::fKNN
Int_t fKNN
Definition: AtPRA.h:41
AtPATTERN::AtPRA::fClusterRadius
Double_t fClusterRadius
Definition: AtPRA.h:48
AtPATTERN::AtPRA::GetSign
constexpr int GetSign(T num, std::true_type is_signed)
Definition: AtPRA.h:84
AtPATTERN::AtPRA::SetStdDevMulkNN
void SetStdDevMulkNN(Double_t stdDevMul)
Definition: AtPRA.h:62
AtPATTERN::AtPRA::SetClusterRadius
void SetClusterRadius(Double_t clusterRadius)
Definition: AtPRA.h:65
AtPATTERN::AtPRA::~AtPRA
virtual ~AtPRA()=default
AtPATTERN::AtPRA::fTrackCand
std::vector< AtTrack > fTrackCand
Definition: AtPRA.h:34
AtPATTERN::AtPRA::fTrackTransformer
std::unique_ptr< AtTools::AtTrackTransformer > fTrackTransformer
Definition: AtPRA.h:47
AtTrack
Definition: AtTrack.h:25
AtPATTERN::AtPRA::SetClusterDistance
void SetClusterDistance(Double_t clusterDistance)
Definition: AtPRA.h:66
AtPATTERN::AtPRA::SetTrackInitialParameters
void SetTrackInitialParameters(AtTrack &track)
Set initial parameters for HC.
Definition: AtPRA.cxx:35
AtPATTERN::AtPRA::fMeanDistance
Float_t fMeanDistance
Definition: AtPRA.h:40
AtDigiPar
Definition: AtDigiPar.h:14
AtTrack.h
AtPATTERN::AtPRA::GetTrackCand
virtual std::vector< AtTrack > GetTrackCand() const
Definition: AtPRA.h:55
AtPATTERN::AtPRA::fMinHits
Int_t fMinHits
Definition: AtPRA.h:39
AtPATTERN::AtPRA::PruneTrack
void PruneTrack(AtTrack &track)
Definition: AtPRA.cxx:203
AtTrackTransformer.h
AtPATTERN::AtPRA::fMaxHits
Int_t fMaxHits
Definition: AtPRA.h:38
AtPATTERN::AtPRA::FindTracks
virtual std::unique_ptr< AtPatternEvent > FindTracks(AtEvent &event)=0
AtPATTERN::AtPRA::GetSign
constexpr int GetSign(T num)
Definition: AtPRA.h:94
AtHit
Point in space with charge.
Definition: AtHit.h:27