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  if (searchForLHE_) {
17  lheEventProduct_ = iPSet.getParameter<edm::InputTag>("lheEventProduct");
18  }
19  xBjorkenHistory.clear();
20 
21  if (searchForLHE_)
22  lheEventProductToken_ = consumes<LHEEventProduct>(lheEventProduct_);
23  else
24  generatedCollectionToken_ = consumes<HepMCProduct>(generatedCollection_);
25 }
26 
28 
31  DQMHelper dqm(&i);
32  i.setCurrentFolder("Generator/DuplicationCheck");
33 
35  xBjorkenME = dqm.book1dHisto("xBjorkenME", "x Bjorken ratio", 1000000, 0., 1.);
36 }
37 
39  double bjorken = 0;
40 
41  double weight = 1.;
42 
43  if (searchForLHE_) {
45  iEvent.getByToken(lheEventProductToken_, evt);
46 
47  const lhef::HEPEUP hepeup_ = evt->hepeup();
48 
49  const std::vector<lhef::HEPEUP::FiveVector> pup_ = hepeup_.PUP;
50 
51  double pz1 = (pup_[0])[3];
52  double pz2 = (pup_[1])[3];
53  bjorken += (pz1 / (pz1 + pz2));
54  } else {
55  //change teh weight in this case
57 
59  iEvent.getByToken(generatedCollectionToken_, evt);
60 
61  const HepMC::PdfInfo *pdf = evt->GetEvent()->pdf_info();
62  if (pdf) {
63  bjorken = ((pdf->x1()) / ((pdf->x1()) + (pdf->x2())));
64  }
65  }
66 
67  xBjorkenHistory.insert(std::pair<double, edm::EventID>(bjorken, iEvent.id()));
68 
69  xBjorkenME->Fill(bjorken, weight);
70 
71 } //analyze
72 
74  associationMap::iterator itr;
75  associationMap::iterator lastElement;
76 
77  theObjects.clear();
78 
79  // locate an iterator to the first pair object associated with key
80  itr = mMap.find(key);
81  if (itr == mMap.end())
82  return; // no elements associated with key, so return immediately
83 
84  // get an iterator to the element that is one past the last element associated with key
85  lastElement = mMap.upper_bound(key);
86 
87  // for each element in the sequence [itr, lastElement)
88  for (; itr != lastElement; ++itr)
89  theObjects.push_back(itr);
90 }
91 
92 /* no corresponding function available in MultiThreaded version
93 void DuplicationChecker::endJob()
94 {
95 
96  itemList theObjects;
97  theObjects.reserve(10);
98 
99  for (associationMap::iterator it = xBjorkenHistory.begin(); it != xBjorkenHistory.end(); it++) {
100  double theKey = (*it).first;
101 
102  findValuesAssociatedWithKey(xBjorkenHistory, theKey, theObjects);
103 
104  if (theObjects.size() > 1) {
105  edm::LogWarning("DuplicatedEventFound") << "Duplicated events found with xBjorken = " << std::fixed << std::setw(16) << std::setprecision(14) << theKey;
106  for (unsigned int i = 0; i < theObjects.size(); i++) {
107  edm::LogPrint("DuplicatedEventList") << "Event = " << (*theObjects[i]).second;
108  }
109  }
110 
111  theObjects.clear();
112 
113  }
114 
115 }
116 */
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< LHEEventProduct > lheEventProductToken_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::InputTag generatedCollection_
struct HEPEUP_ hepeup_
edm::EDGetTokenT< edm::HepMCProduct > generatedCollectionToken_
Definition: weight.py:1
std::multimap< double, edm::EventID > associationMap
std::vector< associationMap::iterator > itemList
associationMap xBjorkenHistory
void Fill(long long x)
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
int iEvent
Definition: GenABIO.cc:224
WeightManager wmanager_
~DuplicationChecker() override
MonitorElement * xBjorkenME
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
edm::InputTag lheEventProduct_
DuplicationChecker(const edm::ParameterSet &)
HLT enums.
double weight(const edm::Event &)
Definition: DQMStore.h:18
void findValuesAssociatedWithKey(associationMap &mMap, double &key, itemList &theObjects)
Definition: Run.h:45
const lhef::HEPEUP & hepeup() const