CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
40 
43 
46 
49 
50 //
51 // constants, enums and typedefs
52 //
53 
54 //
55 // static data member definitions
56 //
57 
58 //
59 // constructors and destructor
60 //
61 //using namespace std;
62 
63 namespace cms {
65  edm::ProducesCollector producesCollector,
67  : Mtag_(iConfig.getParameter<edm::InputTag>("vtxTag")),
68  fallbackMtag_(iConfig.getParameter<edm::InputTag>("vtxFallbackTag")),
69  saveVtxTimes_(iConfig.getParameter<bool>("saveVtxTimes")) {
70  edm::LogInfo("PixelDigitizer ") << "Enter the Pixel Digitizer";
71 
72  const std::string alias("PileupVertexAccum");
73 
74  producesCollector.produces<PileupVertexContent>().setBranchAlias(alias);
75 
78  }
79 
81 
82  //
83  // member functions
84  //
85 
87  // Make sure that the first crossing processed starts indexing the minbias events from zero.
88 
89  pT_Hats_.clear();
90  z_posns_.clear();
91  t_posns_.clear();
92  }
93 
95  // don't do anything for hard-scatter signal events
96  }
97 
99  edm::EventSetup const& iSetup,
100  edm::StreamID const& streamID) {
102  iEvent.getByLabel(Mtag_, MCevt);
103  if (MCevt.whyFailed()) {
104  iEvent.getByLabel(fallbackMtag_, MCevt);
105  }
106 
107  const HepMC::GenEvent* myGenEvent = MCevt->GetEvent();
108 
109  double pthat = myGenEvent->event_scale();
110  float pt_hat = float(pthat);
111 
112  pT_Hats_.push_back(pt_hat);
113 
114  HepMC::GenEvent::vertex_const_iterator viter;
115  HepMC::GenEvent::vertex_const_iterator vbegin = myGenEvent->vertices_begin();
116  HepMC::GenEvent::vertex_const_iterator vend = myGenEvent->vertices_end();
117 
118  // for production point, pick first vertex
119  viter = vbegin;
120 
121  if (viter != vend) {
122  // The origin vertex (turn it to cm's from GenEvent mm's)
123  HepMC::GenVertex* v = *viter;
124  float zpos = v->position().z() * 0.1;
125 
126  z_posns_.push_back(zpos);
127 
128  if (saveVtxTimes_) {
129  float tpos = v->position().t() / 299792458e-6; // turn from mm to ns
130  t_posns_.push_back(tpos);
131  }
132  }
133 
134  // delete myGenEvent;
135  }
136 
137  // ------------ method called to produce write the data ------------
139  std::unique_ptr<PileupVertexContent> PUVtxC(new PileupVertexContent(pT_Hats_, z_posns_, t_posns_));
140 
141  // write output to event
142  iEvent.put(std::move(PUVtxC));
143  }
144 
145 } // namespace cms
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
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
def move
Definition: eostools.py:511
Log< level::Info, false > LogInfo
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
std::shared_ptr< cms::Exception > whyFailed() const
Definition: HandleBase.h:91