Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:35

0001 // -*- C++ -*-
0002 //
0003 // Package:    PileUp
0004 // Class:      MCVerticesWeight
0005 //
0006 /**\class MCVerticesWeight MCVerticesWeight.cc Validation/RecoVertex/MCVerticesWeight.cc
0007 
0008  Description: 
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Venturi
0015 //         Created:  Tue Oct 21 20:55:22 CEST 2008
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 #include <string>
0022 
0023 // user include files
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/global/EDFilter.h"
0026 
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 #include "FWCore/Framework/interface/ESHandle.h"
0030 #include "FWCore/Framework/interface/ESWatcher.h"
0031 
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 
0036 #include "FWCore/Utilities/interface/InputTag.h"
0037 
0038 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
0039 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0040 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
0041 
0042 #include "Validation/RecoVertex/interface/VertexWeighter.h"
0043 
0044 //
0045 // class declaration
0046 //
0047 
0048 class MCVerticesWeight : public edm::global::EDFilter<> {
0049 public:
0050   explicit MCVerticesWeight(const edm::ParameterSet&);
0051   ~MCVerticesWeight() override;
0052 
0053 private:
0054   bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0055 
0056   // ----------member data ---------------------------
0057 
0058   edm::EDGetTokenT<std::vector<PileupSummaryInfo> > m_vecPileupSummaryInfoToken;
0059   edm::EDGetTokenT<edm::HepMCProduct> m_hepMCProductToken;
0060   const VertexWeighter m_weighter;
0061 };
0062 
0063 //
0064 // constructors and destructor
0065 //
0066 MCVerticesWeight::MCVerticesWeight(const edm::ParameterSet& iConfig)
0067     : m_vecPileupSummaryInfoToken(
0068           consumes<std::vector<PileupSummaryInfo> >(iConfig.getParameter<edm::InputTag>("pileupSummaryCollection"))),
0069       m_hepMCProductToken(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("mcTruthCollection"))),
0070       m_weighter(iConfig.getParameter<edm::ParameterSet>("weighterConfig")) {
0071   produces<double>();
0072 }
0073 
0074 MCVerticesWeight::~MCVerticesWeight() {}
0075 
0076 //
0077 // member functions
0078 //
0079 
0080 // ------------ method called on each new Event  ------------
0081 bool MCVerticesWeight::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0082   bool selected = true;
0083 
0084   double computed_weight(1);
0085 
0086   edm::Handle<std::vector<PileupSummaryInfo> > pileupinfos;
0087   iEvent.getByToken(m_vecPileupSummaryInfoToken, pileupinfos);
0088 
0089   // look for the intime PileupSummaryInfo
0090 
0091   std::vector<PileupSummaryInfo>::const_iterator pileupinfo;
0092   for (pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end(); ++pileupinfo) {
0093     if (pileupinfo->getBunchCrossing() == 0)
0094       break;
0095   }
0096 
0097   //
0098   if (pileupinfo->getBunchCrossing() != 0) {
0099     edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing();
0100   } else {
0101     //     pileupinfo->getPU_NumInteractions();
0102 
0103     const std::vector<float>& zpositions = pileupinfo->getPU_zpositions();
0104 
0105     //     for(std::vector<float>::const_iterator zpos = zpositions.begin() ; zpos != zpositions.end() ; ++zpos) {
0106 
0107     //     }
0108 
0109     // main interaction part
0110 
0111     edm::Handle<edm::HepMCProduct> EvtHandle;
0112     iEvent.getByToken(m_hepMCProductToken, EvtHandle);
0113 
0114     const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
0115 
0116     // get the first vertex
0117 
0118     double zmain = 0.0;
0119     if (Evt->vertices_begin() != Evt->vertices_end()) {
0120       zmain = (*Evt->vertices_begin())->point3d().z() / 10.;
0121     }
0122 
0123     //
0124 
0125     computed_weight = m_weighter.weight(zpositions, zmain);
0126   }
0127 
0128   std::unique_ptr<double> weight(new double(computed_weight));
0129 
0130   iEvent.put(std::move(weight));
0131 
0132   //
0133 
0134   return selected;
0135 }
0136 
0137 //define this as a plug-in
0138 DEFINE_FWK_MODULE(MCVerticesWeight);