Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    METProducers
0004 // Class:      PFMETProducer
0005 //
0006 //
0007 
0008 //____________________________________________________________________________||
0009 #include "RecoMET/METProducers/interface/PFMETProducer.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/ParameterSet/interface/EmptyGroupDescription.h"
0012 
0013 //____________________________________________________________________________||
0014 namespace cms {
0015 
0016   //____________________________________________________________________________||
0017   PFMETProducer::PFMETProducer(const edm::ParameterSet& iConfig)
0018       : src_(iConfig.getParameter<edm::InputTag>("src")),
0019         inputToken_(consumes<edm::View<reco::Candidate>>(src_)),
0020         calculateSignificance_(iConfig.getParameter<bool>("calculateSignificance")),
0021         globalThreshold_(iConfig.getParameter<double>("globalThreshold")),
0022         applyWeight_(iConfig.getParameter<bool>("applyWeight")),
0023         weights_(nullptr) {
0024     if (applyWeight_) {
0025       edm::InputTag srcWeights = iConfig.getParameter<edm::InputTag>("srcWeights");
0026       if (srcWeights.label().empty())
0027         throw cms::Exception("InvalidInput") << "applyWeight set to True, but no weights given in PFMETProducer\n";
0028       if (srcWeights.label() == src_.label())
0029         edm::LogWarning("PFMETProducer")
0030             << "Particle and weights collection have the same label. You may be applying the same weights twice.\n";
0031       weightsToken_ = consumes<edm::ValueMap<float>>(srcWeights);
0032     }
0033     if (calculateSignificance_) {
0034       metSigAlgo_ = new metsig::METSignificance(iConfig);
0035 
0036       jetToken_ = mayConsume<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("srcJets"));
0037       std::vector<edm::InputTag> srcLeptonsTags = iConfig.getParameter<std::vector<edm::InputTag>>("srcLeptons");
0038       for (std::vector<edm::InputTag>::const_iterator it = srcLeptonsTags.begin(); it != srcLeptonsTags.end(); it++) {
0039         lepTokens_.push_back(mayConsume<edm::View<reco::Candidate>>(*it));
0040       }
0041 
0042       jetSFToken_ = esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("srcJetSF")));
0043       jetResPtToken_ = esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("srcJetResPt")));
0044       jetResPhiToken_ = esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("srcJetResPhi")));
0045       rhoToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("srcRho"));
0046     }
0047 
0048     std::string alias = iConfig.getParameter<std::string>("alias");
0049 
0050     produces<reco::PFMETCollection>().setBranchAlias(alias);
0051   }
0052 
0053   //____________________________________________________________________________||
0054   void PFMETProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0055     edm::Handle<edm::View<reco::Candidate>> input;
0056     event.getByToken(inputToken_, input);
0057 
0058     if (applyWeight_)
0059       weights_ = &event.get(weightsToken_);
0060 
0061     METAlgo algo;
0062     CommonMETData commonMETdata;
0063     commonMETdata = algo.run(*input.product(), globalThreshold_, weights_);
0064 
0065     const math::XYZTLorentzVector p4(commonMETdata.mex, commonMETdata.mey, 0.0, commonMETdata.met);
0066     const math::XYZPoint vtx(0.0, 0.0, 0.0);
0067 
0068     PFSpecificAlgo pf;
0069     SpecificPFMETData specific;
0070     specific = pf.run(*input.product(), weights_);
0071 
0072     reco::PFMET pfmet(specific, commonMETdata.sumet, p4, vtx);
0073     pfmet.setIsWeighted(applyWeight_);
0074 
0075     if (calculateSignificance_) {
0076       reco::METCovMatrix sigcov = getMETCovMatrix(event, setup, input);
0077       pfmet.setSignificanceMatrix(sigcov);
0078     }
0079 
0080     auto pfmetcoll = std::make_unique<reco::PFMETCollection>();
0081 
0082     pfmetcoll->push_back(pfmet);
0083     event.put(std::move(pfmetcoll));
0084   }
0085 
0086   reco::METCovMatrix PFMETProducer::getMETCovMatrix(const edm::Event& event,
0087                                                     const edm::EventSetup& setup,
0088                                                     const edm::Handle<edm::View<reco::Candidate>>& candInput) const {
0089     // leptons
0090     std::vector<edm::Handle<reco::CandidateView>> leptons;
0091     for (std::vector<edm::EDGetTokenT<edm::View<reco::Candidate>>>::const_iterator srcLeptons_i = lepTokens_.begin();
0092          srcLeptons_i != lepTokens_.end();
0093          ++srcLeptons_i) {
0094       edm::Handle<reco::CandidateView> leptons_i;
0095       event.getByToken(*srcLeptons_i, leptons_i);
0096       leptons.push_back(leptons_i);
0097       /*
0098       for ( reco::CandidateView::const_iterator lepton = leptons_i->begin();
0099         lepton != leptons_i->end(); ++lepton ) {
0100         leptons.push_back(*lepton);
0101       }
0102      */
0103     }
0104 
0105     // jets
0106     edm::Handle<edm::View<reco::Jet>> inputJets;
0107     event.getByToken(jetToken_, inputJets);
0108 
0109     JME::JetResolution resPtObj(setup.getData(jetResPtToken_));
0110     JME::JetResolution resPhiObj(setup.getData(jetResPhiToken_));
0111     JME::JetResolutionScaleFactor resSFObj(setup.getData(jetSFToken_));
0112 
0113     edm::Handle<double> rho;
0114     event.getByToken(rhoToken_, rho);
0115 
0116     //Compute the covariance matrix and fill it
0117     double sumPtUnclustered = 0;
0118     reco::METCovMatrix cov = metSigAlgo_->getCovariance(*inputJets,
0119                                                         leptons,
0120                                                         candInput,
0121                                                         *rho,
0122                                                         resPtObj,
0123                                                         resPhiObj,
0124                                                         resSFObj,
0125                                                         event.isRealData(),
0126                                                         sumPtUnclustered,
0127                                                         weights_);
0128 
0129     return cov;
0130   }
0131 
0132   void PFMETProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0133     edm::ParameterSetDescription desc;
0134     desc.add<edm::InputTag>("src", edm::InputTag("particleFlow"));
0135     desc.add<double>("globalThreshold", 0.);
0136     desc.add<std::string>("alias", "@module_label");
0137     desc.add<bool>("calculateSignificance", false);
0138     desc.addOptional<edm::InputTag>("srcJets");
0139     desc.addOptional<std::vector<edm::InputTag>>("srcLeptons");
0140     desc.addOptional<std::string>("srcJetSF");
0141     desc.addOptional<std::string>("srcJetResPt");
0142     desc.addOptional<std::string>("srcJetResPhi");
0143     desc.addOptional<edm::InputTag>("srcRho");
0144     edm::ParameterSetDescription params;
0145     params.setAllowAnything();  // FIXME: This still needs to be defined in METSignficance
0146     desc.addOptional<edm::ParameterSetDescription>("parameters", params);
0147     edm::ParameterSetDescription desc1 = desc;
0148     edm::ParameterSetDescription desc2 = desc;
0149     desc1.add<bool>("applyWeight", false);
0150     desc1.add<edm::InputTag>("srcWeights", edm::InputTag(""));
0151     descriptions.add("pfMet", desc1);
0152     desc2.add<bool>("applyWeight", true);
0153     desc2.add<edm::InputTag>("srcWeights", edm::InputTag("puppiNoLep"));
0154     descriptions.add("pfMetPuppi", desc2);
0155   }
0156 
0157   //____________________________________________________________________________||
0158   DEFINE_FWK_MODULE(PFMETProducer);
0159 }  // namespace cms
0160 
0161 //____________________________________________________________________________||