CMS 3D CMS Logo

GroupedCkfTrajectoryBuilder.h
Go to the documentation of this file.
1 #ifndef GroupedCkfTrajectoryBuilder_H
2 #define GroupedCkfTrajectoryBuilder_H
3 
7 
10 
12 
13 #include <vector>
14 
15 
17 
18 
25 
26  public:
29 
32 
34  // virtual void setEvent(const edm::Event& event) const;
35 
37  TrajectoryContainer trajectories(const TrajectorySeed&) const override;
38 
40  void trajectories(const TrajectorySeed&, TrajectoryContainer &ret) const override;
41 
44 
46  void trajectories(const TrajectorySeed&, TrajectoryContainer &ret, const TrackingRegion&) const;
47 
49  // also new interface returning the start Trajectory...
52  unsigned int& nCandPerSeed,
53  const TrajectoryFilter*) const override;
54 
55 
56 
64  TrajectoryContainer& result) const override ;
65 
66  // same as above using the precomputed startingTraj..
67  void rebuildTrajectories(TempTrajectory const & startingTraj, const TrajectorySeed&,
68  TrajectoryContainer& result) const override ;
69 
70 
71  // Access to lower level components
72  const TrajectoryStateUpdator& updator() const {return *theUpdator;}
74 
75  // PropagationDirection direction() const {return theDirection;}
76 
78  double chiSquareCut() {return theChiSquareCut;}
79 
81  int maxCand() {return theMaxCand;}
82 
83 
85  float lostHitPenalty() {return theLostHitPenalty;}
86 
87  // /** Tells whether an intermediary cleaning stage should take place during TB. */
88  // bool intermediateCleaning() {return theIntermediateCleaning;}
89 
91  double ptCut() {return theptCut;}
92 
94  double mass() {return theMass;}
95 
96 protected:
97  void setEvent_(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
98 
99  virtual void analyseSeed(const TrajectorySeed& seed) const{}
100 
101  virtual void analyseMeasurements( const std::vector<TM>& meas,
102  const Trajectory& traj) const{}
103  virtual void analyseResult( const TrajectoryContainer& result) const {}
104 
105 private :
106 // /// no copy constructor
107 // GroupedCkfTrajectoryBuilder (const GroupedCkfTrajectoryBuilder&) = default;
108 //
109 // /// no assignment operator
110 // GroupedCkfTrajectoryBuilder& operator= (const GroupedCkfTrajectoryBuilder&) dso_internal;
111 
112 
113  inline bool tkxor(bool a, bool b) const dso_internal {return (a||b) && !(a&&b);}
114  // to be ported later
115 
116  bool advanceOneLayer( const TrajectorySeed& seed,
117  TempTrajectory& traj,
118  const TrajectoryFilter* regionalCondition,
119  const Propagator* propagator,
120  bool inOut,
121  TempTrajectoryContainer& newCand,
122  TempTrajectoryContainer& result) const dso_internal;
123 
124  unsigned int groupedLimitedCandidates( const TrajectorySeed& seed,
125  TempTrajectory const& startingTraj,
126  const TrajectoryFilter* regionalCondition,
127  const Propagator* propagator,
128  bool inOut,
129  TempTrajectoryContainer& result) const dso_internal;
130 
132  void rebuildSeedingRegion (const TrajectorySeed&seed,
133  TempTrajectory const& startingTraj,
134  TempTrajectoryContainer& result) const dso_internal;
135 
136  //** try to find additional hits in seeding region for a candidate
137  //* (returns number of trajectories added) *
138  int rebuildSeedingRegion (const TrajectorySeed&seed,
139  const std::vector<const TrackingRecHit*>& seedHits,
140  TempTrajectory& candidate,
141  TempTrajectoryContainer& result) const dso_internal;
142 
143  // ** Backward fit of trajectory candidate except seed. Fit result is returned. invalid if fit failed
144  // * remaining hits are returned remainingHits.
145  TempTrajectory backwardFit (TempTrajectory& candidate, unsigned int nSeed,
146  const TrajectoryFitter& fitter,
147  std::vector<const TrackingRecHit*>& remainingHits) const dso_internal;
148 
150  bool verifyHits (TempTrajectory::DataContainer::const_iterator rbegin,
151  size_t maxDepth,
152  const std::vector<const TrackingRecHit*>& hits) const dso_internal;
153 
155  void groupedIntermediaryClean(TempTrajectoryContainer& theTrajectories) const dso_internal;
156 
157 
160  if ( dir==alongMomentum ) return oppositeToMomentum;
161  if ( dir==oppositeToMomentum ) return alongMomentum;
162  return dir;
163  }
164 
165 
166 private:
168 
169  // typedef deque< const TrajectoryFilter*> StopCondContainer;
170  // StopCondContainer theStopConditions;
171 
174  double theptCut;
176  double theMass;
187 
188  bool theLockHits;
204 
206 
207 // mutable TempTrajectoryContainer work_; // Better here than alloc every time
208  enum work_MaxSize_Size_ { work_MaxSize_ = 50 }; // if it grows above this number, it is forced to resize to half this amount when cleared
209 };
210 
211 #endif
virtual void analyseMeasurements(const std::vector< TM > &meas, const Trajectory &traj) const
virtual void rebuildSeedingRegion(const TrajectorySeed &, TrajectoryContainer &result) const
~GroupedCkfTrajectoryBuilder() override
destructor
const TrajectoryStateUpdator * theUpdator
virtual void analyseSeed(const TrajectorySeed &seed) const
double mass()
Mass hypothesis used for propagation.
PropagationDirection
virtual TempTrajectory buildTrajectories(const TrajectorySeed &seed, TrajectoryContainer &ret, unsigned int &nCandPerSeed, const TrajectoryFilter *) const
virtual TrajectoryContainer trajectories(const TrajectorySeed &) const =0
virtual void setEvent_(const edm::Event &iEvent, const edm::EventSetup &iSetup)=0
int iEvent
Definition: GenABIO.cc:230
virtual void analyseResult(const TrajectoryContainer &result) const
static PropagationDirection oppositeDirection(PropagationDirection dir)
change of propagation direction
virtual void rebuildTrajectories(TempTrajectory const &startingTraj, const TrajectorySeed &seed, TrajectoryContainer &result) const
const Chi2MeasurementEstimatorBase & estimator() const
std::vector< TempTrajectory > TempTrajectoryContainer
double b
Definition: hdecay.h:120
const TrajectoryStateUpdator & updator() const
#define dso_internal
double a
Definition: hdecay.h:121
bool tkxor(bool a, bool b) const
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< Trajectory > TrajectoryContainer
const Chi2MeasurementEstimatorBase * theEstimator