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