Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "JetMETCorrections/Type1MET/plugins/Type0PFMETcorrInputProducer.h"
0002 
0003 #include "DataFormats/Common/interface/AssociationMap.h"
0004 #include "DataFormats/Common/interface/OneToManyWithQuality.h"
0005 #include "DataFormats/VertexReco/interface/Vertex.h"
0006 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0007 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "DataFormats/METReco/interface/CorrMETData.h"
0010 
0011 #include <TMath.h>
0012 
0013 Type0PFMETcorrInputProducer::Type0PFMETcorrInputProducer(const edm::ParameterSet& cfg)
0014     : moduleLabel_(cfg.getParameter<std::string>("@module_label")), correction_(nullptr) {
0015   pfCandidateToVertexAssociationsToken_ =
0016       consumes<PFCandToVertexAssMap>(cfg.getParameter<edm::InputTag>("srcPFCandidateToVertexAssociations"));
0017   hardScatterVertexToken_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("srcHardScatterVertex"));
0018 
0019   edm::ParameterSet cfgCorrection_function = cfg.getParameter<edm::ParameterSet>("correction");
0020   std::string corrFunctionName = std::string(moduleLabel_).append("correction");
0021   std::string corrFunctionFormula = cfgCorrection_function.getParameter<std::string>("formula");
0022   correction_ = new TFormula(corrFunctionName.data(), corrFunctionFormula.data());
0023   int numParameter = correction_->GetNpar();
0024   for (int iParameter = 0; iParameter < numParameter; ++iParameter) {
0025     std::string parName = Form("par%i", iParameter);
0026     double parValue = cfgCorrection_function.getParameter<double>(parName);
0027     correction_->SetParameter(iParameter, parValue);
0028   }
0029 
0030   minDz_ = cfg.getParameter<double>("minDz");
0031 
0032   produces<CorrMETData>();
0033 }
0034 
0035 Type0PFMETcorrInputProducer::~Type0PFMETcorrInputProducer() { delete correction_; }
0036 
0037 void Type0PFMETcorrInputProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
0038   edm::Handle<reco::VertexCollection> hardScatterVertex;
0039   evt.getByToken(hardScatterVertexToken_, hardScatterVertex);
0040 
0041   edm::Handle<PFCandToVertexAssMap> pfCandidateToVertexAssociations;
0042   evt.getByToken(pfCandidateToVertexAssociationsToken_, pfCandidateToVertexAssociations);
0043 
0044   std::unique_ptr<CorrMETData> pfMEtCorrection(new CorrMETData());
0045 
0046   for (PFCandToVertexAssMap::const_iterator pfCandidateToVertexAssociation = pfCandidateToVertexAssociations->begin();
0047        pfCandidateToVertexAssociation != pfCandidateToVertexAssociations->end();
0048        ++pfCandidateToVertexAssociation) {
0049     reco::VertexRef vertex = pfCandidateToVertexAssociation->key;
0050     const PFCandQualityPairVector& pfCandidates_vertex = pfCandidateToVertexAssociation->val;
0051 
0052     bool isHardScatterVertex = false;
0053     for (reco::VertexCollection::const_iterator hardScatterVertex_i = hardScatterVertex->begin();
0054          hardScatterVertex_i != hardScatterVertex->end();
0055          ++hardScatterVertex_i) {
0056       if (TMath::Abs(vertex->position().z() - hardScatterVertex_i->position().z()) < minDz_) {
0057         isHardScatterVertex = true;
0058         break;
0059       }
0060     }
0061 
0062     if (!isHardScatterVertex) {
0063       reco::Candidate::LorentzVector sumChargedPFCandP4_vertex;
0064       for (PFCandQualityPairVector::const_iterator pfCandidate_vertex = pfCandidates_vertex.begin();
0065            pfCandidate_vertex != pfCandidates_vertex.end();
0066            ++pfCandidate_vertex) {
0067         const reco::PFCandidate& pfCandidate = (*pfCandidate_vertex->first);
0068         if (pfCandidate.particleId() == reco::PFCandidate::h || pfCandidate.particleId() == reco::PFCandidate::e ||
0069             pfCandidate.particleId() == reco::PFCandidate::mu) {
0070           sumChargedPFCandP4_vertex += pfCandidate.p4();
0071         }
0072       }
0073 
0074       double pt = sumChargedPFCandP4_vertex.pt();
0075       double phi = sumChargedPFCandP4_vertex.phi();
0076       double ptCorr = correction_->Eval(pt);
0077       double pxCorr = TMath::Cos(phi) * ptCorr;
0078       double pyCorr = TMath::Sin(phi) * ptCorr;
0079 
0080       pfMEtCorrection->mex += pxCorr;
0081       pfMEtCorrection->mey += pyCorr;
0082     }
0083   }
0084 
0085   evt.put(std::move(pfMEtCorrection));
0086 }
0087 
0088 #include "FWCore/Framework/interface/MakerMacros.h"
0089 
0090 DEFINE_FWK_MODULE(Type0PFMETcorrInputProducer);