CMS 3D CMS Logo

HLTSmartSinglet.cc
Go to the documentation of this file.
1 
12 
14 
17 
19 
21 
22 //
23 // constructors and destructor
24 //
25 template<typename T>
27  inputTag_ (iConfig.template getParameter<edm::InputTag>("inputTag")),
28  inputToken_ (consumes<std::vector<T> >(inputTag_)),
29  triggerType_ (iConfig.template getParameter<int>("triggerType")),
30  cut_ (iConfig.template getParameter<std::string> ("cut" )),
31  min_N_ (iConfig.template getParameter<int> ("MinN" )),
32  select_ (cut_ )
33 {
34  LogDebug("") << "Input/tyre/cut/ncut : "
35  << inputTag_.encode() << " "
36  << triggerType_ << " "
37  << cut_<< " "
38  << min_N_ ;
39 }
40 
41 template<typename T>
43 
44 template<typename T>
45 void
49  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltCollection"));
50  desc.add<int>("triggerType",0);
51  desc.add<std::string>("cut","1>0");
52  desc.add<int>("MinN",1);
53  descriptions.add(defaultModuleLabel<HLTSmartSinglet<T>>(), desc);
54 }
55 
56 //
57 // member functions
58 //
59 
60 // ------------ method called to produce the data ------------
61 template<typename T>
62 bool
64 {
65  using namespace std;
66  using namespace edm;
67  using namespace reco;
68  using namespace trigger;
69 
70  typedef vector<T> TCollection;
71  typedef Ref<TCollection> TRef;
72 
73  // All HLT filters must create and fill an HLT filter object,
74  // recording any reconstructed physics objects satisfying (or not)
75  // this HLT filter, and place it in the Event.
76 
77  // The filter object
78  if (saveTags()) filterproduct.addCollectionTag(inputTag_);
79 
80  // Ref to Candidate object to be recorded in filter object
81  TRef ref;
82 
83  // get hold of collection of objects
85  iEvent.getByToken (inputToken_,objects);
86 
87  // look at all objects, check cuts and add to filter object
88  int n(0);
89  typename TCollection::const_iterator i ( objects->begin() );
90  for (; i!=objects->end(); i++) {
91  if (select_(*i)) {
92  n++;
93  ref=TRef(objects,distance(objects->begin(),i));
94  filterproduct.addObject(triggerType_,ref);
95  }
96  }
97 
98  // filter decision
99  bool accept(n>=min_N_);
100 
101  return accept;
102 }
#define LogDebug(id)
std::string defaultModuleLabel()
HLTSmartSinglet(const edm::ParameterSet &)
virtual bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
edm::InputTag inputTag_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
std::string encode() const
Definition: InputTag.cc:166
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
int iEvent
Definition: GenABIO.cc:230
std::string cut_
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:520
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
StringCutObjectSelector< T, true > select_
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
edm::EDGetTokenT< std::vector< T > > inputToken_
HLT enums.
long double T
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)