CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/GeneratorInterface/LHEInterface/plugins/LHEProducer.cc

Go to the documentation of this file.
00001 #include <string>
00002 #include <memory>
00003 
00004 #include <boost/shared_ptr.hpp>
00005 #include <sigc++/signal.h>
00006 
00007 #include <HepMC/GenEvent.h>
00008 #include <HepMC/SimpleVector.h>
00009 
00010 #include "FWCore/Framework/interface/EDFilter.h"
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/Run.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 
00017 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00018 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00019 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00020 
00021 #include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
00022 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
00023 
00024 #include "GeneratorInterface/LHEInterface/interface/LHERunInfo.h"
00025 #include "GeneratorInterface/LHEInterface/interface/LHEEvent.h"
00026 #include "GeneratorInterface/LHEInterface/interface/Hadronisation.h"
00027 #include "GeneratorInterface/LHEInterface/interface/JetMatching.h"
00028 #include "GeneratorInterface/LHEInterface/interface/JetMatchingMLM.h"
00029 #include "GeneratorInterface/LHEInterface/interface/BranchingRatios.h"
00030 
00031 using namespace lhef;
00032 
00033 class LHEProducer : public edm::EDFilter {
00034     public:
00035         explicit LHEProducer(const edm::ParameterSet &params);
00036         virtual ~LHEProducer();
00037 
00038     protected:
00039   virtual void beginJob();
00040         virtual void endJob();
00041         virtual bool beginRun(edm::Run &run, const edm::EventSetup &es);
00042         virtual bool endRun(edm::Run &run, const edm::EventSetup &es);
00043         virtual bool filter(edm::Event &event, const edm::EventSetup &es);
00044 
00045     private:
00046         double matching(const HepMC::GenEvent *event, bool shower) const;
00047 
00048         bool showeredEvent(const boost::shared_ptr<HepMC::GenEvent> &event);
00049         void onInit();
00050         void onBeforeHadronisation();
00051 
00052         unsigned int                    eventsToPrint;
00053         std::vector<int>                removeResonances;
00054         std::auto_ptr<Hadronisation>    hadronisation;
00055         std::auto_ptr<JetMatching>      jetMatching;
00056 
00057         double                          extCrossSect;
00058         double                          extFilterEff;
00059         bool                            matchSummary;
00060 
00061         boost::shared_ptr<LHEEvent>     partonLevel;
00062         boost::shared_ptr<LHERunInfo>   runInfo;
00063         unsigned int                    index;
00064         bool                            matchingDone;
00065         double                          weight;
00066         BranchingRatios                 branchingRatios;
00067 };
00068 
00069 LHEProducer::LHEProducer(const edm::ParameterSet &params) :
00070         eventsToPrint(params.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
00071         extCrossSect(params.getUntrackedParameter<double>("crossSection", -1.0)),
00072         extFilterEff(params.getUntrackedParameter<double>("filterEfficiency", -1.0)),
00073         matchSummary(false)
00074 {
00075         hadronisation = Hadronisation::create(
00076                 params.getParameter<edm::ParameterSet>("hadronisation"));
00077 
00078         if (params.exists("removeResonances"))
00079                 removeResonances =
00080                         params.getParameter<std::vector<int> >(
00081                                                         "removeResonances");
00082 
00083         std::set<std::string> matchingCapabilities;
00084         if (params.exists("jetMatching")) {
00085                 edm::ParameterSet jetParams =
00086                         params.getUntrackedParameter<edm::ParameterSet>(
00087                                                                 "jetMatching");
00088                 jetMatching = JetMatching::create(jetParams);
00089 
00090                 matchingCapabilities = jetMatching->capabilities();
00091                 hadronisation->matchingCapabilities(matchingCapabilities);
00092         }
00093 
00094         produces<edm::HepMCProduct>();
00095         produces<GenEventInfoProduct>();
00096         produces<GenRunInfoProduct, edm::InRun>();
00097 
00098         if (jetMatching.get()) {
00099                 if (params.getUntrackedParameter<bool>(
00100                                         "preferShowerVetoCallback", true))
00101                         hadronisation->onShoweredEvent().connect(
00102                                 sigc::mem_fun(*this,
00103                                               &LHEProducer::showeredEvent));
00104                 hadronisation->onInit().connect(
00105                                 sigc::mem_fun(*this, &LHEProducer::onInit));
00106                 hadronisation->onBeforeHadronisation().connect(
00107                         sigc::mem_fun(*this,
00108                                       &LHEProducer::onBeforeHadronisation));
00109 
00110                 matchSummary = matchingCapabilities.count("matchSummary");
00111                 if (matchSummary) {
00112                         produces< std::vector<double> >("matchDeltaR");
00113                         produces< std::vector<double> >("matchDeltaPRel");
00114                 }
00115         }
00116 
00117         // force total branching ratio for QCD/QED to 1
00118         for(int i = 0; i < 6; i++)
00119                 branchingRatios.set(i);
00120         for(int i = 9; i < 23; i++)
00121                 branchingRatios.set(i);
00122 }
00123 
00124 LHEProducer::~LHEProducer()
00125 {
00126 }
00127 
00128 void LHEProducer::beginJob()
00129 {
00130         hadronisation->init();
00131 }
00132 
00133 void LHEProducer::endJob()
00134 {
00135         hadronisation.reset();
00136         jetMatching.reset();
00137 }
00138 
00139 bool LHEProducer::beginRun(edm::Run &run, const edm::EventSetup &es)
00140 {
00141         edm::Handle<LHERunInfoProduct> product;
00142         run.getByLabel("source", product);
00143 
00144         runInfo.reset(new LHERunInfo(*product));
00145         index = 0;
00146 
00147         return true;
00148 }
00149 
00150 bool LHEProducer::endRun(edm::Run &run, const edm::EventSetup &es)
00151 {
00152         hadronisation->statistics();
00153 
00154         LHERunInfo::XSec crossSection;
00155         if (runInfo) {
00156                 crossSection = runInfo->xsec();
00157                 runInfo->statistics();
00158         }
00159 
00160         std::auto_ptr<GenRunInfoProduct> runInfo(new GenRunInfoProduct);
00161 
00162         runInfo->setInternalXSec(
00163                         GenRunInfoProduct::XSec(crossSection.value,
00164                                                 crossSection.error));
00165         runInfo->setExternalXSecLO(extCrossSect);
00166         runInfo->setFilterEfficiency(extFilterEff);
00167 
00168         run.put(runInfo);
00169 
00170         runInfo.reset();
00171 
00172         return true;
00173 }
00174 
00175 bool LHEProducer::filter(edm::Event &event, const edm::EventSetup &es)
00176 {
00177         std::auto_ptr<edm::HepMCProduct> result(new edm::HepMCProduct);
00178 
00179         edm::Handle<LHEEventProduct> product;
00180         event.getByLabel("source", product);
00181 
00182         partonLevel.reset(new LHEEvent(runInfo, *product));
00183         if (!removeResonances.empty())
00184                 partonLevel->removeResonances(removeResonances);
00185 
00186         if (product->pdf())
00187                 partonLevel->setPDF(
00188                         std::auto_ptr<LHEEvent::PDF>(
00189                                 new LHEEvent::PDF(*product->pdf())));
00190 
00191         hadronisation->setEvent(partonLevel);
00192 
00193         double br = branchingRatios.getFactor(hadronisation.get(),
00194                                               partonLevel);
00195 
00196         matchingDone = false;
00197         weight = 1.0;
00198         std::auto_ptr<HepMC::GenEvent> hadronLevel(hadronisation->hadronize());
00199 
00200         if (!hadronLevel.get()) {
00201                 if (matchingDone) {
00202                         if (weight == 0.0)
00203                                 partonLevel->count(LHERunInfo::kSelected, br);
00204                         else
00205                                 partonLevel->count(LHERunInfo::kKilled,
00206                                                    br, weight);
00207                 } else
00208                         partonLevel->count(LHERunInfo::kTried, br);
00209         }
00210 
00211         if (!matchingDone && jetMatching.get() && hadronLevel.get())
00212                 weight = matching(hadronLevel.get(), false);
00213 
00214         if (weight == 0.0) {
00215                 edm::LogInfo("Generator|LHEInterface")
00216                         << "Event got rejected by the "
00217                            "jet matching." << std::endl;
00218 
00219                 if (hadronLevel.get()) {
00220                         partonLevel->count(LHERunInfo::kSelected);
00221                         hadronLevel.reset();
00222                 }
00223         }
00224 
00225         if (!hadronLevel.get()) {
00226                 event.put(result);
00227                 std::auto_ptr<GenEventInfoProduct> info(
00228                                                 new GenEventInfoProduct);
00229                 event.put(info);
00230                 return false;
00231         }
00232 
00233         partonLevel->count(LHERunInfo::kAccepted, br, weight);
00234 
00235         hadronLevel->set_event_number(++index);
00236 
00237         if (eventsToPrint) {
00238                 eventsToPrint--;
00239                 hadronLevel->print();
00240         }
00241 
00242         std::auto_ptr<GenEventInfoProduct> info(
00243                                 new GenEventInfoProduct(hadronLevel.get()));
00244         result->addHepMCData(hadronLevel.release());
00245         event.put(result);
00246         event.put(info);
00247 
00248         if (jetMatching.get() && matchSummary) {
00249                 std::auto_ptr< std::vector<double> > matchDeltaR(
00250                                                 new std::vector<double>);
00251                 std::auto_ptr< std::vector<double> > matchDeltaPRel(
00252                                                 new std::vector<double>);
00253 
00254                 typedef std::vector<JetMatching::JetPartonMatch> Matches;
00255                 Matches matches = jetMatching->getMatchSummary();
00256 
00257                 for(Matches::const_iterator iter = matches.begin();
00258                     iter != matches.end(); ++iter) {
00259                         if (!iter->isMatch())
00260                                 continue;
00261 
00262                         matchDeltaR->push_back(iter->delta);
00263                         matchDeltaPRel->push_back(iter->jet.rho() /
00264                                                   iter->parton.rho() - 1.0);
00265                 }
00266 
00267                 event.put(matchDeltaR, "matchDeltaR");
00268                 event.put(matchDeltaPRel, "matchDeltaPRel");
00269         }
00270 
00271         return true;
00272 }
00273 
00274 double LHEProducer::matching(const HepMC::GenEvent *event, bool shower) const
00275 {
00276         if (!jetMatching.get())
00277                 return 1.0;
00278 
00279         return jetMatching->match(partonLevel->asHepMCEvent().get(),
00280                                   event, shower);
00281 }
00282 
00283 bool LHEProducer::showeredEvent(const boost::shared_ptr<HepMC::GenEvent> &event)
00284 {
00285         weight = matching(event.get(), true);
00286         matchingDone = true;
00287         return weight == 0.0;
00288 }
00289 
00290 void LHEProducer::onInit()
00291 {
00292         jetMatching->init(runInfo);
00293 }
00294 
00295 void LHEProducer::onBeforeHadronisation()
00296 {
00297         jetMatching->beforeHadronisation(partonLevel);
00298 }
00299 
00300 DEFINE_FWK_MODULE(LHEProducer);
00301 
00302 DEFINE_LHE_JETMATCHING_PLUGIN(JetMatchingMLM);