CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/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: 2012/10/16 14:49:03 $
00006  *  $Revision: 1.4 $
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   if (searchForLHE_) {
00020     lheEventProduct_ = iPSet.getParameter<edm::InputTag>("lheEventProduct");
00021   }
00022   dbe = 0;
00023   dbe = edm::Service<DQMStore>().operator->();
00024 
00025   xBjorkenHistory.clear();
00026 }
00027 
00028 DuplicationChecker::~DuplicationChecker() 
00029 {
00030   xBjorkenHistory.clear();
00031 }
00032 
00033 void DuplicationChecker::beginJob()
00034 {
00035   if(dbe){
00037         dbe->setCurrentFolder("Generator/DuplicationCheck");
00038         
00040         xBjorkenME = dbe->book1D("xBjorkenME", "x Bjorken ratio", 1000000, 0., 1.);
00041   }
00042 }
00043 
00044 void DuplicationChecker::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup)
00045 {
00046     
00047   double bjorken = 0;
00048  
00049   double weight = 1.;
00050 
00051   if (searchForLHE_) {
00052 
00053     Handle<LHEEventProduct> evt;
00054     iEvent.getByLabel(lheEventProduct_, evt);
00055 
00056     const lhef::HEPEUP hepeup_ = evt->hepeup();
00057 
00058     const std::vector<lhef::HEPEUP::FiveVector> pup_ = hepeup_.PUP;
00059 
00060     double pz1=(pup_[0])[3];
00061     double pz2=(pup_[1])[3];
00062     bjorken+=(pz1/(pz1+pz2));
00063   }
00064   else {
00065     //change teh weight in this case
00066     weight = _wmanager.weight(iEvent);
00067 
00068     edm::Handle<HepMCProduct> evt;
00069     iEvent.getByLabel(generatedCollection_, evt);
00070 
00071     const HepMC::PdfInfo *pdf = evt->GetEvent()->pdf_info();    
00072     if(pdf){
00073       bjorken = ((pdf->x1())/((pdf->x1())+(pdf->x2())));
00074     }
00075 
00076   }
00077 
00078   xBjorkenHistory.insert(std::pair<double,edm::EventID>(bjorken,iEvent.id()));
00079 
00080   xBjorkenME->Fill(bjorken,weight);
00081 
00082 }//analyze
00083 
00084 void DuplicationChecker::findValuesAssociatedWithKey(associationMap &mMap, double &key, itemList &theObjects)
00085 {
00086   associationMap::iterator itr;
00087   associationMap::iterator lastElement;
00088         
00089   theObjects.clear();
00090 
00091   // locate an iterator to the first pair object associated with key
00092   itr = mMap.find(key);
00093   if (itr == mMap.end())
00094     return; // no elements associated with key, so return immediately
00095 
00096   // get an iterator to the element that is one past the last element associated with key
00097   lastElement = mMap.upper_bound(key);
00098 
00099   // for each element in the sequence [itr, lastElement)
00100   for ( ; itr != lastElement; ++itr)
00101     theObjects.push_back(itr);
00102 }  
00103 
00104 void DuplicationChecker::endJob()
00105 {
00106 
00107   itemList theObjects;
00108   theObjects.reserve(10);
00109 
00110   for (associationMap::iterator it = xBjorkenHistory.begin(); it != xBjorkenHistory.end(); it++) {
00111     double theKey = (*it).first;
00112 
00113     findValuesAssociatedWithKey(xBjorkenHistory, theKey, theObjects);
00114 
00115     if (theObjects.size() > 1) {
00116       edm::LogWarning("DuplicatedEventFound") << "Duplicated events found with xBjorken = " << std::fixed << std::setw(16) << std::setprecision(14) << theKey; 
00117       for (unsigned int i = 0; i < theObjects.size(); i++) {
00118         edm::LogPrint("DuplicatedEventList") << "Event = " << (*theObjects[i]).second;
00119       }
00120     }
00121 
00122     theObjects.clear();
00123  
00124   }
00125 
00126 }