CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/Validation/EventGenerator/plugins/DuplicationChecker.cc

Go to the documentation of this file.
00001 /*class DuplicationChecker
00002  *  
00003  *  Class to monitor duplication of events
00004  *
00005  *  $Date: 2011/12/29 10:53:11 $
00006  *  $Revision: 1.2 $
00007  *
00008  */
00009  
00010 #include "Validation/EventGenerator/interface/DuplicationChecker.h"
00011 
00012 using namespace edm;
00013 
00014 DuplicationChecker::DuplicationChecker(const edm::ParameterSet& iPSet):
00015   _wmanager(iPSet),
00016   generatedCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
00017   searchForLHE_(iPSet.getParameter<bool>("searchForLHE"))
00018 {    
00019   dbe = 0;
00020   dbe = edm::Service<DQMStore>().operator->();
00021 
00022   xBjorkenHistory.clear();
00023 }
00024 
00025 DuplicationChecker::~DuplicationChecker() 
00026 {
00027   xBjorkenHistory.clear();
00028 }
00029 
00030 void DuplicationChecker::beginJob()
00031 {
00032   if(dbe){
00034         dbe->setCurrentFolder("Generator/DuplicationCheck");
00035         
00037         xBjorkenME = dbe->book1D("xBjorkenME", "x Bjorken ratio", 1000000, 0., 1.);
00038   }
00039 }
00040 
00041 void DuplicationChecker::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00042 {
00043     
00044   double bjorken = 0;
00045   
00046   double weight = _wmanager.weight(iEvent);
00047 
00048   if (searchForLHE_) {
00049 
00050     Handle<LHEEventProduct> evt;
00051     iEvent.getByType( evt );
00052 
00053     const lhef::HEPEUP hepeup_ = evt->hepeup();
00054 
00055     const std::vector<lhef::HEPEUP::FiveVector> pup_ = hepeup_.PUP;
00056 
00057     double pz1=(pup_[0])[3];
00058     double pz2=(pup_[1])[3];
00059     bjorken+=(pz1/(pz1+pz2));
00060   }
00061   else {
00062 
00063     edm::Handle<HepMCProduct> evt;
00064     iEvent.getByLabel(generatedCollection_, evt);
00065 
00066     const HepMC::PdfInfo *pdf = evt->GetEvent()->pdf_info();    
00067     if(pdf){
00068       bjorken = ((pdf->x1())/((pdf->x1())+(pdf->x2())));
00069     }
00070 
00071   }
00072 
00073   xBjorkenHistory.insert(std::pair<double,edm::EventID>(bjorken,iEvent.id()));
00074 
00075   xBjorkenME->Fill(bjorken,weight);
00076 
00077 }//analyze
00078 
00079 void DuplicationChecker::findValuesAssociatedWithKey(associationMap &mMap, double &key, itemList &theObjects)
00080 {
00081   associationMap::iterator itr;
00082   associationMap::iterator lastElement;
00083         
00084   theObjects.clear();
00085 
00086   // locate an iterator to the first pair object associated with key
00087   itr = mMap.find(key);
00088   if (itr == mMap.end())
00089     return; // no elements associated with key, so return immediately
00090 
00091   // get an iterator to the element that is one past the last element associated with key
00092   lastElement = mMap.upper_bound(key);
00093 
00094   // for each element in the sequence [itr, lastElement)
00095   for ( ; itr != lastElement; ++itr)
00096     theObjects.push_back(itr);
00097 }  
00098 
00099 void DuplicationChecker::endJob()
00100 {
00101 
00102   itemList theObjects;
00103   theObjects.reserve(10);
00104 
00105   for (associationMap::iterator it = xBjorkenHistory.begin(); it != xBjorkenHistory.end(); it++) {
00106     double theKey = (*it).first;
00107 
00108     findValuesAssociatedWithKey(xBjorkenHistory, theKey, theObjects);
00109 
00110     if (theObjects.size() > 1) {
00111       edm::LogWarning("DuplicatedEventFound") << "Duplicated events found with xBjorken = " << std::fixed << std::setw(16) << std::setprecision(14) << theKey; 
00112       for (unsigned int i = 0; i < theObjects.size(); i++) {
00113         edm::LogPrint("DuplicatedEventList") << "Event = " << (*theObjects[i]).second;
00114       }
00115     }
00116 
00117     theObjects.clear();
00118  
00119   }
00120 
00121 }