Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:45

0001 #ifndef RecoMET_METProducers_MinMETProducerT_h
0002 #define RecoMET_METProducers_MinMETProducerT_h
0003 
0004 /** \class MinMETProducerT
0005  *
0006  * Produce MET object representing the minimum missing transverse energy
0007  * of set of MET objects given as input
0008  *
0009  * NOTE: class is templated to that it works with reco::CaloMET as well as with reco::PFMET objects as input
0010  *
0011  * \author Christian Veelken, LLR
0012  *
0013  *
0014  *
0015  */
0016 
0017 #include "FWCore/Framework/interface/stream/EDProducer.h"
0018 #include "FWCore/Framework/interface/Event.h"
0019 #include "FWCore/Framework/interface/EventSetup.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/Utilities/interface/InputTag.h"
0022 #include "FWCore/Utilities/interface/Exception.h"
0023 #include "FWCore/Framework/interface/ConsumesCollector.h"
0024 
0025 #include <vector>
0026 
0027 template <typename T>
0028 class MinMETProducerT : public edm::stream::EDProducer<> {
0029   typedef std::vector<T> METCollection;
0030 
0031 public:
0032   explicit MinMETProducerT(const edm::ParameterSet& cfg)
0033       : moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
0034     src_ = cfg.getParameter<vInputTag>("src");
0035 
0036     for (vInputTag::const_iterator src_i = src_.begin(); src_i != src_.end(); ++src_i) {
0037       src_token_.push_back(consumes<METCollection>(*src_i));
0038     }
0039     produces<METCollection>();
0040   }
0041   ~MinMETProducerT() override {}
0042 
0043 private:
0044   void produce(edm::Event& evt, const edm::EventSetup& es) override {
0045     auto outputMETs = std::make_unique<METCollection>();
0046 
0047     // check that all MET collections given as input have the same number of entries
0048     int numMEtObjects = -1;
0049     //    for ( vInputTag::const_iterator src_i = src_.begin();
0050     //        src_i != src_.end(); ++src_i ) {
0051     for (typename vInputToken::const_iterator src_i = src_token_.begin(); src_i != src_token_.end(); ++src_i) {
0052       edm::Handle<METCollection> inputMETs;
0053       //      evt.getByLabel(*src_i, inputMETs);
0054       evt.getByToken(*src_i, inputMETs);
0055       if (numMEtObjects == -1)
0056         numMEtObjects = inputMETs->size();
0057       else if (numMEtObjects != (int)inputMETs->size())
0058         throw cms::Exception("MinMETProducer::produce") << "Mismatch in number of input MET objects !!\n";
0059     }
0060 
0061     for (int iMEtObject = 0; iMEtObject < numMEtObjects; ++iMEtObject) {
0062       const T* minMET = nullptr;
0063       //      for ( vInputTag::const_iterator src_i = src_.begin();
0064       //        src_i != src_.end(); ++src_i ) {
0065       for (typename vInputToken::const_iterator src_i = src_token_.begin(); src_i != src_token_.end(); ++src_i) {
0066         edm::Handle<METCollection> inputMETs;
0067         //  evt.getByLabel(*src_i, inputMETs);
0068         evt.getByToken(*src_i, inputMETs);
0069         const T& inputMET = inputMETs->at(iMEtObject);
0070         if (minMET == nullptr || inputMET.pt() < minMET->pt())
0071           minMET = &inputMET;
0072       }
0073       assert(minMET);
0074       outputMETs->push_back(T(*minMET));
0075     }
0076 
0077     evt.put(std::move(outputMETs));
0078   }
0079 
0080   std::string moduleLabel_;
0081 
0082   typedef std::vector<edm::InputTag> vInputTag;
0083   vInputTag src_;
0084   typedef std::vector<edm::EDGetTokenT<METCollection> > vInputToken;
0085   vInputToken src_token_;
0086 };
0087 
0088 #include "DataFormats/METReco/interface/CaloMET.h"
0089 #include "DataFormats/METReco/interface/PFMET.h"
0090 
0091 namespace reco {
0092   typedef MinMETProducerT<reco::CaloMET> MinCaloMETProducer;
0093   typedef MinMETProducerT<reco::PFMET> MinPFMETProducer;
0094 }  // namespace reco
0095 
0096 #endif