CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/GeneratorInterface/RivetInterface/plugins/RivetAnalyzer.cc

Go to the documentation of this file.
00001 #include "GeneratorInterface/RivetInterface/interface/RivetAnalyzer.h"
00002 
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/Framework/interface/MakerMacros.h"
00006 
00007 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00008 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 
00011 #include "Rivet/AnalysisHandler.hh"
00012 #include "Rivet/Analysis.hh"
00013 #include "Rivet/RivetAIDA.hh"
00014 #include "LWH/AIManagedObject.h"
00015 
00016 #include <string>
00017 #include <vector>
00018 #include <iostream>
00019 #include <cstdlib>
00020 #include <cstring>
00021 
00022 using namespace Rivet;
00023 using namespace edm;
00024 
00025 RivetAnalyzer::RivetAnalyzer(const edm::ParameterSet& pset) : 
00026 _analysisHandler(),
00027 _isFirstEvent(true),
00028 _outFileName(pset.getParameter<std::string>("OutputFile")),
00029 //decide whether to finlaize tthe plots or not.
00030 //deciding not to finalize them can be useful for further harvesting of many jobs
00031 _doFinalize(pset.getParameter<bool>("DoFinalize")),
00032 _produceDQM(pset.getParameter<bool>("ProduceDQMOutput"))
00033 {
00034   //retrive the analysis name from paarmeter set
00035   std::vector<std::string> analysisNames = pset.getParameter<std::vector<std::string> >("AnalysisNames");
00036   
00037   _hepmcCollection = pset.getParameter<edm::InputTag>("HepMCCollection");
00038 
00039   _useExternalWeight = pset.getParameter<bool>("UseExternalWeight");
00040   if (_useExternalWeight) {
00041     if (!pset.exists("GenEventInfoCollection")){
00042       throw cms::Exception("RivetAnalyzer") << "when using an external event weight you have to specify the GenEventInfoProduct collection from which the weight has to be taken " ; 
00043     }
00044     _genEventInfoCollection = pset.getParameter<edm::InputTag>("GenEventInfoCollection");
00045   }
00046 
00047   //get the analyses
00048   _analysisHandler.addAnalyses(analysisNames);
00049 
00050   //go through the analyses and check those that need the cross section
00051   const std::set< AnaHandle, AnaHandleLess > & analyses = _analysisHandler.analyses();
00052 
00053   std::set< AnaHandle, AnaHandleLess >::const_iterator ibeg = analyses.begin();
00054   std::set< AnaHandle, AnaHandleLess >::const_iterator iend = analyses.end();
00055   std::set< AnaHandle, AnaHandleLess >::const_iterator iana; 
00056   double xsection = -1.;
00057   xsection = pset.getParameter<double>("CrossSection");
00058   for (iana = ibeg; iana != iend; ++iana){
00059     if ((*iana)->needsCrossSection())
00060       (*iana)->setCrossSection(xsection);
00061   }
00062   if (_produceDQM){
00063     // book stuff needed for DQM
00064     dbe = 0;
00065     dbe = edm::Service<DQMStore>().operator->();
00066     dbe->setVerbose(50);
00067   }  
00068 
00069 }
00070 
00071 RivetAnalyzer::~RivetAnalyzer(){
00072 }
00073 
00074 void RivetAnalyzer::beginJob(){
00075   //set the environment, very ugly but rivet is monolithic when it comes to paths
00076   char * cmsswbase    = getenv("CMSSW_BASE");
00077   char * cmsswrelease = getenv("CMSSW_RELEASE_BASE");
00078   std::string rivetref, rivetinfo;
00079   rivetref = "RIVET_REF_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
00080   rivetinfo = "RIVET_INFO_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
00081   putenv(strdup(rivetref.c_str()));
00082   putenv(strdup(rivetinfo.c_str()));
00083 }
00084 
00085 void RivetAnalyzer::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup){
00086   return;
00087 }
00088 
00089 void RivetAnalyzer::analyze(const edm::Event& iEvent,const edm::EventSetup& iSetup){
00090   
00091   //get the hepmc product from the event
00092   edm::Handle<HepMCProduct> evt;
00093   iEvent.getByLabel(_hepmcCollection, evt);
00094 
00095   // get HepMC GenEvent
00096   const HepMC::GenEvent *myGenEvent = evt->GetEvent();
00097   //if you want to use an external weight we have to clene the GenEvent and change the weight  
00098   if ( _useExternalWeight ){
00099     HepMC::GenEvent * tmpGenEvtPtr = new HepMC::GenEvent( *(evt->GetEvent()) );
00100     if (tmpGenEvtPtr->weights().size() == 0) {
00101       throw cms::Exception("RivetAnalyzer") << "Original weight container has 0 size ";
00102     }
00103     if (tmpGenEvtPtr->weights().size() > 1) {
00104       edm::LogWarning("RivetAnalyzer") << "Original event weight size is " << tmpGenEvtPtr->weights().size() << ". Will change only the first one ";  
00105     }
00106     edm::Handle<GenEventInfoProduct> genEventInfoProduct;
00107     iEvent.getByLabel(_genEventInfoCollection, genEventInfoProduct);
00108     tmpGenEvtPtr->weights()[0] = genEventInfoProduct->weight();
00109     myGenEvent = tmpGenEvtPtr; 
00110   }
00111     
00112 
00113   //aaply the beams initialization on the first event
00114   if (_isFirstEvent){
00115     _analysisHandler.init(*myGenEvent);
00116     _isFirstEvent = false;
00117   }
00118 
00119   //run the analysis
00120   _analysisHandler.analyze(*myGenEvent);
00121 
00122   //if we have cloned the GenEvent, we delete it
00123   if ( _useExternalWeight ) 
00124     delete myGenEvent;
00125 }
00126 
00127 
00128 void RivetAnalyzer::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){
00129   if (_doFinalize)
00130     _analysisHandler.finalize();
00131   else {
00132     //if we don't finalize we just want to do the transformation from histograms to DPS
00133     normalizeTree(_analysisHandler.tree());
00134   }
00135   _analysisHandler.writeData(_outFileName);
00136 
00137   return;
00138 }
00139 
00140 void RivetAnalyzer::endJob(){
00141 }
00142 
00143 void RivetAnalyzer::normalizeTree(AIDA::ITree& tree)    {
00144   using namespace AIDA;
00145   std::vector<string> analyses = _analysisHandler.analysisNames();
00146   tree.ls(".", true);
00147   const string tmpdir = "/RivetNormalizeTmp";
00148   tree.mkdir(tmpdir);
00149   foreach (const string& analysis, analyses) {
00150     if (_produceDQM){
00151       dbe->setCurrentFolder(("Rivet/"+analysis).c_str());
00152       //global variables that are always present
00153       //sumOfWeights
00154       TH1F nevent("nEvt", "n analyzed Events", 1, 0., 1.);
00155       nevent.SetBinContent(1,_analysisHandler.sumOfWeights());
00156       _mes.push_back(dbe->book1D("nEvt",&nevent));
00157     }  
00158     //cross section
00159     //TH1F xsection("xSection", "Cross Section", 1, 0., 1.);
00160     //xsection.SetBinContent(1,_analysisHandler.crossSection());
00161     //_mes.push_back(dbe->book1D("xSection",&xsection)); 
00162     //now loop over the histograms
00163     const vector<string> paths = tree.listObjectNames("/"+analysis, true); // args set recursive listing
00164     std::cout << "Number of objects in AIDA tree for analysis " << analysis << " = " << paths.size() << std::endl;
00165     foreach (const string& path, paths) {
00166       IManagedObject* hobj = tree.find(path);
00167       if (hobj) {
00168         // Weird seg fault on SLC4 when trying to dyn cast an IProfile ptr to a IHistogram
00169         // Fix by attempting to cast to IProfile first, only try IHistogram if it fails.
00170         IHistogram1D* histo = 0;
00171         IProfile1D* prof = dynamic_cast<IProfile1D*>(hobj);
00172         if (!prof) histo = dynamic_cast<IHistogram1D*>(hobj);
00173 
00174         std::cout << "Converting histo " << path << " to DPS" << std::endl;
00175         tree.mv(path, tmpdir);
00176         const size_t lastslash = path.find_last_of("/");
00177         const string basename = path.substr(lastslash+1, path.length() - (lastslash+1));
00178         const string tmppath = tmpdir + "/" + basename;
00179 
00180         // If it's a normal histo:
00181         if (histo) {
00182           IHistogram1D* tmphisto = dynamic_cast<IHistogram1D*>(tree.find(tmppath));
00183           if (tmphisto) {
00184             _analysisHandler.datapointsetFactory().create(path, *tmphisto);
00185           }
00186           //now convert to root and then ME
00187           TH1F* h = aida2root<IHistogram1D, TH1F>(histo, basename);
00188           if (_produceDQM)
00189             _mes.push_back(dbe->book1D(h->GetName(), h));
00190           delete h;
00191           tree.rm(tmppath);
00192         }
00193         // If it's a profile histo:
00194         else if (prof) {
00195           IProfile1D* tmpprof = dynamic_cast<IProfile1D*>(tree.find(tmppath));
00196           if (tmpprof) {
00197             _analysisHandler.datapointsetFactory().create(path, *tmpprof);
00198           }
00199           //now convert to root and then ME
00200           TProfile* p = aida2root<IProfile1D, TProfile>(prof, basename);
00201           if (_produceDQM)
00202             _mes.push_back(dbe->bookProfile(p->GetName(), p));
00203           delete p;
00204           tree.rm(tmppath);
00205         }
00206       }
00207     }
00208   }
00209   tree.rmdir(tmpdir);  
00210 }
00211 
00212 
00213 DEFINE_FWK_MODULE(RivetAnalyzer);