CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DuplicationChecker.cc
Go to the documentation of this file.
1 /*class DuplicationChecker
2  *
3  * Class to monitor duplication of events
4  *
5  * $Date: 2010/05/19 13:05:17 $
6  * $Revision: 1.4 $
7  *
8  */
9 
11 
12 using namespace edm;
13 
15  generatedCollection_(iPSet.getParameter<edm::InputTag>("generatedCollection")),
16  searchForLHE_(iPSet.getParameter<bool>("searchForLHE"))
17 {
18  dbe = 0;
20 
21  xBjorkenHistory.clear();
22 }
23 
25 {
26  xBjorkenHistory.clear();
27 }
28 
30 {
31  if(dbe){
33  dbe->setCurrentFolder("Generator/DuplicationCheck");
34 
36  xBjorkenME = dbe->book1D("xBjorkenME", "x Bjorken ratio", 1000000, 0., 1.);
37  }
38 }
39 
41 {
42 
43  double bjorken = 0;
44 
45  if (searchForLHE_) {
46 
48  iEvent.getByType( evt );
49 
50  const lhef::HEPEUP hepeup_ = evt->hepeup();
51 
52  const std::vector<lhef::HEPEUP::FiveVector> pup_ = hepeup_.PUP;
53 
54  double pz1=(pup_[0])[3];
55  double pz2=(pup_[1])[3];
56  bjorken+=(pz1/(pz1+pz2));
57  }
58  else {
59 
61  iEvent.getByLabel(generatedCollection_, evt);
62 
63  const HepMC::PdfInfo *pdf = evt->GetEvent()->pdf_info();
64  if(pdf){
65  bjorken = ((pdf->x1())/((pdf->x1())+(pdf->x2())));
66  }
67 
68  }
69 
70  xBjorkenHistory.insert(std::pair<double,edm::EventID>(bjorken,iEvent.id()));
71 
72  xBjorkenME->Fill(bjorken);
73 
74 }//analyze
75 
77 {
78  associationMap::iterator itr;
79  associationMap::iterator lastElement;
80 
81  theObjects.clear();
82 
83  // locate an iterator to the first pair object associated with key
84  itr = mMap.find(key);
85  if (itr == mMap.end())
86  return; // no elements associated with key, so return immediately
87 
88  // get an iterator to the element that is one past the last element associated with key
89  lastElement = mMap.upper_bound(key);
90 
91  // for each element in the sequence [itr, lastElement)
92  for ( ; itr != lastElement; ++itr)
93  theObjects.push_back(itr);
94 }
95 
97 {
98 
99  itemList theObjects;
100  theObjects.reserve(10);
101 
102  for (associationMap::iterator it = xBjorkenHistory.begin(); it != xBjorkenHistory.end(); it++) {
103  double theKey = (*it).first;
104 
105  findValuesAssociatedWithKey(xBjorkenHistory, theKey, theObjects);
106 
107  if (theObjects.size() > 1) {
108  edm::LogWarning("DuplicatedEventFound") << "Duplicated events found with xBjorken = " << std::fixed << std::setw(16) << std::setprecision(14) << theKey;
109  for (unsigned int i = 0; i < theObjects.size(); i++) {
110  edm::LogPrint("DuplicatedEventList") << "Event = " << (*theObjects[i]).second;
111  }
112  }
113 
114  theObjects.clear();
115 
116  }
117 
118 }
int i
Definition: DBlmapReader.cc:9
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
edm::InputTag generatedCollection_
struct HEPEUP_ hepeup_
bool getByType(Handle< PROD > &result) const
Definition: Event.h:397
std::vector< associationMap::iterator > itemList
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup)
associationMap xBjorkenHistory
void Fill(long long x)
U second(std::pair< T, U > const &p)
int iEvent
Definition: GenABIO.cc:243
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::multimap< double, edm::EventID > associationMap
MonitorElement * xBjorkenME
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
DuplicationChecker(const edm::ParameterSet &)
edm::EventID id() const
Definition: EventBase.h:56
list key
Definition: combine.py:13
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
void findValuesAssociatedWithKey(associationMap &mMap, double &key, itemList &theObjects)