CMS 3D CMS Logo

MCVerticesWeight.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PileUp
4 // Class: MCVerticesWeight
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Tue Oct 21 20:55:22 CEST 2008
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <string>
22 
23 // user include files
26 
31 
33 
35 
37 
41 
43 
44 //
45 // class declaration
46 //
47 
49 public:
50  explicit MCVerticesWeight(const edm::ParameterSet&);
51  ~MCVerticesWeight() override;
52 
53 private:
54  bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
55 
56  // ----------member data ---------------------------
57 
61 };
62 
63 //
64 // constructors and destructor
65 //
67  : m_vecPileupSummaryInfoToken(
68  consumes<std::vector<PileupSummaryInfo> >(iConfig.getParameter<edm::InputTag>("pileupSummaryCollection"))),
69  m_hepMCProductToken(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("mcTruthCollection"))),
70  m_weighter(iConfig.getParameter<edm::ParameterSet>("weighterConfig")) {
71  produces<double>();
72 }
73 
75 
76 //
77 // member functions
78 //
79 
80 // ------------ method called on each new Event ------------
82  bool selected = true;
83 
84  double computed_weight(1);
85 
87  iEvent.getByToken(m_vecPileupSummaryInfoToken, pileupinfos);
88 
89  // look for the intime PileupSummaryInfo
90 
91  std::vector<PileupSummaryInfo>::const_iterator pileupinfo;
92  for (pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end(); ++pileupinfo) {
93  if (pileupinfo->getBunchCrossing() == 0)
94  break;
95  }
96 
97  //
98  if (pileupinfo->getBunchCrossing() != 0) {
99  edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing();
100  } else {
101  // pileupinfo->getPU_NumInteractions();
102 
103  const std::vector<float>& zpositions = pileupinfo->getPU_zpositions();
104 
105  // for(std::vector<float>::const_iterator zpos = zpositions.begin() ; zpos != zpositions.end() ; ++zpos) {
106 
107  // }
108 
109  // main interaction part
110 
112  iEvent.getByToken(m_hepMCProductToken, EvtHandle);
113 
114  const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
115 
116  // get the first vertex
117 
118  double zmain = 0.0;
119  if (Evt->vertices_begin() != Evt->vertices_end()) {
120  zmain = (*Evt->vertices_begin())->point3d().z() / 10.;
121  }
122 
123  //
124 
125  computed_weight = m_weighter.weight(zpositions, zmain);
126  }
127 
128  std::unique_ptr<double> weight(new double(computed_weight));
129 
130  iEvent.put(std::move(weight));
131 
132  //
133 
134  return selected;
135 }
136 
137 //define this as a plug-in
const double weight(const std::vector< float > &zpositions, const float &zmain) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Definition: weight.py:1
Log< level::Error, false > LogError
MCVerticesWeight(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > m_vecPileupSummaryInfoToken
bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int iEvent
Definition: GenABIO.cc:224
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
const VertexWeighter m_weighter
edm::EDGetTokenT< edm::HepMCProduct > m_hepMCProductToken
~MCVerticesWeight() override
HLT enums.
def move(src, dest)
Definition: eostools.py:511