CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
PileupInformation Class Reference

#include <PileupInformation.h>

Inheritance diagram for PileupInformation:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 PileupInformation (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

typedef std::map
< EncodedEventId, unsigned int > 
EncodedEventIdToIndex
 
typedef std::map< int, int > myindex
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &)
 

Private Attributes

edm::ParameterSet conf_
 
double distanceCut_
 
myindex event_index_
 
std::string MessageCategory_
 
std::vector< int > ntrks_highpT
 
std::vector< int > ntrks_lowpT
 
edm::InputTag PileupInfoLabel_
 
double pTcut_1_
 
double pTcut_2_
 
std::string simHitLabel_
 
std::auto_ptr< MixCollection
< SimTrack > > 
simTracks_
 
std::auto_ptr< MixCollection
< SimVertex > > 
simVertexes_
 
std::vector< float > sumpT_highpT
 
std::vector< float > sumpT_lowpT
 
edm::InputTag trackingTruth_
 
double volumeRadius_
 
double volumeZ_
 
std::vector< float > zpositions
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 29 of file PileupInformation.h.

Member Typedef Documentation

typedef std::map<EncodedEventId, unsigned int> PileupInformation::EncodedEventIdToIndex
private

Definition at line 42 of file PileupInformation.h.

typedef std::map< int, int > PileupInformation::myindex
private

Definition at line 43 of file PileupInformation.h.

Constructor & Destructor Documentation

PileupInformation::PileupInformation ( const edm::ParameterSet config)
explicit

Definition at line 20 of file PileupInformation.cc.

References distanceCut_, edm::ParameterSet::getParameter(), MessageCategory_, PileupInfoLabel_, pTcut_1_, pTcut_2_, simHitLabel_, AlCaHLTBitMon_QueryRunRegistry::string, trackingTruth_, volumeRadius_, and volumeZ_.

21 {
22  // Initialize global parameters
23 
24  pTcut_1_ = 0.1;
25  pTcut_2_ = 0.5; // defaults
26  distanceCut_ = config.getParameter<double>("vertexDistanceCut");
27  volumeRadius_ = config.getParameter<double>("volumeRadius");
28  volumeZ_ = config.getParameter<double>("volumeZ");
29  pTcut_1_ = config.getParameter<double>("pTcut_1");
30  pTcut_2_ = config.getParameter<double>("pTcut_2");
31 
32  PileupInfoLabel_ = config.getParameter<std::string>("PileupMixingLabel");
33 
34  trackingTruth_ = config.getParameter<std::string>("TrackingParticlesLabel");
35 
36  simHitLabel_ = config.getParameter<std::string>("simHitLabel");
37 
38  MessageCategory_ = "PileupInformation";
39 
40  edm::LogInfo (MessageCategory_) << "Setting up PileupInformation";
41  edm::LogInfo (MessageCategory_) << "Vertex distance cut set to " << distanceCut_ << " mm";
42  edm::LogInfo (MessageCategory_) << "Volume radius set to " << volumeRadius_ << " mm";
43  edm::LogInfo (MessageCategory_) << "Volume Z set to " << volumeZ_ << " mm";
44  edm::LogInfo (MessageCategory_) << "Lower pT Threshold set to " << pTcut_1_ << " GeV";
45  edm::LogInfo (MessageCategory_) << "Upper pT Threshold set to " << pTcut_2_ << " GeV";
46 
47 
48  produces< std::vector<PileupSummaryInfo> >();
49  //produces<PileupSummaryInfo>();
50 }
T getParameter(std::string const &) const
std::string MessageCategory_
edm::InputTag PileupInfoLabel_
edm::InputTag trackingTruth_
std::string simHitLabel_

Member Function Documentation

void PileupInformation::produce ( edm::Event event,
const edm::EventSetup setup 
)
privatevirtual

Implements edm::EDProducer.

Definition at line 53 of file PileupInformation.cc.

References EncodedEventId::bunchCrossing(), event_index_, edm::Event::getByLabel(), PileupMixingContent::getMix_bunchCrossing(), PileupMixingContent::getMix_Ninteractions(), PileupMixingContent::getMix_TrueInteractions(), cuy::ib, getHLTprescales::index, ntrks_highpT, ntrks_lowpT, PileupInfoLabel_, edm::Handle< T >::product(), pTcut_1_, pTcut_2_, simHitLabel_, simVertexes_, mathSSE::sqrt(), sumpT_highpT, sumpT_lowpT, trackingTruth_, and zpositions.

54 {
55 
56  std::auto_ptr<std::vector<PileupSummaryInfo> > PSIVector(new std::vector<PileupSummaryInfo>);
57 
58  edm::Handle< PileupMixingContent > MixingPileup; // Get True pileup information from MixingModule
59  event.getByLabel(PileupInfoLabel_, MixingPileup);
60 
61  std::vector<int> BunchCrossings;
62  std::vector<int> Interactions_Xing;
63  std::vector<float> TrueInteractions_Xing;
64 
65  const PileupMixingContent* MixInfo = MixingPileup.product();
66 
67  if(MixInfo) { // extract information - way easier than counting vertices
68 
69  const std::vector<int> bunchCrossing = MixInfo->getMix_bunchCrossing();
70  const std::vector<int> interactions = MixInfo->getMix_Ninteractions();
71  const std::vector<float> TrueInteractions = MixInfo->getMix_TrueInteractions();
72 
73 
74  for(int ib=0; ib<(int)bunchCrossing.size(); ++ib){
75  // std::cout << " bcr, nint " << bunchCrossing[ib] << " " << interactions[ib] << std::endl;
76  BunchCrossings.push_back(bunchCrossing[ib]);
77  Interactions_Xing.push_back(interactions[ib]);
78  TrueInteractions_Xing.push_back(TrueInteractions[ib]);
79  }
80  }
81  else{ //If no Mixing Truth, work with SimVertices (probably should throw an exception, but...)
82 
83  // Collect all the simvertex from the crossing frame
85  if( event.getByLabel("mix", simHitLabel_, cfSimVertexes) ) {
86 
87  // Create a mix collection from one simvertex collection
88  simVertexes_ = std::auto_ptr<MixCollection<SimVertex> >( new MixCollection<SimVertex>(cfSimVertexes.product()) );
89 
90  int index = 0;
91  // Solution to the problem of not having vertexId
92  // bool FirstL = true;
93  EncodedEventIdToIndex vertexId;
94  EncodedEventId oldEventId;
95  int oldBX = -1000;
96 
97  std::vector<int> BunchCrossings2;
98  std::list<int> Interactions_Xing2;
99 
100 
101  // Loop for finding repeated vertexId (vertexId problem hack)
102  for (MixCollection<SimVertex>::MixItr iterator = simVertexes_->begin(); iterator != simVertexes_->end(); ++iterator, ++index)
103  {
104  // std::cout << " SimVtx eventid, vertexid " << iterator->eventId().event() << " " << iterator->eventId().bunchCrossing() << std::endl;
105  if (!index || iterator->eventId() != oldEventId)
106  {
107  if(iterator->eventId().bunchCrossing()==0 && iterator->eventId().event()==0){
108  continue;
109  }
110  if(iterator->eventId().bunchCrossing() != oldBX) {
111  BunchCrossings2.push_back(iterator->eventId().bunchCrossing());
112  Interactions_Xing2.push_back(iterator->eventId().event());
113  oldBX = iterator->eventId().bunchCrossing();
114  }
115  else { Interactions_Xing2.pop_back();
116  Interactions_Xing2.push_back(iterator->eventId().event());
117  }
118 
119 
120  oldEventId = iterator->eventId();
121  continue;
122  }
123 
124  }
125 
126  std::vector<int>::iterator viter;
127  std::list<int>::iterator liter = Interactions_Xing2.begin();
128 
129  for(viter = BunchCrossings2.begin(); viter != BunchCrossings2.end(); ++viter, ++liter){
130  //std::cout << " bcr, nint from VTX " << (*viter) << " " << (*liter) << std::endl;
131  BunchCrossings.push_back((*viter));
132  Interactions_Xing.push_back((*liter));
133  TrueInteractions_Xing.push_back(-1.); // no idea what the true number is
134  }
135 
136  } // end of "did we find vertices?"
137  } // end of look at SimVertices
138 
139 
140  //Now, get information on valid particles that look like they could be in the tracking volume
141 
142 
143  zpositions.clear();
144  sumpT_lowpT.clear();
145  sumpT_highpT.clear();
146  ntrks_lowpT.clear();
147  ntrks_highpT.clear();
148  event_index_.clear();
149 
150  int lastEvent = 0; // zero is the true MC hard-scatter event
151 
152  // int lastBunchCrossing = 0; // 0 is the true bunch crossing, should always come first.
153 
156 
157  bool HaveTrackingParticles = false;
158 
159  TrackingVertexCollection::const_iterator iVtx;
160  TrackingVertexCollection::const_iterator iVtxTest;
161  TrackingParticleCollection::const_iterator iTrackTest;
162 
163 
164  if(event.getByLabel(trackingTruth_, mergedPH) && event.getByLabel(trackingTruth_, mergedVH)) {
165 
166  HaveTrackingParticles = true;
167 
168  iVtxTest = mergedVH->begin();
169  iTrackTest = mergedPH->begin();
170  }
171 
172  int nminb_vtx = 0;
173  // bool First = true;
174  // bool flag_new = false;
175 
176  std::vector<int>::iterator BXIter;
177  std::vector<int>::iterator InteractionsIter = Interactions_Xing.begin();
178  std::vector<float>::iterator TInteractionsIter = TrueInteractions_Xing.begin();
179 
180  // loop over the bunch crossings and interactions we have extracted
181 
182  for( BXIter = BunchCrossings.begin(); BXIter != BunchCrossings.end(); ++BXIter, ++InteractionsIter, ++TInteractionsIter) {
183 
184  //std::cout << "looking for BX: " << (*BXIter) << std::endl;
185 
186  if(HaveTrackingParticles) { // leave open the option of not producing TrackingParticles and just keeping # interactions
187 
188  for (iVtx = iVtxTest; iVtx != mergedVH->end(); ++iVtx) {
189 
190  if(iVtx->eventId().bunchCrossing() == (*BXIter) ) { // found first vertex in this bunch crossing
191 
192  if(iVtx->eventId().event() != lastEvent) {
193 
194  //std::cout << "BX,event " << iVtx->eventId().bunchCrossing() << " " << iVtx->eventId().event() << std::endl;
195 
196  float zpos = 0.;
197  zpos = iVtx->position().z();
198  zpositions.push_back(zpos); //save z position of each vertex
199  sumpT_lowpT.push_back(0.);
200  sumpT_highpT.push_back(0.);
201  ntrks_lowpT.push_back(0);
202  ntrks_highpT.push_back(0);
203 
204  lastEvent = iVtx->eventId().event();
205  iVtxTest = --iVtx; // just for security
206 
207  // turns out events aren't sequential... save map of indices
208 
209  event_index_.insert(myindex::value_type(lastEvent,nminb_vtx));
210 
211  ++nminb_vtx;
212 
213  continue;
214  }
215  }
216  }
217 
218  // next loop over tracks to get information
219 
220  for (TrackingParticleCollection::const_iterator iTrack = iTrackTest; iTrack != mergedPH->end(); ++iTrack)
221  {
222  bool FoundTrk = false;
223 
224  float zpos=0.;
225 
226  if(iTrack->eventId().bunchCrossing() == (*BXIter) && iTrack->eventId().event() > 0 )
227  {
228  FoundTrk = true;
229  int correct_index = event_index_[iTrack->eventId().event()];
230 
231  //std::cout << " track index, correct index " << iTrack->eventId().event() << " " << correct_index << std::endl;
232 
233  zpos = zpositions[correct_index];
234  if(iTrack->matchedHit()>0) {
235  if(fabs(iTrack->parentVertex()->position().z()-zpos)<0.1) { //make sure track really comes from this vertex
236  //std::cout << *iTrack << std::endl;
237  float Tpx = iTrack->p4().px();
238  float Tpy = iTrack->p4().py();
239  float TpT = sqrt(Tpx*Tpx + Tpy*Tpy);
240  if( TpT>pTcut_1_ ) {
241  sumpT_lowpT[correct_index]+=TpT;
242  ++ntrks_lowpT[correct_index];
243  }
244  if( TpT>pTcut_2_ ){
245  sumpT_highpT[correct_index]+=TpT;
246  ++ntrks_highpT[correct_index];
247  }
248  }
249  }
250  }
251  else{
252  if(FoundTrk) {
253 
254  iTrackTest = --iTrack; // reset so we can start over next time
255  --iTrackTest; // just to be sure
256  break;
257  }
258 
259  }
260 
261  } // end of track loop
262 
263  } // end of check that we have TrackingParticles to begin with...
264 
265 
266  // now that we have all of the track information for a given bunch crossing,
267  // make PileupSummary for this one and move on
268 
269  // std::cout << "Making PSI for bunch " << lastBunchCrossing << std::endl;
270 
271  if(!HaveTrackingParticles) { // stick in one value so we don't have empty vectors
272 
273  zpositions.push_back(-999.);
274  sumpT_lowpT.push_back(0.);
275  sumpT_highpT.push_back(0.);
276  ntrks_lowpT.push_back(0);
277  ntrks_highpT.push_back(0);
278 
279  }
280 
282  (*InteractionsIter),
283  zpositions,
284  sumpT_lowpT,
285  sumpT_highpT,
286  ntrks_lowpT,
287  ntrks_highpT,
288  (*BXIter),
289  (*TInteractionsIter)
290  );
291 
292  //std::cout << " " << std::endl;
293  //std::cout << "Adding Bunch Crossing, nint " << (*BXIter) << " " << (*InteractionsIter) << std::endl;
294 
295  //for(int iv = 0; iv<(*InteractionsIter); ++iv){
296 
297  // std::cout << "Z position " << zpositions[iv] << std::endl;
298  // std::cout << "ntrks_lowpT " << ntrks_lowpT[iv] << std::endl;
299  // std::cout << "sumpT_lowpT " << sumpT_lowpT[iv] << std::endl;
300  // std::cout << "ntrks_highpT " << ntrks_highpT[iv] << std::endl;
301  // std::cout << "sumpT_highpT " << sumpT_highpT[iv] << std::endl;
302  //}
303 
304  PSIVector->push_back(PSI_bunch);
305 
306  // if(HaveTrackingParticles) lastBunchCrossing = iVtx->eventId().bunchCrossing();
307 
308  event_index_.clear();
309  zpositions.clear();
310  sumpT_lowpT.clear();
311  sumpT_highpT.clear();
312  ntrks_lowpT.clear();
313  ntrks_highpT.clear();
314  nminb_vtx = 0;
315  lastEvent=0;
316 
317 
318  } // end of loop over bunch crossings
319 
320  // put our vector of PileupSummaryInfo objects into the event.
321 
322  event.put(PSIVector);
323 
324 
325 }
int ib
Definition: cuy.py:660
const std::vector< float > & getMix_TrueInteractions() const
const std::vector< int > & getMix_bunchCrossing() const
std::map< EncodedEventId, unsigned int > EncodedEventIdToIndex
const std::vector< int > & getMix_Ninteractions() const
T sqrt(T t)
Definition: SSEVec.h:48
int bunchCrossing() const
get the detector field from this detid
edm::InputTag PileupInfoLabel_
std::vector< float > sumpT_highpT
Container::value_type value_type
std::vector< float > zpositions
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
std::auto_ptr< MixCollection< SimVertex > > simVertexes_
std::vector< int > ntrks_lowpT
T const * product() const
Definition: Handle.h:81
std::vector< float > sumpT_lowpT
edm::InputTag trackingTruth_
std::vector< int > ntrks_highpT
std::string simHitLabel_

Member Data Documentation

edm::ParameterSet PileupInformation::conf_
private

Definition at line 40 of file PileupInformation.h.

double PileupInformation::distanceCut_
private

Definition at line 53 of file PileupInformation.h.

Referenced by PileupInformation().

myindex PileupInformation::event_index_
private

Definition at line 44 of file PileupInformation.h.

Referenced by produce().

std::string PileupInformation::MessageCategory_
private

Definition at line 62 of file PileupInformation.h.

Referenced by PileupInformation().

std::vector<int> PileupInformation::ntrks_highpT
private

Definition at line 50 of file PileupInformation.h.

Referenced by produce().

std::vector<int> PileupInformation::ntrks_lowpT
private

Definition at line 49 of file PileupInformation.h.

Referenced by produce().

edm::InputTag PileupInformation::PileupInfoLabel_
private

Definition at line 60 of file PileupInformation.h.

Referenced by PileupInformation(), and produce().

double PileupInformation::pTcut_1_
private

Definition at line 56 of file PileupInformation.h.

Referenced by PileupInformation(), and produce().

double PileupInformation::pTcut_2_
private

Definition at line 57 of file PileupInformation.h.

Referenced by PileupInformation(), and produce().

std::string PileupInformation::simHitLabel_
private

Definition at line 63 of file PileupInformation.h.

Referenced by PileupInformation(), and produce().

std::auto_ptr<MixCollection<SimTrack> > PileupInformation::simTracks_
private

Definition at line 64 of file PileupInformation.h.

std::auto_ptr<MixCollection<SimVertex> > PileupInformation::simVertexes_
private

Definition at line 65 of file PileupInformation.h.

Referenced by produce().

std::vector<float> PileupInformation::sumpT_highpT
private

Definition at line 48 of file PileupInformation.h.

Referenced by produce().

std::vector<float> PileupInformation::sumpT_lowpT
private

Definition at line 47 of file PileupInformation.h.

Referenced by produce().

edm::InputTag PileupInformation::trackingTruth_
private

Definition at line 59 of file PileupInformation.h.

Referenced by PileupInformation(), and produce().

double PileupInformation::volumeRadius_
private

Definition at line 54 of file PileupInformation.h.

Referenced by PileupInformation().

double PileupInformation::volumeZ_
private

Definition at line 55 of file PileupInformation.h.

Referenced by PileupInformation().

std::vector<float> PileupInformation::zpositions
private

Definition at line 46 of file PileupInformation.h.

Referenced by produce().