CMS 3D CMS Logo

PileupInformation.cc
Go to the documentation of this file.
1 // File: PileupInformation.cc
2 // Description: adds pileup information object to event
3 // Author: Mike Hildreth
4 //
5 // Adds a vector of PileupSummaryInfo objects to the event.
6 // One for each bunch crossing.
7 //
8 //--------------------------------------------
16 
18 
19 
21 {
22  // Initialize global parameters
23 
24  pTcut_1_ = 0.1;
25  pTcut_2_ = 0.5; // defaults
26  isPreMixed_ = config.getParameter<bool>("isPreMixed");
27 
28  if ( !isPreMixed_ ) {
29  distanceCut_ = config.getParameter<double>("vertexDistanceCut");
30  volumeRadius_ = config.getParameter<double>("volumeRadius");
31  volumeZ_ = config.getParameter<double>("volumeZ");
32  pTcut_1_ = config.getParameter<double>("pTcut_1");
33  pTcut_2_ = config.getParameter<double>("pTcut_2");
34 
35  PileupInfoLabel_ = consumes<PileupMixingContent>(config.getParameter<edm::InputTag>("PileupMixingLabel"));
36 
37  PileupVtxLabel_ = consumes<PileupVertexContent>(config.getParameter<edm::InputTag>("PileupMixingLabel"));
38 
39  LookAtTrackingTruth_ = config.getUntrackedParameter<bool>("doTrackTruth");
40 
41  trackingTruthT_ = mayConsume<TrackingParticleCollection>(config.getParameter<edm::InputTag>("TrackingParticlesLabel"));
42  trackingTruthV_ = mayConsume<TrackingVertexCollection>(config.getParameter<edm::InputTag>("TrackingParticlesLabel"));
43 
44  saveVtxTimes_ = config.getParameter<bool>("saveVtxTimes");
45 
46  MessageCategory_ = "PileupInformation";
47 
48  edm::LogInfo (MessageCategory_) << "Setting up PileupInformation";
49  edm::LogInfo (MessageCategory_) << "Vertex distance cut set to " << distanceCut_ << " mm";
50  edm::LogInfo (MessageCategory_) << "Volume radius set to " << volumeRadius_ << " mm";
51  edm::LogInfo (MessageCategory_) << "Volume Z set to " << volumeZ_ << " mm";
52  edm::LogInfo (MessageCategory_) << "Lower pT Threshold set to " << pTcut_1_ << " GeV";
53  edm::LogInfo (MessageCategory_) << "Upper pT Threshold set to " << pTcut_2_ << " GeV";
54  }
55  else{
56  pileupSummaryToken_=consumes<std::vector<PileupSummaryInfo> >(config.getParameter<edm::InputTag>("PileupSummaryInfoInputTag"));
57  bunchSpacingToken_=consumes<int>(config.getParameter<edm::InputTag>("BunchSpacingInputTag"));
58  }
59 
60  produces< std::vector<PileupSummaryInfo> >();
61  produces<int>("bunchSpacing");
62 
63  //produces<PileupSummaryInfo>();
64 }
65 
66 
68 {
69 
70  std::unique_ptr<std::vector<PileupSummaryInfo> > PSIVector(new std::vector<PileupSummaryInfo>);
71 
72  if ( isPreMixed_ ) {
74  event.getByToken(pileupSummaryToken_,psiInput);
75 
76  std::vector<PileupSummaryInfo>::const_iterator PSiter;
77 
78  for(PSiter = psiInput.product()->begin(); PSiter != psiInput.product()->end(); PSiter++){
79 
80  PSIVector->push_back(*PSiter);
81  }
82 
83 
84  edm::Handle< int> bsInput;
85  event.getByToken(bunchSpacingToken_,bsInput);
86  int bunchSpacing=*(bsInput.product());
87 
88  event.put(std::move(PSIVector));
89 
90  //add bunch spacing to the event as a seperate integer for use by downstream modules
91  std::unique_ptr<int> bunchSpacingP(new int(bunchSpacing));
92  event.put(std::move(bunchSpacingP),"bunchSpacing");
93 
94  return;
95  }
96 
97  edm::Handle< PileupMixingContent > MixingPileup; // Get True pileup information from MixingModule
98  event.getByToken(PileupInfoLabel_, MixingPileup);
99 
100  std::vector<int> BunchCrossings;
101  std::vector<int> Interactions_Xing;
102  std::vector<float> TrueInteractions_Xing;
103  std::vector< std::vector<edm::EventID> > eventInfoList_Xing;
104 
105  int bunchSpacing;
106 
107  const PileupMixingContent* MixInfo = MixingPileup.product();
108 
109  if(MixInfo) { // extract information - way easier than counting vertices
110 
111  const std::vector<int>& bunchCrossing = MixInfo->getMix_bunchCrossing();
112  const std::vector<int>& interactions = MixInfo->getMix_Ninteractions();
113  const std::vector<float>& TrueInteractions = MixInfo->getMix_TrueInteractions();
114  const std::vector<edm::EventID> eventInfoList= MixInfo->getMix_eventInfo();
115 
116  bunchSpacing = MixInfo->getMix_bunchSpacing();
117  unsigned int totalIntPU=0;
118 
119  for(int ib=0; ib<(int)bunchCrossing.size(); ++ib){
120  // std::cout << " bcr, nint " << bunchCrossing[ib] << " " << interactions[ib] << std::endl;
121  BunchCrossings.push_back(bunchCrossing[ib]);
122  Interactions_Xing.push_back(interactions[ib]);
123  TrueInteractions_Xing.push_back(TrueInteractions[ib]);
124 
125  std::vector<edm::EventID> eventInfos;
126  eventInfos.reserve( interactions[ib] );
127  for ( int pu=0; pu< interactions[ib]; pu++) {
128  eventInfos.push_back(eventInfoList[totalIntPU+pu]);
129  }
130  totalIntPU+=(interactions[ib]);
131  eventInfoList_Xing.push_back(eventInfos);
132 
133  }
134  }
135  else{ // have to throw an exception..
136 
137  throw cms::Exception("PileupInformation") << " PileupMixingContent is missing from the event.\n"
138  "There must be some breakdown in the Simulation Chain.\n"
139  "You must run the MixingModule before calling this routine.";
140 
141  }
142 
143  // store information from pileup vertices, if it's in the event. Have to loop on interactions again.
144 
145  edm::Handle< PileupVertexContent > MixingPileupVtx; // Get True pileup information from MixingModule
146  event.getByToken(PileupVtxLabel_, MixingPileupVtx);
147 
148  const PileupVertexContent* MixVtxInfo = MixingPileupVtx.product();
149 
150  std::vector< std::vector<float> > ptHatList_Xing;
151  std::vector< std::vector<float> > zPosList_Xing;
152  std::vector< std::vector<float> > tPosList_Xing;
153 
154  bool Have_pThats = false;
155 
156  if(MixVtxInfo) { // extract information - way easier than counting vertices
157 
158 
159  Have_pThats = true;
160 
161  const std::vector<int>& bunchCrossing = MixInfo->getMix_bunchCrossing();
162  const std::vector<int>& interactions = MixInfo->getMix_Ninteractions();
163 
164  const std::vector<float>& PtHatInput = MixVtxInfo->getMix_pT_hats();
165  const std::vector<float>& ZposInput = MixVtxInfo->getMix_z_Vtxs();
166  const std::vector<float>& TposInput = MixVtxInfo->getMix_t_Vtxs();
167 
168  // store information from pileup vertices, if it's in the event:
169 
170  unsigned int totalIntPU=0;
171 
172  for(int ib=0; ib<(int)bunchCrossing.size(); ++ib){
173  // std::cout << " bcr, nint " << bunchCrossing[ib] << " " << interactions[ib] << std::endl;
174 
175  std::vector<float> zposBX, tposBX;
176  std::vector<float> pthatBX;
177  zposBX.reserve( interactions[ib] );
178  if (saveVtxTimes_) tposBX.reserve( interactions[ib] );
179  pthatBX.reserve( interactions[ib] );
180  for ( int pu=0; pu< interactions[ib]; pu++) {
181  zposBX.push_back(ZposInput[totalIntPU+pu]);
182  if (saveVtxTimes_) tposBX.push_back(TposInput[totalIntPU+pu]);
183  pthatBX.push_back(PtHatInput[totalIntPU+pu]);
184  }
185  totalIntPU+=(interactions[ib]);
186  zPosList_Xing.push_back(zposBX);
187  tPosList_Xing.push_back(tposBX);
188  ptHatList_Xing.push_back(pthatBX);
189 
190  }
191  } // end of VertexInfo block
192 
193 
194  //Now, get information on valid particles that look like they could be in the tracking volume
195 
196 
197  zpositions.clear();
198  sumpT_lowpT.clear();
199  sumpT_highpT.clear();
200  ntrks_lowpT.clear();
201  ntrks_highpT.clear();
202  event_index_.clear();
203 
204  int lastEvent = 0; // zero is the true MC hard-scatter event
205 
206  // int lastBunchCrossing = 0; // 0 is the true bunch crossing, should always come first.
207 
208  bool HaveTrackingParticles = false;
209 
212 
213  TrackingVertexCollection::const_iterator iVtx;
214  TrackingVertexCollection::const_iterator iVtxTest;
215  TrackingParticleCollection::const_iterator iTrackTest;
216 
217  if( LookAtTrackingTruth_ ){
218 
219  if(event.getByToken(trackingTruthT_, mergedPH) && event.getByToken(trackingTruthV_, mergedVH)) {
220 
221  HaveTrackingParticles = true;
222 
223  iVtxTest = mergedVH->begin();
224  iTrackTest = mergedPH->begin();
225  }
226 
227  }
228 
229  int nminb_vtx = 0;
230  // bool First = true;
231  // bool flag_new = false;
232 
233  std::vector<int>::iterator BXIter;
234  std::vector<int>::iterator InteractionsIter = Interactions_Xing.begin();
235  std::vector<float>::iterator TInteractionsIter = TrueInteractions_Xing.begin();
236  std::vector< std::vector<edm::EventID> >::iterator TEventInfoIter = eventInfoList_Xing.begin();
237 
238  std::vector< std::vector<float> >::iterator zPosIter;
239  std::vector< std::vector<float> >::iterator tPosIter;
240  std::vector< std::vector<float> >::iterator pThatIter;
241 
242  if(Have_pThats) {
243  zPosIter = zPosList_Xing.begin();
244  tPosIter = tPosList_Xing.begin();
245  pThatIter = ptHatList_Xing.begin();
246  }
247 
248  // loop over the bunch crossings and interactions we have extracted
249 
250  for( BXIter = BunchCrossings.begin(); BXIter != BunchCrossings.end(); ++BXIter, ++InteractionsIter, ++TInteractionsIter, ++TEventInfoIter) {
251 
252  //std::cout << "looking for BX: " << (*BXIter) << std::endl;
253 
254  if(HaveTrackingParticles) { // leave open the option of not producing TrackingParticles and just keeping # interactions
255 
256  for (iVtx = iVtxTest; iVtx != mergedVH->end(); ++iVtx) {
257 
258  if(iVtx->eventId().bunchCrossing() == (*BXIter) ) { // found first vertex in this bunch crossing
259 
260  if(iVtx->eventId().event() != lastEvent) {
261 
262  //std::cout << "BX,event " << iVtx->eventId().bunchCrossing() << " " << iVtx->eventId().event() << std::endl;
263 
264  float zpos = 0.;
265  zpos = iVtx->position().z();
266  zpositions.push_back(zpos); //save z position of each vertex
267  sumpT_lowpT.push_back(0.);
268  sumpT_highpT.push_back(0.);
269  ntrks_lowpT.push_back(0);
270  ntrks_highpT.push_back(0);
271 
272  lastEvent = iVtx->eventId().event();
273  iVtxTest = --iVtx; // just for security
274 
275  // turns out events aren't sequential... save map of indices
276 
277  event_index_.insert(myindex::value_type(lastEvent,nminb_vtx));
278 
279  ++nminb_vtx;
280 
281  continue;
282  }
283  }
284  }
285 
286  // next loop over tracks to get information
287 
288  for (TrackingParticleCollection::const_iterator iTrack = iTrackTest; iTrack != mergedPH->end(); ++iTrack)
289  {
290  bool FoundTrk = false;
291 
292  float zpos=0.;
293 
294  if(iTrack->eventId().bunchCrossing() == (*BXIter) && iTrack->eventId().event() > 0 )
295  {
296  FoundTrk = true;
297  int correct_index = event_index_[iTrack->eventId().event()];
298 
299  //std::cout << " track index, correct index " << iTrack->eventId().event() << " " << correct_index << std::endl;
300 
301  zpos = zpositions[correct_index];
302  if(iTrack->matchedHit()>0) {
303  if(fabs(iTrack->parentVertex()->position().z()-zpos)<0.1) { //make sure track really comes from this vertex
304  //std::cout << *iTrack << std::endl;
305  float Tpx = iTrack->p4().px();
306  float Tpy = iTrack->p4().py();
307  float TpT = sqrt(Tpx*Tpx + Tpy*Tpy);
308  if( TpT>pTcut_1_ ) {
309  sumpT_lowpT[correct_index]+=TpT;
310  ++ntrks_lowpT[correct_index];
311  }
312  if( TpT>pTcut_2_ ){
313  sumpT_highpT[correct_index]+=TpT;
314  ++ntrks_highpT[correct_index];
315  }
316  }
317  }
318  }
319  else{
320  if(FoundTrk) {
321 
322  iTrackTest = --iTrack; // reset so we can start over next time
323  --iTrackTest; // just to be sure
324  break;
325  }
326 
327  }
328 
329  } // end of track loop
330 
331  } // end of check that we have TrackingParticles to begin with...
332 
333 
334  // now that we have all of the track information for a given bunch crossing,
335  // make PileupSummary for this one and move on
336 
337  // std::cout << "Making PSI for bunch " << lastBunchCrossing << std::endl;
338 
339  if(!HaveTrackingParticles) { // stick in one value so we don't have empty vectors
340 
341  zpositions.push_back(-999.);
342  sumpT_lowpT.push_back(0.);
343  sumpT_highpT.push_back(0.);
344  ntrks_lowpT.push_back(0);
345  ntrks_highpT.push_back(0);
346 
347  }
348 
349  if(Have_pThats) {
350 
352  (*InteractionsIter),
353  (*zPosIter),
354  (*tPosIter),
355  sumpT_lowpT,
356  sumpT_highpT,
357  ntrks_lowpT,
358  ntrks_highpT,
359  (*TEventInfoIter),
360  (*pThatIter),
361  (*BXIter),
362  (*TInteractionsIter),
363  bunchSpacing
364  );
365  PSIVector->push_back(PSI_bunch);
366 
367  zPosIter++;
368  tPosIter++;
369  pThatIter++;
370  }
371  else{
372 
373  std::vector<float> zposZeros( (*TEventInfoIter).size(), 0);
374  std::vector<float> tposZeros( (*TEventInfoIter).size(), 0);
375  std::vector<float> pThatZeros( (*TEventInfoIter).size(), 0);
376 
378  (*InteractionsIter),
379  zposZeros,
380  tposZeros,
381  sumpT_lowpT,
382  sumpT_highpT,
383  ntrks_lowpT,
384  ntrks_highpT,
385  (*TEventInfoIter),
386  pThatZeros,
387  (*BXIter),
388  (*TInteractionsIter),
389  bunchSpacing
390  );
391 
392  PSIVector->push_back(PSI_bunch);
393  }
394  //std::cout << " " << std::endl;
395  //std::cout << "Adding Bunch Crossing, nint " << (*BXIter) << " " << (*InteractionsIter) << std::endl;
396 
397  //for(int iv = 0; iv<(*InteractionsIter); ++iv){
398 
399  // std::cout << "Z position " << zpositions[iv] << std::endl;
400  // std::cout << "ntrks_lowpT " << ntrks_lowpT[iv] << std::endl;
401  // std::cout << "sumpT_lowpT " << sumpT_lowpT[iv] << std::endl;
402  // std::cout << "ntrks_highpT " << ntrks_highpT[iv] << std::endl;
403  // std::cout << "sumpT_highpT " << sumpT_highpT[iv] << std::endl;
404  //std::cout << iv << " " << PSI_bunch.getPU_EventID()[iv] << std::endl;
405  //}
406 
407 
408 
409  // if(HaveTrackingParticles) lastBunchCrossing = iVtx->eventId().bunchCrossing();
410 
411  event_index_.clear();
412  zpositions.clear();
413  sumpT_lowpT.clear();
414  sumpT_highpT.clear();
415  ntrks_lowpT.clear();
416  ntrks_highpT.clear();
417  nminb_vtx = 0;
418  lastEvent=0;
419 
420 
421  } // end of loop over bunch crossings
422 
423  // put our vector of PileupSummaryInfo objects into the event.
424 
425  event.put(std::move(PSIVector));
426 
427  //add bunch spacing to the event as a seperate integer for use by downstream modules
428  std::unique_ptr<int> bunchSpacingP(new int(bunchSpacing));
429  event.put(std::move(bunchSpacingP),"bunchSpacing");
430 
431 }
432 
433 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< PileupMixingContent > PileupInfoLabel_
edm::EDGetTokenT< std::vector< PileupSummaryInfo > > pileupSummaryToken_
std::string MessageCategory_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
const std::vector< float > & getMix_TrueInteractions() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const std::vector< int > & getMix_bunchCrossing() const
const std::vector< float > & getMix_pT_hats() const
edm::EDGetTokenT< int > bunchSpacingToken_
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
Definition: config.py:1
void produce(edm::Event &, const edm::EventSetup &) override
PileupInformation(const edm::ParameterSet &)
const std::vector< int > & getMix_Ninteractions() const
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< float > sumpT_highpT
const std::vector< edm::EventID > getMix_eventInfo() const
std::vector< float > zpositions
T const * product() const
Definition: Handle.h:81
std::vector< int > ntrks_lowpT
edm::EDGetTokenT< TrackingParticleCollection > trackingTruthT_
edm::EDGetTokenT< PileupVertexContent > PileupVtxLabel_
const std::vector< float > & getMix_z_Vtxs() const
std::vector< float > sumpT_lowpT
const std::vector< float > & getMix_t_Vtxs() const
const int & getMix_bunchSpacing() const
edm::EDGetTokenT< TrackingVertexCollection > trackingTruthV_
std::vector< int > ntrks_highpT
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1
ib
Definition: cuy.py:661