CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CkfBaseTrajectoryFilter.h
Go to the documentation of this file.
1 #ifndef CkfBaseTrajectoryFilter_H
2 #define CkfBaseTrajectoryFilter_H
3 
5 
14 
15 
17 public:
18 
20  //define the filters by default in the BaseCkfTrajectoryBuilder
29  {}
30 
31  void setEvent(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
32  theChargeSignificanceTrajectoryFilter->setEvent(iEvent, iSetup);
33  theMaxLostHitsTrajectoryFilter->setEvent(iEvent, iSetup);
34  theMaxConsecLostHitsTrajectoryFilter->setEvent(iEvent, iSetup);
35  theMinPtTrajectoryFilter->setEvent(iEvent, iSetup);
36  theMaxHitsTrajectoryFilter->setEvent(iEvent, iSetup);
37  theMinHitsTrajectoryFilter->setEvent(iEvent, iSetup);
38  theLostHitsFractionTrajectoryFilter->setEvent(iEvent, iSetup);
39  theLooperTrajectoryFilter->setEvent(iEvent, iSetup);
40  }
41 
42  virtual bool qualityFilter( const Trajectory& traj) const {return QF<Trajectory>(traj);}
43  virtual bool qualityFilter( const TempTrajectory& traj) const {return QF<TempTrajectory>(traj);}
44 
45  virtual bool toBeContinued( Trajectory& traj) const {return TBC<Trajectory>(traj);}
46  virtual bool toBeContinued( TempTrajectory& traj) const {return TBC<TempTrajectory>(traj);}
47 
48  virtual std::string name() const { return "CkfBaseTrajectoryFilter";}
49 
50 protected:
51 
52  template <class T> bool QF(const T& traj) const{
53  if (!theChargeSignificanceTrajectoryFilter->qualityFilter(traj)) return false;
54  if (!theMinHitsTrajectoryFilter->qualityFilter(traj)) return false;
55  if (!theMinPtTrajectoryFilter->qualityFilter(traj)) return false;
56  if (!theLooperTrajectoryFilter->qualityFilter(traj)) return false;
57  return true;}
58 
59  template <class T> bool TBC(T& traj) const{
60  if (!theMaxHitsTrajectoryFilter->toBeContinued(traj)) return false;
61  if (!theMaxLostHitsTrajectoryFilter->toBeContinued(traj)) return false;
62  if (!theMaxConsecLostHitsTrajectoryFilter->toBeContinued(traj)) return false;
63  if (!theLostHitsFractionTrajectoryFilter->toBeContinued(traj)) return false;
64  if (!theMinPtTrajectoryFilter->toBeContinued(traj)) return false;
65  if (!theChargeSignificanceTrajectoryFilter->toBeContinued(traj)) return false;
66  if (!theLooperTrajectoryFilter->toBeContinued(traj)) return false;
67  return true;}
68 
69 
70 
71  std::unique_ptr<ChargeSignificanceTrajectoryFilter> theChargeSignificanceTrajectoryFilter;
72  std::unique_ptr<MaxConsecLostHitsTrajectoryFilter> theMaxConsecLostHitsTrajectoryFilter;
73  std::unique_ptr<MaxHitsTrajectoryFilter> theMaxHitsTrajectoryFilter;
74  std::unique_ptr<MaxLostHitsTrajectoryFilter> theMaxLostHitsTrajectoryFilter;
75  std::unique_ptr<LostHitsFractionTrajectoryFilter> theLostHitsFractionTrajectoryFilter;
76  std::unique_ptr<MinHitsTrajectoryFilter> theMinHitsTrajectoryFilter;
77  std::unique_ptr<MinPtTrajectoryFilter> theMinPtTrajectoryFilter;
78  std::unique_ptr<LooperTrajectoryFilter> theLooperTrajectoryFilter;
79 };
80 
81 #endif
std::unique_ptr< LooperTrajectoryFilter > theLooperTrajectoryFilter
CkfBaseTrajectoryFilter(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
int iEvent
Definition: GenABIO.cc:230
std::unique_ptr< MaxLostHitsTrajectoryFilter > theMaxLostHitsTrajectoryFilter
virtual bool toBeContinued(Trajectory &traj) const
std::unique_ptr< ChargeSignificanceTrajectoryFilter > theChargeSignificanceTrajectoryFilter
bool QF(const T &traj) const
std::unique_ptr< MinPtTrajectoryFilter > theMinPtTrajectoryFilter
std::unique_ptr< MinHitsTrajectoryFilter > theMinHitsTrajectoryFilter
std::unique_ptr< LostHitsFractionTrajectoryFilter > theLostHitsFractionTrajectoryFilter
virtual std::string name() const
std::unique_ptr< MaxConsecLostHitsTrajectoryFilter > theMaxConsecLostHitsTrajectoryFilter
void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
std::unique_ptr< MaxHitsTrajectoryFilter > theMaxHitsTrajectoryFilter
virtual bool qualityFilter(const Trajectory &traj) const
virtual bool toBeContinued(TempTrajectory &traj) const
long double T
virtual bool qualityFilter(const TempTrajectory &traj) const