CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MinMETProducerT.h
Go to the documentation of this file.
1 #ifndef RecoMET_METProducers_MinMETProducerT_h
2 #define RecoMET_METProducers_MinMETProducerT_h
3 
24 
25 #include <vector>
26 
27 template <typename T>
29 {
30  typedef std::vector<T> METCollection;
31 
32  public:
33 
34  explicit MinMETProducerT(const edm::ParameterSet& cfg) : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
35  {
36  src_ = cfg.getParameter<vInputTag>("src");
37 
38  for ( vInputTag::const_iterator src_i = src_.begin();
39  src_i != src_.end(); ++src_i ) {
40  src_token_.push_back( consumes<METCollection>(*src_i) );
41  }
42  produces<METCollection>();
43  }
45 
46  private:
47 
48  void produce(edm::Event& evt, const edm::EventSetup& es) override
49  {
50  std::auto_ptr<METCollection> outputMETs(new METCollection());
51 
52  // check that all MET collections given as input have the same number of entries
53  int numMEtObjects = -1;
54  // for ( vInputTag::const_iterator src_i = src_.begin();
55  // src_i != src_.end(); ++src_i ) {
56  for ( typename vInputToken::const_iterator src_i = src_token_.begin();
57  src_i != src_token_.end(); ++src_i ) {
59  // evt.getByLabel(*src_i, inputMETs);
60  evt.getByToken(*src_i, inputMETs);
61  if ( numMEtObjects == -1 ) numMEtObjects = inputMETs->size();
62  else if ( numMEtObjects != (int)inputMETs->size() )
63  throw cms::Exception("MinMETProducer::produce")
64  << "Mismatch in number of input MET objects !!\n";
65  }
66 
67  for ( int iMEtObject = 0; iMEtObject < numMEtObjects; ++iMEtObject ) {
68  const T* minMET = 0;
69  // for ( vInputTag::const_iterator src_i = src_.begin();
70  // src_i != src_.end(); ++src_i ) {
71  for ( typename vInputToken::const_iterator src_i = src_token_.begin();
72  src_i != src_token_.end(); ++src_i ) {
74  // evt.getByLabel(*src_i, inputMETs);
75  evt.getByToken(*src_i, inputMETs);
76  const T& inputMET = inputMETs->at(iMEtObject);
77  if ( minMET == 0 || inputMET.pt() < minMET->pt() ) minMET = &inputMET;
78  }
79  assert(minMET);
80  outputMETs->push_back(T(*minMET));
81  }
82 
83  evt.put(outputMETs);
84  }
85 
87 
88  typedef std::vector<edm::InputTag> vInputTag;
90  typedef std::vector<edm::EDGetTokenT<METCollection> > vInputToken;
92 };
93 
96 
97 namespace reco
98 {
101 }
102 
103 #endif
104 
105 
106 
107 
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:259
std::vector< edm::InputTag > vInputTag
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
assert(m_qm.get())
vInputToken src_token_
MinMETProducerT< reco::PFMET > MinPFMETProducer
std::vector< T > METCollection
std::string moduleLabel_
std::vector< edm::EDGetTokenT< METCollection > > vInputToken
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
void produce(edm::Event &evt, const edm::EventSetup &es) override
MinMETProducerT< reco::CaloMET > MinCaloMETProducer
MinMETProducerT(const edm::ParameterSet &cfg)
long double T