CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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  typedef std::vector<T> METCollection;
30 
31 public:
33  : moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
34  src_ = cfg.getParameter<vInputTag>("src");
35 
36  for (vInputTag::const_iterator src_i = src_.begin(); src_i != src_.end(); ++src_i) {
37  src_token_.push_back(consumes<METCollection>(*src_i));
38  }
39  produces<METCollection>();
40  }
41  ~MinMETProducerT() override {}
42 
43 private:
44  void produce(edm::Event& evt, const edm::EventSetup& es) override {
45  auto outputMETs = std::make_unique<METCollection>();
46 
47  // check that all MET collections given as input have the same number of entries
48  int numMEtObjects = -1;
49  // for ( vInputTag::const_iterator src_i = src_.begin();
50  // src_i != src_.end(); ++src_i ) {
51  for (typename vInputToken::const_iterator src_i = src_token_.begin(); src_i != src_token_.end(); ++src_i) {
53  // evt.getByLabel(*src_i, inputMETs);
54  evt.getByToken(*src_i, inputMETs);
55  if (numMEtObjects == -1)
56  numMEtObjects = inputMETs->size();
57  else if (numMEtObjects != (int)inputMETs->size())
58  throw cms::Exception("MinMETProducer::produce") << "Mismatch in number of input MET objects !!\n";
59  }
60 
61  for (int iMEtObject = 0; iMEtObject < numMEtObjects; ++iMEtObject) {
62  const T* minMET = nullptr;
63  // for ( vInputTag::const_iterator src_i = src_.begin();
64  // src_i != src_.end(); ++src_i ) {
65  for (typename vInputToken::const_iterator src_i = src_token_.begin(); src_i != src_token_.end(); ++src_i) {
67  // evt.getByLabel(*src_i, inputMETs);
68  evt.getByToken(*src_i, inputMETs);
69  const T& inputMET = inputMETs->at(iMEtObject);
70  if (minMET == nullptr || inputMET.pt() < minMET->pt())
71  minMET = &inputMET;
72  }
73  assert(minMET);
74  outputMETs->push_back(T(*minMET));
75  }
76 
77  evt.put(std::move(outputMETs));
78  }
79 
81 
82  typedef std::vector<edm::InputTag> vInputTag;
84  typedef std::vector<edm::EDGetTokenT<METCollection> > vInputToken;
86 };
87 
90 
91 namespace reco {
94 } // namespace reco
95 
96 #endif
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple cfg
Definition: looper.py:296
std::vector< edm::InputTag > vInputTag
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
vInputToken src_token_
MinMETProducerT< reco::PFMET > MinPFMETProducer
std::vector< T > METCollection
assert(be >=bs)
std::string moduleLabel_
std::vector< edm::EDGetTokenT< METCollection > > vInputToken
~MinMETProducerT() override
def move
Definition: eostools.py:511
void produce(edm::Event &evt, const edm::EventSetup &es) override
MinMETProducerT< reco::CaloMET > MinCaloMETProducer
MinMETProducerT(const edm::ParameterSet &cfg)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
long double T