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