CMS 3D CMS Logo

DuplicationChecker.cc
Go to the documentation of this file.
1 /*class DuplicationChecker
2  *
3  * Class to monitor duplication of events
4  *
5  *
6  */
7 
10 using namespace edm;
11 
13  wmanager_(iPSet,consumesCollector()),
14  generatedCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
15  searchForLHE_(iPSet.getParameter<bool>("searchForLHE"))
16 {
17  if (searchForLHE_) {
18  lheEventProduct_ = iPSet.getParameter<edm::InputTag>("lheEventProduct");
19  }
20  xBjorkenHistory.clear();
21 
22  if (searchForLHE_) lheEventProductToken_=consumes<LHEEventProduct>(lheEventProduct_);
23  else generatedCollectionToken_=consumes<HepMCProduct>(generatedCollection_);
24 
25 }
26 
28 {
29  xBjorkenHistory.clear();
30 }
31 
34  DQMHelper dqm(&i); i.setCurrentFolder("Generator/DuplicationCheck");
35 
37  xBjorkenME = dqm.book1dHisto("xBjorkenME", "x Bjorken ratio", 1000000, 0., 1.);
38 }
39 
41 {
42 
43  double bjorken = 0;
44 
45  double weight = 1.;
46 
47  if (searchForLHE_) {
48 
50  iEvent.getByToken(lheEventProductToken_, evt);
51 
52  const lhef::HEPEUP hepeup_ = evt->hepeup();
53 
54  const std::vector<lhef::HEPEUP::FiveVector> pup_ = hepeup_.PUP;
55 
56  double pz1=(pup_[0])[3];
57  double pz2=(pup_[1])[3];
58  bjorken+=(pz1/(pz1+pz2));
59  }
60  else {
61  //change teh weight in this case
62  weight = wmanager_.weight(iEvent);
63 
66 
67  const HepMC::PdfInfo *pdf = evt->GetEvent()->pdf_info();
68  if(pdf){
69  bjorken = ((pdf->x1())/((pdf->x1())+(pdf->x2())));
70  }
71 
72  }
73 
74  xBjorkenHistory.insert(std::pair<double,edm::EventID>(bjorken,iEvent.id()));
75 
76  xBjorkenME->Fill(bjorken,weight);
77 
78 }//analyze
79 
81 {
82  associationMap::iterator itr;
83  associationMap::iterator lastElement;
84 
85  theObjects.clear();
86 
87  // locate an iterator to the first pair object associated with key
88  itr = mMap.find(key);
89  if (itr == mMap.end())
90  return; // no elements associated with key, so return immediately
91 
92  // get an iterator to the element that is one past the last element associated with key
93  lastElement = mMap.upper_bound(key);
94 
95  // for each element in the sequence [itr, lastElement)
96  for ( ; itr != lastElement; ++itr)
97  theObjects.push_back(itr);
98 }
99 
100 /* no corresponding function available in MultiThreaded version
101 void DuplicationChecker::endJob()
102 {
103 
104  itemList theObjects;
105  theObjects.reserve(10);
106 
107  for (associationMap::iterator it = xBjorkenHistory.begin(); it != xBjorkenHistory.end(); it++) {
108  double theKey = (*it).first;
109 
110  findValuesAssociatedWithKey(xBjorkenHistory, theKey, theObjects);
111 
112  if (theObjects.size() > 1) {
113  edm::LogWarning("DuplicatedEventFound") << "Duplicated events found with xBjorken = " << std::fixed << std::setw(16) << std::setprecision(14) << theKey;
114  for (unsigned int i = 0; i < theObjects.size(); i++) {
115  edm::LogPrint("DuplicatedEventList") << "Event = " << (*theObjects[i]).second;
116  }
117  }
118 
119  theObjects.clear();
120 
121  }
122 
123 }
124 */
T getParameter(std::string const &) const
const lhef::HEPEUP & hepeup() const
edm::EDGetTokenT< LHEEventProduct > lheEventProductToken_
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
edm::InputTag generatedCollection_
struct HEPEUP_ hepeup_
edm::EDGetTokenT< edm::HepMCProduct > generatedCollectionToken_
Definition: weight.py:1
std::vector< associationMap::iterator > itemList
MonitorElement * book1dHisto(std::string name, std::string title, int n, double xmin, double xmax, std::string xaxis, std::string yaxis)
Definition: DQMHelper.cc:12
associationMap xBjorkenHistory
void Fill(long long x)
virtual void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
int iEvent
Definition: GenABIO.cc:230
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
WeightManager wmanager_
std::multimap< double, edm::EventID > associationMap
MonitorElement * xBjorkenME
edm::InputTag lheEventProduct_
DuplicationChecker(const edm::ParameterSet &)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
double weight(const edm::Event &)
void findValuesAssociatedWithKey(associationMap &mMap, double &key, itemList &theObjects)
Definition: Run.h:43