CMS 3D CMS Logo

PileupVertexAccumulator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PileupVertexAccumulator
4 // Class: PileupVertexAccumulator
5 //
13 //
14 // Original Author: Mike Hildreth - Notre Dame
15 // Created: Wed Jan 21 05:14:48 CET 2015
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <set>
22 
23 // user include files
25 
31 
36 
37 // user include files
39 
42 
45 
48 
49 //
50 // constants, enums and typedefs
51 //
52 
53 //
54 // static data member definitions
55 //
56 
57 //
58 // constructors and destructor
59 //
60 //using namespace std;
61 
62 namespace cms {
64  edm::ProducesCollector producesCollector,
66  : Mtag_(iConfig.getParameter<edm::InputTag>("vtxTag")),
67  fallbackMtag_(iConfig.getParameter<edm::InputTag>("vtxFallbackTag")),
68  saveVtxTimes_(iConfig.getParameter<bool>("saveVtxTimes")) {
69  edm::LogInfo("PixelDigitizer ") << "Enter the Pixel Digitizer";
70 
71  const std::string alias("PileupVertexAccum");
72 
73  producesCollector.produces<PileupVertexContent>().setBranchAlias(alias);
74 
77  }
78 
80 
81  //
82  // member functions
83  //
84 
86  // Make sure that the first crossing processed starts indexing the minbias events from zero.
87 
88  pT_Hats_.clear();
89  z_posns_.clear();
90  t_posns_.clear();
91  }
92 
94  // don't do anything for hard-scatter signal events
95  }
96 
98  edm::EventSetup const& iSetup,
99  edm::StreamID const& streamID) {
101  iEvent.getByLabel(Mtag_, MCevt);
102  if (MCevt.whyFailed()) {
103  iEvent.getByLabel(fallbackMtag_, MCevt);
104  }
105 
106  const HepMC::GenEvent* myGenEvent = MCevt->GetEvent();
107 
108  double pthat = myGenEvent->event_scale();
109  float pt_hat = float(pthat);
110 
111  pT_Hats_.push_back(pt_hat);
112 
113  HepMC::GenEvent::vertex_const_iterator viter;
114  HepMC::GenEvent::vertex_const_iterator vbegin = myGenEvent->vertices_begin();
115  HepMC::GenEvent::vertex_const_iterator vend = myGenEvent->vertices_end();
116 
117  // for production point, pick first vertex
118  viter = vbegin;
119 
120  if (viter != vend) {
121  // The origin vertex (turn it to cm's from GenEvent mm's)
122  HepMC::GenVertex* v = *viter;
123  float zpos = v->position().z() * 0.1;
124 
125  z_posns_.push_back(zpos);
126 
127  if (saveVtxTimes_) {
128  float tpos = v->position().t() / 299792458e-6; // turn from mm to ns
129  t_posns_.push_back(tpos);
130  }
131  }
132 
133  // delete myGenEvent;
134  }
135 
136  // ------------ method called to produce write the data ------------
138  std::unique_ptr<PileupVertexContent> PUVtxC(new PileupVertexContent(pT_Hats_, z_posns_, t_posns_));
139 
140  // write output to event
141  iEvent.put(std::move(PUVtxC));
142  }
143 
144 } // namespace cms
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
ProductRegistryHelper::BranchAliasSetterT< ProductType > produces()
PileupVertexAccumulator(const edm::ParameterSet &conf, edm::ProducesCollector, edm::ConsumesCollector &iC)
void finalizeEvent(edm::Event &e, edm::EventSetup const &c) override
void accumulate(edm::Event const &e, edm::EventSetup const &c) override
int iEvent
Definition: GenABIO.cc:224
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:91
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
Namespace of DDCMS conversion namespace.
Log< level::Info, false > LogInfo
HLT enums.
void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
def move(src, dest)
Definition: eostools.py:511