CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/GeneratorInterface/AlpgenInterface/plugins/AlpgenSource.cc

Go to the documentation of this file.
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/ExternalInputSource.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::ExternalInputSource {
00030 public:
00032   AlpgenSource(const edm::ParameterSet &params,
00033                const edm::InputSourceDescription &desc);
00034 
00036   virtual ~AlpgenSource();
00037 
00038 private:
00039   virtual bool produce(edm::Event &event);
00040   virtual void beginRun(edm::Run &run);
00041 
00043   template<typename T>
00044   T getParameter(AlpgenHeader::Parameter index) const;
00046   template<typename T>
00047   T getParameter(AlpgenHeader::Parameter index, const T &defValue) const;
00048 
00051   std::string slhaMassLine(int pdgId, AlpgenHeader::Masses mass,
00052                            const std::string &comment) const;
00053 
00057   unsigned int processID() const;
00058 
00060   bool readAlpgenEvent(lhef::HEPEUP &hepeup);
00061 
00063   std::string                   fileName_;
00064 
00066   std::auto_ptr<std::ifstream>  inputFile_;
00067 
00069   unsigned long                 skipEvents_;
00070 
00072   unsigned long                 nEvent_;
00073 
00075   LHERunInfoProduct::Header     lheAlpgenUnwParHeader;
00077   AlpgenHeader                  header;
00078 
00080   std::string                   extraHeaderFileName_;
00081 
00083   std::string                   extraHeaderName_;
00084 
00086   bool                          writeAlpgenWgtFile;
00087   bool                          writeAlpgenParFile;
00088   bool                          writeExtraHeader;
00089 };
00090 
00091 AlpgenSource::AlpgenSource(const edm::ParameterSet &params,
00092                            const edm::InputSourceDescription &desc) :
00093   edm::ExternalInputSource(params, desc, false), 
00094   skipEvents_(params.getUntrackedParameter<unsigned int>("skipEvents", 0)),
00095   nEvent_(0), lheAlpgenUnwParHeader("AlpgenUnwParFile"),
00096   extraHeaderFileName_(params.getUntrackedParameter<std::string>("extraHeaderFileName","")),
00097   extraHeaderName_(params.getUntrackedParameter<std::string>("extraHeaderName","")),
00098   writeAlpgenWgtFile(params.getUntrackedParameter<bool>("writeAlpgenWgtFile", true)),
00099   writeAlpgenParFile(params.getUntrackedParameter<bool>("writeAlpgenParFile", true)),
00100   writeExtraHeader(params.getUntrackedParameter<bool>("writeExtraHeader", false))
00101 {
00102   std::vector<std::string> allFileNames = fileNames();
00103 
00104   // Only one filename
00105   if (allFileNames.size() != 1)
00106     throw cms::Exception("Generator|AlpgenInterface")
00107       << "AlpgenSource needs exactly one file specified "
00108       "for now." << std::endl;
00109 
00110   fileName_ = allFileNames[0];
00111 
00112   // Strip the "file:" prefix 
00113   if (fileName_.find("file:") != 0)
00114     throw cms::Exception("Generator|AlpgenInterface") << "AlpgenSource only supports the file: scheme "
00115       "for now." << std::endl;
00116   fileName_.erase(0, 5);
00117 
00118   // Open the _unw.par file to store additional 
00119   // informations in the LHERunInfoProduct
00120   std::ifstream reader((fileName_ + "_unw.par").c_str());
00121   if (!reader.good())
00122     throw cms::Exception("Generator|AlpgenInterface")
00123       << "AlpgenSource was unable to open the file \""
00124       << fileName_ << "_unw.par\"." << std::endl;
00125 
00126   // A full copy of the _unw.par file in the LHE header.
00127   char buffer[256];
00128   while(reader.getline(buffer, sizeof buffer))
00129     lheAlpgenUnwParHeader.addLine(std::string(buffer) + "\n");
00130 
00131   // Parse that header to setup an Alpgen header,
00132   // which will be used in the production itself.
00133    if (!header.parse(lheAlpgenUnwParHeader.begin(),
00134                     lheAlpgenUnwParHeader.end()))
00135     throw cms::Exception("Generator|AlpgenInterface")
00136       << "AlpgenSource was unable to parse the Alpgen "
00137       << "unweighted parameter file." << std::endl;
00138 
00139   // Declare the products.
00140   produces<LHERunInfoProduct, edm::InRun>();
00141   produces<LHEEventProduct>();
00142 }
00143 
00144 AlpgenSource::~AlpgenSource()
00145 {
00146 }
00147 
00148 std::string AlpgenSource::slhaMassLine(int pdgId, AlpgenHeader::Masses mass,
00149                                        const std::string &comment) const
00150 {
00151   std::ostringstream ss;
00152   ss << std::setw(9) << pdgId << "     " << std::scientific
00153      << std::setprecision(9) << header.masses[mass] << "   # "
00154      << comment << std::endl;
00155   return ss.str();
00156 }
00157 
00158 void AlpgenSource::beginRun(edm::Run &run)
00159 {
00160   // At this point, the lheUnwParHeader has the full contents of the _unw.par
00161   // file. So we can get the HEPRUP information from the LHE header itself.
00162   lhef::HEPRUP heprup;
00163 
00164   // Get basic run information.
00165   // Beam identity.
00166   heprup.IDBMUP.first = 2212;
00167   switch(getParameter<int>(AlpgenHeader::ih2)) {
00168   case 1:
00169     heprup.IDBMUP.second = 2212;
00170     break;
00171   case -1:
00172     heprup.IDBMUP.second = -2212;
00173     break;
00174   default:
00175     throw cms::Exception("Generator|AlpgenInterface")
00176       << "AlpgenSource was unable to understand the ih2 "
00177       << "parameter." << std::endl;
00178   }
00179 
00180   // Beam energy.
00181   heprup.EBMUP.second = heprup.EBMUP.first =
00182     getParameter<double>(AlpgenHeader::ebeam);
00183 
00184   // PDF info. Initially, Alpgen doesn't fill it.
00185   heprup.PDFGUP.first = -1;
00186   heprup.PDFGUP.second = -1;
00187   heprup.PDFSUP.first = -1;
00188   heprup.PDFSUP.second = -1;
00189 
00190   // Unweighted events.
00191   heprup.IDWTUP = 3;
00192 
00193   // Only one process.
00194   heprup.resize(1);
00195 
00196   // Cross section and error.
00197   heprup.XSECUP[0] = header.xsec;
00198   heprup.XERRUP[0] = header.xsecErr;
00199 
00200   // Maximum weight.
00201   heprup.XMAXUP[0] = header.xsec;
00202 
00203   // Process code for Pythia.
00204   heprup.LPRUP[0] = processID();
00205 
00206   // Comments on top.
00207   LHERunInfoProduct::Header comments;
00208   comments.addLine("\n");
00209   comments.addLine("\tExtracted by AlpgenInterface\n");
00210 
00211   // Add SLHA header containing particle masses from Alpgen.
00212   // Pythia6Hadronisation will feed the masses to Pythia automatically.
00213   LHERunInfoProduct::Header slha("slha");
00214   slha.addLine("\n# SLHA header containing masses from Alpgen\n");
00215   slha.addLine("Block MASS      #  Mass spectrum (kinematic masses)\n");
00216   slha.addLine("#       PDG   Mass\n");
00217   slha.addLine(slhaMassLine(4, AlpgenHeader::mc, "charm    pole mass"));
00218   slha.addLine(slhaMassLine(5, AlpgenHeader::mb, "bottom   pole mass"));
00219   slha.addLine(slhaMassLine(6, AlpgenHeader::mt, "top      pole mass"));
00220   slha.addLine(slhaMassLine(23, AlpgenHeader::mz, "Z        mass"));
00221   slha.addLine(slhaMassLine(24, AlpgenHeader::mw, "W        mass"));
00222   slha.addLine(slhaMassLine(25, AlpgenHeader::mh, "H        mass"));
00223 
00224   char buffer[512];
00225 
00226   // We also add the information on weighted events.
00227   LHERunInfoProduct::Header lheAlpgenWgtHeader("AlpgenWgtFile");
00228   if (writeAlpgenWgtFile) {
00229     std::ifstream wgtascii((fileName_+".wgt").c_str());
00230     while(wgtascii.getline(buffer,512)) {
00231       lheAlpgenWgtHeader.addLine(std::string(buffer) + "\n");
00232     }
00233   }
00234 
00235   LHERunInfoProduct::Header lheAlpgenParHeader("AlpgenParFile");
00236   if (writeAlpgenParFile) {
00237     std::ifstream parascii((fileName_+".par").c_str());
00238     while(parascii.getline(buffer,512)) {
00239       lheAlpgenParHeader.addLine(std::string(buffer) + "\n");
00240     }
00241   }
00242 
00243   // If requested by the user, we also add any specific header provided.
00244   // Nota bene: the header is put in the LHERunInfoProduct AS IT WAS GIVEN.
00245   // That means NO CROSS-CHECKS WHATSOEVER. Use with care.
00246   LHERunInfoProduct::Header extraHeader(extraHeaderName_.c_str());
00247   if(writeExtraHeader) {
00248     std::ifstream extraascii(extraHeaderFileName_.c_str());
00249     while(extraascii.getline(buffer,512)) {
00250       extraHeader.addLine(std::string(buffer) + "\n");
00251     }
00252   }
00253 
00254   // Build the final Run info object. Backwards-compatible order.
00255   std::auto_ptr<LHERunInfoProduct> runInfo(new LHERunInfoProduct(heprup));
00256   runInfo->addHeader(comments);
00257   runInfo->addHeader(lheAlpgenUnwParHeader);
00258   if (writeAlpgenWgtFile)
00259     runInfo->addHeader(lheAlpgenWgtHeader);
00260   if (writeAlpgenParFile)
00261     runInfo->addHeader(lheAlpgenParHeader);
00262   runInfo->addHeader(slha);
00263   if(writeExtraHeader)
00264     runInfo->addHeader(extraHeader);
00265   run.put(runInfo);
00266 
00267   // Open the .unw file in the heap, and set the global pointer to it.
00268   inputFile_.reset(new std::ifstream((fileName_ + ".unw").c_str()));
00269   if (!inputFile_->good())
00270     throw cms::Exception("Generator|AlpgenInterface")
00271       << "AlpgenSource was unable to open the file \""
00272       << fileName_ << ".unw\"." << std::endl;
00273 
00274 }
00275 
00276 template<typename T>
00277 T AlpgenSource::getParameter(AlpgenHeader::Parameter index) const
00278 {
00279   std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
00280     header.params.find(index);
00281   if (pos == header.params.end())
00282     throw cms::Exception("Generator|AlpgenInterface")
00283       << "Requested Alpgen parameter \""
00284       << AlpgenHeader::parameterName(index) << "\" "
00285       "not found in Alpgen parameter file." << std::endl;
00286 
00287   return T(pos->second);
00288 }
00289 
00290 template<typename T>
00291 T AlpgenSource::getParameter(AlpgenHeader::Parameter index,
00292                              const T &defValue) const
00293 {
00294   std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
00295     header.params.find(index);
00296   if (pos == header.params.end())
00297     return defValue;
00298   else
00299     return T(pos->second);
00300 }
00301 
00302 unsigned int AlpgenSource::processID() const
00303 {
00304   // return 661; // The original, old thing.
00305   // digits #XYZ: X = ihrd, Y = ihvy, Z = njets
00306   return header.ihrd * 100 +
00307     getParameter<unsigned int>(AlpgenHeader::ihvy, 0) * 10 +
00308     getParameter<unsigned int>(AlpgenHeader::njets, 0);
00309 }
00310 
00311 bool AlpgenSource::readAlpgenEvent(lhef::HEPEUP &hepeup)
00312 {
00313   char buffer[512];
00314   double dummy;
00315   int nPart;
00316   double sWgtRes;
00317   double sQ;
00318 
00319   inputFile_->getline(buffer, sizeof buffer);
00320   if (!inputFile_->good())
00321     return false;
00322 
00323   std::istringstream ls(buffer);
00324 
00325   // Event number and process don't matter (or do they?)
00326   ls >> dummy >> dummy;
00327 
00328   // Number of particles in the record, sample's average weight and Q scale
00329   ls >> nPart >> sWgtRes >> sQ;
00330 
00331   if (ls.bad() || nPart < 1 || nPart > 1000)
00332     return false;
00333 
00334   // Make room for the particles listed in the Alpgen file
00335   hepeup.resize(nPart);
00336 
00337   // Scales, weights and process ID.
00338   hepeup.SCALUP = sQ;
00339   hepeup.XWGTUP = sWgtRes;
00340   hepeup.IDPRUP = processID();
00341 
00342   // Incoming lines
00343   for(int i = 0; i != 2; ++i) {
00344     inputFile_->getline(buffer, sizeof buffer);
00345     std::istringstream ls(buffer);
00346     int flavour; ls >> flavour;
00347     int colour1; ls >> colour1;
00348     int colour2; ls >> colour2;
00349     double zMomentum; ls >> zMomentum;
00350 
00351     if (inputFile_->bad())
00352       return false;
00353 
00354     // Setting the HEPEUP of the incoming lines.
00355     hepeup.IDUP[i] = flavour;
00356     hepeup.ISTUP[i] = -1;
00357     hepeup.MOTHUP[i].first = 0;
00358     hepeup.MOTHUP[i].second = 0;
00359     hepeup.PUP[i][0] = 0.;
00360     hepeup.PUP[i][1] = 0.;
00361     hepeup.PUP[i][2] = zMomentum;
00362     hepeup.PUP[i][3] = std::fabs(zMomentum);
00363     hepeup.PUP[i][4] = 0.;
00364     if (colour1) colour1 += 500;
00365     if (colour2) colour2 += 500;
00366     hepeup.ICOLUP[i].first = colour1;
00367     hepeup.ICOLUP[i].second = colour2;
00368   }
00369 
00370   // Outgoing lines
00371   for(int i = 2; i != nPart; ++i) {
00372     inputFile_->getline(buffer, sizeof buffer);
00373     std::istringstream ls(buffer);
00374     int flavour; ls >> flavour;
00375     int colour1; ls >> colour1;
00376     int colour2; ls >> colour2;
00377     double px, py, pz, mass; 
00378     ls >> px >> py >> pz >> mass;
00379     double energy = std::sqrt(px*px + py*py + pz*pz + mass*mass);
00380 
00381     if (inputFile_->bad())
00382       return false;
00383 
00384     // Setting the HEPEUP of the outgoing lines.
00385     hepeup.IDUP[i] = flavour;
00386     hepeup.ISTUP[i] = 1;
00387     hepeup.MOTHUP[i].first = 1;
00388     hepeup.MOTHUP[i].second = 2;
00389     hepeup.PUP[i][0] = px;
00390     hepeup.PUP[i][1] = py;
00391     hepeup.PUP[i][2] = pz;
00392     hepeup.PUP[i][3] = energy;
00393     hepeup.PUP[i][4] = mass;
00394     if (colour1) colour1 += 500;
00395     if (colour2) colour2 += 500;
00396     hepeup.ICOLUP[i].first = colour1;
00397     hepeup.ICOLUP[i].second = colour2;
00398   }
00399 
00400   return true;
00401 }
00402 
00403 bool AlpgenSource::produce(edm::Event &event)
00404 {
00405   // The LHE Event Record
00406   lhef::HEPEUP hepeup;
00407 
00408   // Read the .unw file until it is over.
00409   for(;;) {
00410     if (!readAlpgenEvent(hepeup)) {
00411       if (inputFile_->eof())
00412         return false;
00413 
00414       throw cms::Exception("Generator|AlpgenInterface")
00415         << "AlpgenSource is not able read event no. "
00416         << nEvent_ << std::endl;
00417     }
00418 
00419     nEvent_++;
00420     if (skipEvents_ > 0)
00421       skipEvents_--;
00422     else
00423       break;
00424   }
00425 
00426   // Here are the Alpgen routines for filling up the rest
00427   // of the LHE Event Record. The .unw file has the information
00428   // in a compressed way, e.g. it doesn't list the W boson - 
00429   // one has to reconstruct it from the e nu pair.
00430 
00431   switch(header.ihrd) {
00432   case 1:
00433   case 2:
00434   case 3:
00435   case 4:
00436   case 10:
00437   case 14:
00438   case 15:
00439     alpgen::fixEventWZ(hepeup);
00440     break;
00441   case 5:
00442     alpgen::fixEventMultiBoson(hepeup);
00443     break;
00444   case 6:
00445     alpgen::fixEventTTbar(hepeup);
00446     break;
00447   case 8:
00448     alpgen::fixEventHiggsTTbar(hepeup);
00449     break;
00450   case 13:
00451     alpgen::fixEventSingleTop(hepeup, header.masses[AlpgenHeader::mb],
00452                               int(header.params[AlpgenHeader::itopprc]));
00453     break;    
00454   case 7:
00455   case 9:
00456   case 11:
00457   case 12:
00458   case 16:
00459     // No fixes needed.
00460     break;
00461 
00462   default: 
00463     throw cms::Exception("Generator|AlpgenInterface") 
00464       << "Unrecognized IHRD process code" << std::endl;
00465   }
00466 
00467   // Create the LHEEventProduct and put it into the Event.
00468   std::auto_ptr<LHEEventProduct> lheEvent(new LHEEventProduct(hepeup));
00469   event.put(lheEvent);
00470 
00471   return true;
00472 }
00473 
00474 DEFINE_FWK_INPUT_SOURCE(AlpgenSource);