Go to the documentation of this file.00001 #include <algorithm>
00002 #include <functional>
00003 #include <iostream>
00004 #include <string>
00005 #include <memory>
00006
00007 #include <boost/bind.hpp>
00008 #include <boost/ptr_container/ptr_deque.hpp>
00009
00010 #include "FWCore/Framework/interface/InputSourceMacros.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 "DataFormats/Common/interface/OrphanHandle.h"
00018
00019 #include "SimDataFormats/GeneratorProducts/interface/LesHouches.h"
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/LHEReader.h"
00026
00027 #include "LHESource.h"
00028
00029 using namespace lhef;
00030
00031 LHESource::LHESource(const edm::ParameterSet ¶ms,
00032 const edm::InputSourceDescription &desc) :
00033 ProducerSourceFromFiles(params, desc, false),
00034 reader(new LHEReader(fileNames(), params.getUntrackedParameter<unsigned int>("skipEvents", 0))),
00035 wasMerged(false)
00036 {
00037 produces<LHEEventProduct>();
00038 produces<LHERunInfoProduct, edm::InRun>();
00039 }
00040
00041 LHESource::~LHESource()
00042 {
00043 }
00044
00045 void LHESource::endJob()
00046 {
00047 reader.reset();
00048 }
00049
00050 void LHESource::nextEvent()
00051 {
00052 if (partonLevel)
00053 return;
00054
00055 bool newFileOpened = false;
00056 partonLevel = reader->next(&newFileOpened);
00057 if(newFileOpened) incrementFileIndex();
00058 if (!partonLevel)
00059 return;
00060
00061 boost::shared_ptr<LHERunInfo> runInfoThis = partonLevel->getRunInfo();
00062 if (runInfoThis != runInfoLast) {
00063 runInfo = runInfoThis;
00064 runInfoLast = runInfoThis;
00065 }
00066 }
00067
00068 void LHESource::beginRun(edm::Run &run)
00069 {
00070 nextEvent();
00071 if (runInfoLast) {
00072 runInfo = runInfoLast;
00073
00074 std::auto_ptr<LHERunInfoProduct> product(
00075 new LHERunInfoProduct(*runInfo->getHEPRUP()));
00076 std::for_each(runInfo->getHeaders().begin(),
00077 runInfo->getHeaders().end(),
00078 boost::bind(
00079 &LHERunInfoProduct::addHeader,
00080 product.get(), _1));
00081 std::for_each(runInfo->getComments().begin(),
00082 runInfo->getComments().end(),
00083 boost::bind(&LHERunInfoProduct::addComment,
00084 product.get(), _1));
00085
00086
00087 runInfoProducts.push_back(new LHERunInfoProduct(*product));
00088 wasMerged = false;
00089
00090 run.put(product);
00091
00092 runInfo.reset();
00093 }
00094 }
00095
00096 void LHESource::endRun(edm::Run &run)
00097 {
00098 if (!runInfoProducts.empty()) {
00099 std::auto_ptr<LHERunInfoProduct> product(
00100 runInfoProducts.pop_front().release());
00101 run.put(product);
00102 }
00103 }
00104
00105 bool LHESource::setRunAndEventInfo(edm::EventID&, edm::TimeValue_t&)
00106 {
00107 nextEvent();
00108 if (!partonLevel)
00109 return false;
00110 return true;
00111 }
00112
00113 void LHESource::produce(edm::Event &event)
00114 {
00115 std::auto_ptr<LHEEventProduct> product(
00116 new LHEEventProduct(*partonLevel->getHEPEUP()));
00117 if (partonLevel->getPDF())
00118 product->setPDF(*partonLevel->getPDF());
00119 std::for_each(partonLevel->getComments().begin(),
00120 partonLevel->getComments().end(),
00121 boost::bind(&LHEEventProduct::addComment,
00122 product.get(), _1));
00123 event.put(product);
00124
00125 if (runInfo) {
00126 std::auto_ptr<LHERunInfoProduct> product(
00127 new LHERunInfoProduct(*runInfo->getHEPRUP()));
00128 std::for_each(runInfo->getHeaders().begin(),
00129 runInfo->getHeaders().end(),
00130 boost::bind(
00131 &LHERunInfoProduct::addHeader,
00132 product.get(), _1));
00133 std::for_each(runInfo->getComments().begin(),
00134 runInfo->getComments().end(),
00135 boost::bind(&LHERunInfoProduct::addComment,
00136 product.get(), _1));
00137
00138 if (!runInfoProducts.empty()) {
00139 runInfoProducts.front().mergeProduct(*product);
00140 if (!wasMerged) {
00141 runInfoProducts.pop_front();
00142 runInfoProducts.push_front(product);
00143 wasMerged = true;
00144 }
00145 }
00146
00147 runInfo.reset();
00148 }
00149
00150 partonLevel.reset();
00151 }
00152
00153 DEFINE_FWK_INPUT_SOURCE(LHESource);