00001 #include <iostream>
00002 #include <iomanip>
00003 #include <sstream>
00004 #include <fstream>
00005 #include <string>
00006 #include <memory>
00007 #include <cmath>
00008
00009 #include <boost/algorithm/string/classification.hpp>
00010 #include <boost/algorithm/string/split.hpp>
00011
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/Run.h"
00014 #include "FWCore/Framework/interface/MakerMacros.h"
00015 #include "FWCore/Framework/interface/InputSourceMacros.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/Utilities/interface/Exception.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 #include "FWCore/Sources/interface/ProducerSourceFromFiles.h"
00020
00021 #include "SimDataFormats/GeneratorProducts/interface/LesHouches.h"
00022 #include "SimDataFormats/GeneratorProducts/interface/LHECommonBlocks.h"
00023 #include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
00024 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
00025
00026 #include "GeneratorInterface/AlpgenInterface/interface/AlpgenHeader.h"
00027 #include "GeneratorInterface/AlpgenInterface/interface/AlpgenEventRecordFixes.h"
00028
00029 class AlpgenSource : public edm::ProducerSourceFromFiles {
00030 public:
00032 AlpgenSource(const edm::ParameterSet ¶ms,
00033 const edm::InputSourceDescription &desc);
00034
00036 virtual ~AlpgenSource();
00037
00038 private:
00039 virtual bool setRunAndEventInfo(edm::EventID&, edm::TimeValue_t&);
00040 virtual void produce(edm::Event &event);
00041 virtual void beginRun(edm::Run &run);
00042
00044 template<typename T>
00045 T getParameter(AlpgenHeader::Parameter index) const;
00047 template<typename T>
00048 T getParameter(AlpgenHeader::Parameter index, const T &defValue) const;
00049
00052 std::string slhaMassLine(int pdgId, AlpgenHeader::Masses mass,
00053 const std::string &comment) const;
00054
00058 unsigned int processID() const;
00059
00061 bool readAlpgenEvent(lhef::HEPEUP &hepeup);
00062
00064 std::string fileName_;
00065
00067 std::unique_ptr<std::ifstream> inputFile_;
00068
00070 unsigned long skipEvents_;
00071
00073 unsigned long nEvent_;
00074
00076 LHERunInfoProduct::Header lheAlpgenUnwParHeader;
00078 AlpgenHeader header;
00079
00080 std::unique_ptr<lhef::HEPEUP> hepeup_;
00081
00083 std::string extraHeaderFileName_;
00084
00086 std::string extraHeaderName_;
00087
00089 bool writeAlpgenWgtFile;
00090 bool writeAlpgenParFile;
00091 bool writeExtraHeader;
00092 };
00093
00094 AlpgenSource::AlpgenSource(const edm::ParameterSet ¶ms,
00095 const edm::InputSourceDescription &desc) :
00096 edm::ProducerSourceFromFiles(params, desc, false),
00097 skipEvents_(params.getUntrackedParameter<unsigned int>("skipEvents", 0)),
00098 nEvent_(0), lheAlpgenUnwParHeader("AlpgenUnwParFile"),
00099 extraHeaderFileName_(params.getUntrackedParameter<std::string>("extraHeaderFileName","")),
00100 extraHeaderName_(params.getUntrackedParameter<std::string>("extraHeaderName","")),
00101 writeAlpgenWgtFile(params.getUntrackedParameter<bool>("writeAlpgenWgtFile", true)),
00102 writeAlpgenParFile(params.getUntrackedParameter<bool>("writeAlpgenParFile", true)),
00103 writeExtraHeader(params.getUntrackedParameter<bool>("writeExtraHeader", false))
00104 {
00105 std::vector<std::string> allFileNames = fileNames();
00106
00107
00108 if (allFileNames.size() != 1)
00109 throw cms::Exception("Generator|AlpgenInterface")
00110 << "AlpgenSource needs exactly one file specified "
00111 "for now." << std::endl;
00112
00113 fileName_ = allFileNames[0];
00114
00115
00116 if (fileName_.find("file:") != 0)
00117 throw cms::Exception("Generator|AlpgenInterface") << "AlpgenSource only supports the file: scheme "
00118 "for now." << std::endl;
00119 fileName_.erase(0, 5);
00120
00121
00122
00123 std::ifstream reader((fileName_ + "_unw.par").c_str());
00124 if (!reader.good())
00125 throw cms::Exception("Generator|AlpgenInterface")
00126 << "AlpgenSource was unable to open the file \""
00127 << fileName_ << "_unw.par\"." << std::endl;
00128
00129
00130 char buffer[256];
00131 while(reader.getline(buffer, sizeof buffer))
00132 lheAlpgenUnwParHeader.addLine(std::string(buffer) + "\n");
00133
00134
00135
00136 if (!header.parse(lheAlpgenUnwParHeader.begin(),
00137 lheAlpgenUnwParHeader.end()))
00138 throw cms::Exception("Generator|AlpgenInterface")
00139 << "AlpgenSource was unable to parse the Alpgen "
00140 << "unweighted parameter file." << std::endl;
00141
00142
00143 produces<LHERunInfoProduct, edm::InRun>();
00144 produces<LHEEventProduct>();
00145 }
00146
00147 AlpgenSource::~AlpgenSource()
00148 {
00149 }
00150
00151 std::string AlpgenSource::slhaMassLine(int pdgId, AlpgenHeader::Masses mass,
00152 const std::string &comment) const
00153 {
00154 std::ostringstream ss;
00155 ss << std::setw(9) << pdgId << " " << std::scientific
00156 << std::setprecision(9) << header.masses[mass] << " # "
00157 << comment << std::endl;
00158 return ss.str();
00159 }
00160
00161 void AlpgenSource::beginRun(edm::Run &run)
00162 {
00163
00164
00165 lhef::HEPRUP heprup;
00166
00167
00168
00169 heprup.IDBMUP.first = 2212;
00170 switch(getParameter<int>(AlpgenHeader::ih2)) {
00171 case 1:
00172 heprup.IDBMUP.second = 2212;
00173 break;
00174 case -1:
00175 heprup.IDBMUP.second = -2212;
00176 break;
00177 default:
00178 throw cms::Exception("Generator|AlpgenInterface")
00179 << "AlpgenSource was unable to understand the ih2 "
00180 << "parameter." << std::endl;
00181 }
00182
00183
00184 heprup.EBMUP.second = heprup.EBMUP.first =
00185 getParameter<double>(AlpgenHeader::ebeam);
00186
00187
00188 heprup.PDFGUP.first = -1;
00189 heprup.PDFGUP.second = -1;
00190 heprup.PDFSUP.first = -1;
00191 heprup.PDFSUP.second = -1;
00192
00193
00194 heprup.IDWTUP = 3;
00195
00196
00197 heprup.resize(1);
00198
00199
00200 heprup.XSECUP[0] = header.xsec;
00201 heprup.XERRUP[0] = header.xsecErr;
00202
00203
00204 heprup.XMAXUP[0] = header.xsec;
00205
00206
00207 heprup.LPRUP[0] = processID();
00208
00209
00210 LHERunInfoProduct::Header comments;
00211 comments.addLine("\n");
00212 comments.addLine("\tExtracted by AlpgenInterface\n");
00213
00214
00215
00216 LHERunInfoProduct::Header slha("slha");
00217 slha.addLine("\n# SLHA header containing masses from Alpgen\n");
00218 slha.addLine("Block MASS # Mass spectrum (kinematic masses)\n");
00219 slha.addLine("# PDG Mass\n");
00220 slha.addLine(slhaMassLine(4, AlpgenHeader::mc, "charm pole mass"));
00221 slha.addLine(slhaMassLine(5, AlpgenHeader::mb, "bottom pole mass"));
00222 slha.addLine(slhaMassLine(6, AlpgenHeader::mt, "top pole mass"));
00223 slha.addLine(slhaMassLine(23, AlpgenHeader::mz, "Z mass"));
00224 slha.addLine(slhaMassLine(24, AlpgenHeader::mw, "W mass"));
00225 slha.addLine(slhaMassLine(25, AlpgenHeader::mh, "H mass"));
00226
00227 char buffer[512];
00228
00229
00230 LHERunInfoProduct::Header lheAlpgenWgtHeader("AlpgenWgtFile");
00231 if (writeAlpgenWgtFile) {
00232 std::ifstream wgtascii((fileName_+".wgt").c_str());
00233 while(wgtascii.getline(buffer,512)) {
00234 lheAlpgenWgtHeader.addLine(std::string(buffer) + "\n");
00235 }
00236 }
00237
00238 LHERunInfoProduct::Header lheAlpgenParHeader("AlpgenParFile");
00239 if (writeAlpgenParFile) {
00240 std::ifstream parascii((fileName_+".par").c_str());
00241 while(parascii.getline(buffer,512)) {
00242 lheAlpgenParHeader.addLine(std::string(buffer) + "\n");
00243 }
00244 }
00245
00246
00247
00248
00249 LHERunInfoProduct::Header extraHeader(extraHeaderName_.c_str());
00250 if(writeExtraHeader) {
00251 std::ifstream extraascii(extraHeaderFileName_.c_str());
00252 while(extraascii.getline(buffer,512)) {
00253 extraHeader.addLine(std::string(buffer) + "\n");
00254 }
00255 }
00256
00257
00258 std::auto_ptr<LHERunInfoProduct> runInfo(new LHERunInfoProduct(heprup));
00259 runInfo->addHeader(comments);
00260 runInfo->addHeader(lheAlpgenUnwParHeader);
00261 if (writeAlpgenWgtFile)
00262 runInfo->addHeader(lheAlpgenWgtHeader);
00263 if (writeAlpgenParFile)
00264 runInfo->addHeader(lheAlpgenParHeader);
00265 runInfo->addHeader(slha);
00266 if(writeExtraHeader)
00267 runInfo->addHeader(extraHeader);
00268 run.put(runInfo);
00269
00270
00271 inputFile_.reset(new std::ifstream((fileName_ + ".unw").c_str()));
00272 if (!inputFile_->good())
00273 throw cms::Exception("Generator|AlpgenInterface")
00274 << "AlpgenSource was unable to open the file \""
00275 << fileName_ << ".unw\"." << std::endl;
00276
00277 }
00278
00279 template<typename T>
00280 T AlpgenSource::getParameter(AlpgenHeader::Parameter index) const
00281 {
00282 std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
00283 header.params.find(index);
00284 if (pos == header.params.end())
00285 throw cms::Exception("Generator|AlpgenInterface")
00286 << "Requested Alpgen parameter \""
00287 << AlpgenHeader::parameterName(index) << "\" "
00288 "not found in Alpgen parameter file." << std::endl;
00289
00290 return T(pos->second);
00291 }
00292
00293 template<typename T>
00294 T AlpgenSource::getParameter(AlpgenHeader::Parameter index,
00295 const T &defValue) const
00296 {
00297 std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
00298 header.params.find(index);
00299 if (pos == header.params.end())
00300 return defValue;
00301 else
00302 return T(pos->second);
00303 }
00304
00305 unsigned int AlpgenSource::processID() const
00306 {
00307
00308
00309 return header.ihrd * 100 +
00310 getParameter<unsigned int>(AlpgenHeader::ihvy, 0) * 10 +
00311 getParameter<unsigned int>(AlpgenHeader::njets, 0);
00312 }
00313
00314 bool AlpgenSource::readAlpgenEvent(lhef::HEPEUP &hepeup)
00315 {
00316 char buffer[512];
00317 double dummy;
00318 int nPart;
00319 double sWgtRes;
00320 double sQ;
00321
00322 inputFile_->getline(buffer, sizeof buffer);
00323 if (!inputFile_->good())
00324 return false;
00325
00326 std::istringstream ls(buffer);
00327
00328
00329 ls >> dummy >> dummy;
00330
00331
00332 ls >> nPart >> sWgtRes >> sQ;
00333
00334 if (ls.bad() || nPart < 1 || nPart > 1000)
00335 return false;
00336
00337
00338 hepeup.resize(nPart);
00339
00340
00341 hepeup.SCALUP = sQ;
00342 hepeup.XWGTUP = sWgtRes;
00343 hepeup.IDPRUP = processID();
00344
00345
00346 for(int i = 0; i != 2; ++i) {
00347 inputFile_->getline(buffer, sizeof buffer);
00348 std::istringstream ls(buffer);
00349 int flavour; ls >> flavour;
00350 int colour1; ls >> colour1;
00351 int colour2; ls >> colour2;
00352 double zMomentum; ls >> zMomentum;
00353
00354 if (inputFile_->bad())
00355 return false;
00356
00357
00358 hepeup.IDUP[i] = flavour;
00359 hepeup.ISTUP[i] = -1;
00360 hepeup.MOTHUP[i].first = 0;
00361 hepeup.MOTHUP[i].second = 0;
00362 hepeup.PUP[i][0] = 0.;
00363 hepeup.PUP[i][1] = 0.;
00364 hepeup.PUP[i][2] = zMomentum;
00365 hepeup.PUP[i][3] = std::fabs(zMomentum);
00366 hepeup.PUP[i][4] = 0.;
00367 if (colour1) colour1 += 500;
00368 if (colour2) colour2 += 500;
00369 hepeup.ICOLUP[i].first = colour1;
00370 hepeup.ICOLUP[i].second = colour2;
00371 }
00372
00373
00374 for(int i = 2; i != nPart; ++i) {
00375 inputFile_->getline(buffer, sizeof buffer);
00376 std::istringstream ls(buffer);
00377 int flavour; ls >> flavour;
00378 int colour1; ls >> colour1;
00379 int colour2; ls >> colour2;
00380 double px, py, pz, mass;
00381 ls >> px >> py >> pz >> mass;
00382 double energy = std::sqrt(px*px + py*py + pz*pz + mass*mass);
00383
00384 if (inputFile_->bad())
00385 return false;
00386
00387
00388 hepeup.IDUP[i] = flavour;
00389 hepeup.ISTUP[i] = 1;
00390 hepeup.MOTHUP[i].first = 1;
00391 hepeup.MOTHUP[i].second = 2;
00392 hepeup.PUP[i][0] = px;
00393 hepeup.PUP[i][1] = py;
00394 hepeup.PUP[i][2] = pz;
00395 hepeup.PUP[i][3] = energy;
00396 hepeup.PUP[i][4] = mass;
00397 if (colour1) colour1 += 500;
00398 if (colour2) colour2 += 500;
00399 hepeup.ICOLUP[i].first = colour1;
00400 hepeup.ICOLUP[i].second = colour2;
00401 }
00402
00403 return true;
00404 }
00405
00406 bool AlpgenSource::setRunAndEventInfo(edm::EventID&, edm::TimeValue_t&)
00407 {
00408
00409 hepeup_.reset(new lhef::HEPEUP);
00410
00411 lhef::HEPEUP& hepeup = *hepeup_;
00412
00413
00414 for(;;) {
00415 if (!readAlpgenEvent(hepeup)) {
00416 if (inputFile_->eof())
00417 return false;
00418
00419 throw cms::Exception("Generator|AlpgenInterface")
00420 << "AlpgenSource is not able read event no. "
00421 << nEvent_ << std::endl;
00422 }
00423
00424 nEvent_++;
00425 if (skipEvents_ > 0)
00426 skipEvents_--;
00427 else
00428 break;
00429 }
00430 return true;
00431 }
00432
00433 void AlpgenSource::produce(edm::Event &event)
00434 {
00435
00436
00437
00438
00439
00440 lhef::HEPEUP& hepeup = *hepeup_;
00441
00442 switch(header.ihrd) {
00443 case 1:
00444 case 2:
00445 case 3:
00446 case 4:
00447 case 10:
00448 case 14:
00449 case 15:
00450 alpgen::fixEventWZ(hepeup);
00451 break;
00452 case 5:
00453 alpgen::fixEventMultiBoson(hepeup);
00454 break;
00455 case 6:
00456 alpgen::fixEventTTbar(hepeup);
00457 break;
00458 case 8:
00459 alpgen::fixEventHiggsTTbar(hepeup);
00460 break;
00461 case 13:
00462 alpgen::fixEventSingleTop(hepeup, header.masses[AlpgenHeader::mb],
00463 int(header.params[AlpgenHeader::itopprc]));
00464 break;
00465 case 7:
00466 case 9:
00467 case 11:
00468 case 12:
00469 case 16:
00470
00471 break;
00472
00473 default:
00474 throw cms::Exception("Generator|AlpgenInterface")
00475 << "Unrecognized IHRD process code" << std::endl;
00476 }
00477
00478
00479 std::auto_ptr<LHEEventProduct> lheEvent(new LHEEventProduct(hepeup));
00480 event.put(lheEvent);
00481
00482 hepeup_.reset();
00483 }
00484
00485 DEFINE_FWK_INPUT_SOURCE(AlpgenSource);