CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CompositeTrajectoryFilter.h
Go to the documentation of this file.
1 #ifndef CompositeTrajectoryFilter_H
2 #define CompositeTrajectoryFilter_H
3 
8 
16 public:
17 
18  explicit CompositeTrajectoryFilter(){filters.clear();}
20  {
21  //look for VPSet of filters
22  std::vector<edm::ParameterSet> vpset=pset.getParameter<std::vector<edm::ParameterSet> >("filters");
23  for (unsigned int i=0;i!= vpset.size();i++)
24  {filters.emplace_back(TrajectoryFilterFactory::get()->create(vpset[i].getParameter<std::string>("ComponentType"),
25  vpset[i], iC));}
26  }
27 
29 
30  void setEvent(const edm::Event& iEvent, const edm::EventSetup& iSetup) override {
31  for(auto& f: filters) {
32  f->setEvent(iEvent, iSetup);
33  }
34  }
35 
36  virtual bool qualityFilter( const Trajectory& traj) const { return QF<Trajectory>(traj);}
37  virtual bool qualityFilter( const TempTrajectory& traj) const { return QF<TempTrajectory>(traj);}
38 
39  virtual bool toBeContinued( Trajectory& traj) const { return TBC<Trajectory>(traj);}
40  virtual bool toBeContinued( TempTrajectory& traj) const { return TBC<TempTrajectory>(traj);}
41 
42  virtual std::string name() const { std::string rname="CompositeTrajectoryFilter";
43  unsigned int i=0;
44  unsigned int n=filters.size();
45  for (;i<n;i++){ rname+="_"+filters[i]->name();}
46  return rname;
47  }
48 
49 protected:
50  template <class T> bool TBC( T& traj)const{
51  unsigned int i=0;
52  unsigned int n=filters.size();
53  for (;i<n;i++){ if (!filters[i]->toBeContinued(traj)) return false; }
54  return true;}
55 
56  template <class T> bool QF(const T& traj)const{
57  unsigned int i=0;
58  unsigned int n=filters.size();
59  for (;i<n;i++){ if (!filters[i]->qualityFilter(traj)) return false; }
60  return true;}
61  protected:
62  std::vector<std::unique_ptr<TrajectoryFilter> > filters;
63 
64 };
65 
66 #endif
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual bool qualityFilter(const TempTrajectory &traj) const
CompositeTrajectoryFilter(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
virtual bool toBeContinued(TempTrajectory &traj) const
int iEvent
Definition: GenABIO.cc:230
std::vector< std::unique_ptr< TrajectoryFilter > > filters
void setEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
bool QF(const T &traj) const
double f[11][100]
virtual std::string name() const
virtual bool toBeContinued(Trajectory &traj) const
long double T
virtual bool qualityFilter(const Trajectory &traj) const
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:55