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 
25 
26 #include <vector>
27 
28 template <typename T>
30 {
31  typedef std::vector<T> METCollection;
32 
33  public:
34 
35  explicit MinMETProducerT(const edm::ParameterSet& cfg)
36  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
37  {
38  src_ = cfg.getParameter<vInputTag>("src");
39 
40  produces<METCollection>();
41  }
43 
44  private:
45 
46  void produce(edm::Event& evt, const edm::EventSetup& es)
47  {
48  std::auto_ptr<METCollection> outputMETs(new METCollection());
49 
50  // check that all MET collections given as input have the same number of entries
51  int numMEtObjects = -1;
52  for ( vInputTag::const_iterator src_i = src_.begin();
53  src_i != src_.end(); ++src_i ) {
55  evt.getByLabel(*src_i, inputMETs);
56  if ( numMEtObjects == -1 ) numMEtObjects = inputMETs->size();
57  else if ( numMEtObjects != (int)inputMETs->size() )
58  throw cms::Exception("MinMETProducer::produce")
59  << "Mismatch in number of input MET objects !!\n";
60  }
61 
62  for ( int iMEtObject = 0; iMEtObject < numMEtObjects; ++iMEtObject ) {
63  const T* minMET = 0;
64  for ( vInputTag::const_iterator src_i = src_.begin();
65  src_i != src_.end(); ++src_i ) {
67  evt.getByLabel(*src_i, inputMETs);
68  const T& inputMET = inputMETs->at(iMEtObject);
69  if ( minMET == 0 || inputMET.pt() < minMET->pt() ) minMET = &inputMET;
70  }
71  assert(minMET);
72  outputMETs->push_back(T(*minMET));
73  }
74 
75  evt.put(outputMETs);
76  }
77 
78  std::string moduleLabel_;
79 
80  typedef std::vector<edm::InputTag> vInputTag;
82 };
83 
86 
87 namespace reco
88 {
91 }
92 
93 #endif
94 
95 
96 
97 
T getParameter(std::string const &) const
std::vector< edm::InputTag > vInputTag
MinMETProducerT< reco::PFMET > MinPFMETProducer
std::vector< T > METCollection
std::string moduleLabel_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
MinMETProducerT< reco::CaloMET > MinCaloMETProducer
void produce(edm::Event &evt, const edm::EventSetup &es)
MinMETProducerT(const edm::ParameterSet &cfg)
long double T