CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/GeneratorInterface/Herwig6Interface/plugins/Herwig6Hadronizer.cc

Go to the documentation of this file.
00001 #include <cstring>
00002 #include <sstream>
00003 #include <string>
00004 #include <vector>
00005 #include <memory>
00006 #include <map>
00007 #include <set>
00008 
00009 #include <boost/shared_ptr.hpp>
00010 #include <boost/algorithm/string/classification.hpp>
00011 #include <boost/algorithm/string/split.hpp>
00012 
00013 #include <HepMC/GenEvent.h>
00014 #include <HepMC/GenParticle.h>
00015 #include <HepMC/GenVertex.h>
00016 #include <HepMC/PdfInfo.h>
00017 #include <HepMC/HerwigWrapper.h>
00018 #include <HepMC/HEPEVT_Wrapper.h>
00019 #include <HepMC/IO_HERWIG.h>
00020 #include "HepPID/ParticleIDTranslations.hh"
00021 
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 
00024 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00025 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00026 
00027 #include "SimDataFormats/GeneratorProducts/interface/LesHouches.h"
00028 #include "SimDataFormats/GeneratorProducts/interface/LHECommonBlocks.h"
00029 
00030 #include "GeneratorInterface/Core/interface/ParameterCollector.h"
00031 #include "GeneratorInterface/Core/interface/BaseHadronizer.h"
00032 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
00033 #include "GeneratorInterface/Core/interface/HadronizerFilter.h"
00034 
00035 #include "GeneratorInterface/LHEInterface/interface/LHEEvent.h"
00036 
00037 #include "GeneratorInterface/Herwig6Interface/interface/Herwig6Instance.h"
00038 #include "GeneratorInterface/Herwig6Interface/interface/herwig.h"
00039 
00040 #include "DataFormats/Math/interface/LorentzVector.h"
00041 
00042 extern "C" {
00043   void hwuidt_(int *iopt, int *ipdg, int *iwig, char nwig[8]);
00044   double hwualf_(int *mode, double* scale);
00045   double hwuaem_(double* scale);
00046 }
00047 
00048 // helpers
00049 namespace {
00050         static inline bool call_hwmsct()
00051         {
00052                 int result;
00053                 hwmsct(&result);
00054                 return result;
00055         }
00056 
00057         static int pdgToHerwig(int ipdg, char *nwig)
00058         {
00059                 int iopt = 1;
00060                 int iwig = 0;
00061                 hwuidt_(&iopt, &ipdg, &iwig, nwig);
00062                 return ipdg ? iwig : 0;
00063         }
00064 
00065         static bool markStable(int pdgId)
00066         {
00067                 char nwig[9] = "        ";
00068                 if (!pdgToHerwig(pdgId, nwig))
00069                         return false;
00070                 hwusta(nwig, 1);
00071                 return true;
00072         }
00073 }
00074 
00075 class Herwig6Hadronizer : public gen::BaseHadronizer,
00076                           public gen::Herwig6Instance {
00077     public:
00078         Herwig6Hadronizer(const edm::ParameterSet &params);
00079         ~Herwig6Hadronizer();
00080 
00081         void setSLHAFromHeader(const std::vector<std::string> &lines);
00082         bool initialize(const lhef::HEPRUP *heprup);
00083 
00084         bool readSettings( int );
00085         
00086         // bool initializeForInternalPartons() { return initialize(0); }
00087         // bool initializeForExternalPartons() { return initialize(lheRunInfo()->getHEPRUP()); }
00088 
00089         bool initializeForInternalPartons();
00090         bool initializeForExternalPartons() { return initializeForInternalPartons(); }
00091 
00092         bool declareStableParticles( const std::vector<int> &pdgIds);
00093         bool declareSpecialSettings( const std::vector<std::string> );
00094 
00095         void statistics();
00096 
00097 
00098         bool generatePartonsAndHadronize() { return hadronize(); }
00099         bool hadronize();
00100         bool decay();
00101         bool residualDecay();
00102         void finalizeEvent();
00103 
00104         const char *classname() const { return "Herwig6Hadronizer"; }
00105 
00106     private:
00107         void clear();
00108 
00109         int pythiaStatusCode(const HepMC::GenParticle *p) const;
00110         void pythiaStatusCodes();
00111 
00112         virtual void upInit();
00113         virtual void upEvnt();
00114 
00115         HepMC::IO_HERWIG                conv;
00116         bool                            needClear;
00117         bool                            externalPartons;
00118 
00119         gen::ParameterCollector         parameters;
00120         int                             herwigVerbosity;
00121         int                             hepmcVerbosity;
00122         int                             maxEventsToPrint;
00123         bool                            printCards;
00124         bool                            emulatePythiaStatusCodes;
00125         double                          comEnergy;
00126         bool                            useJimmy;
00127         bool                            doMPInteraction;
00128         int                             numTrials;
00129         bool                            fConvertToPDG;  // convert PIDs 
00130 
00131         bool                            readMCatNLOfile;
00132 
00133   // -------------------------------------------------------------------------------
00134         std::string                     particleSpecFileName; //Lars 20/Jul/2011
00135         bool                            readParticleSpecFile;
00136   // -------------------------------------------------------------------------------
00137 
00138 };
00139 
00140 extern "C" {
00141         void hwwarn_(const char *method, int *id);
00142         void setherwpdf_(void);
00143         void mysetpdfpath_(const char *path);
00144 } // extern "C"
00145 
00146 Herwig6Hadronizer::Herwig6Hadronizer(const edm::ParameterSet &params) :
00147         BaseHadronizer(params),
00148         needClear(false),
00149         parameters(params.getParameter<edm::ParameterSet>("HerwigParameters")),
00150         herwigVerbosity(params.getUntrackedParameter<int>("herwigVerbosity", 0)),
00151         hepmcVerbosity(params.getUntrackedParameter<int>("hepmcVerbosity", 0)),
00152         maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
00153         printCards(params.getUntrackedParameter<bool>("printCards", false)),
00154         emulatePythiaStatusCodes(params.getUntrackedParameter<bool>("emulatePythiaStatusCodes", false)),
00155         comEnergy(params.getParameter<double>("comEnergy")),
00156         useJimmy(params.getParameter<bool>("useJimmy")),
00157         doMPInteraction(params.getParameter<bool>("doMPInteraction")),
00158         numTrials(params.getUntrackedParameter<int>("numTrialsMPI", 100)),
00159         readMCatNLOfile(false),
00160 
00161         // added to be able to read external particle spectrum file
00162         particleSpecFileName(params.getUntrackedParameter<std::string>("ParticleSpectrumFileName","")),
00163         readParticleSpecFile(params.getUntrackedParameter<bool>("readParticleSpecFile", false))
00164   
00165 {
00166   
00167   fConvertToPDG = false;
00168   if ( params.exists( "doPDGConvert" ) )
00169     fConvertToPDG = params.getParameter<bool>("doPDGConvert");
00170 }
00171 
00172 Herwig6Hadronizer::~Herwig6Hadronizer()
00173 {
00174         clear();
00175 }
00176 
00177 void Herwig6Hadronizer::clear()
00178 {
00179         if (!needClear)
00180                 return;
00181 
00182         // teminate elementary process
00183         call(hwefin);
00184         if (useJimmy)
00185           call(jmefin);
00186 
00187         needClear = false;
00188 }
00189 
00190 void Herwig6Hadronizer::setSLHAFromHeader(
00191                                 const std::vector<std::string> &lines)
00192 {
00193         std::set<std::string> blocks;
00194         std::string block;
00195         for(std::vector<std::string>::const_iterator iter = lines.begin();
00196             iter != lines.end(); ++iter) {
00197                 std::string line = *iter;
00198                 std::transform(line.begin(), line.end(),
00199                                line.begin(), (int(*)(int))std::toupper);
00200                 std::string::size_type pos = line.find('#');
00201                 if (pos != std::string::npos)
00202                         line.resize(pos);
00203 
00204                 if (line.empty())
00205                         continue;
00206 
00207                 if (!boost::algorithm::is_space()(line[0])) {
00208                         std::vector<std::string> tokens;
00209                         boost::split(tokens, line,
00210                                      boost::algorithm::is_space(),
00211                                      boost::token_compress_on);
00212                         if (!tokens.size())
00213                                 continue;
00214                         block.clear();
00215                         if (tokens.size() < 2)
00216                                 continue;
00217                         if (tokens[0] == "BLOCK")
00218                                 block = tokens[1];
00219                         else if (tokens[0] == "DECAY")
00220                                 block = "DECAY";
00221 
00222                         if (block.empty())
00223                                 continue;
00224 
00225                         if (!blocks.count(block)) {
00226                                 blocks.insert(block);
00227                                 edm::LogWarning("Generator|Herwig6Hadronzier")
00228                                         << "Unsupported SLHA block \"" << block
00229                                         << "\".  It will be ignored."
00230                                         << std::endl;
00231                         }
00232                 }
00233         }
00234 }
00235 
00236 bool Herwig6Hadronizer:: readSettings( int key )
00237 {
00238 
00239    clear();
00240    const lhef::HEPRUP* heprup = lheRunInfo()->getHEPRUP();
00241    externalPartons = ( heprup != 0 );
00242    
00243    if ( key == 0 && externalPartons ) return false;
00244    if ( key > 0 && !externalPartons ) return false;
00245 
00246    std::ostringstream info;
00247    info << "---------------------------------------------------\n"; 
00248    info << "Taking in settinsg for Herwig6Hadronizer for "
00249         << (externalPartons ? "external" : "internal") << " partons\n";
00250    info << "---------------------------------------------------\n";
00251 
00252    info << "   Herwig verbosity level         = " << herwigVerbosity << "\n";
00253    info << "   HepMC verbosity                = " << hepmcVerbosity << "\n";
00254    info << "   Number of events to be printed = " << maxEventsToPrint << "\n";
00255 
00256    // Call hwudat to set up HERWIG block data
00257    hwudat();
00258 
00259    // Setting basic parameters
00260    if (externalPartons) {
00261                 hwproc.PBEAM1 = heprup->EBMUP.first;
00262                 hwproc.PBEAM2 = heprup->EBMUP.second;
00263                 pdgToHerwig(heprup->IDBMUP.first, hwbmch.PART1);
00264                 pdgToHerwig(heprup->IDBMUP.second, hwbmch.PART2);
00265    } 
00266    else 
00267    {
00268                 hwproc.PBEAM1 = 0.5 * comEnergy;
00269                 hwproc.PBEAM2 = 0.5 * comEnergy;
00270                 pdgToHerwig(2212, hwbmch.PART1);
00271                 pdgToHerwig(2212, hwbmch.PART2);
00272    }
00273 
00274    if (useJimmy) 
00275    {
00276                 info << "   HERWIG will be using JIMMY for UE/MI.\n";
00277                 jmparm.MSFLAG = 1;
00278                 if (doMPInteraction)
00279                         info << "   JIMMY trying to generate multiple interactions.\n";
00280    }
00281 
00282 
00283    // set the IPROC already here... needed for VB pairs
00284 
00285    bool iprocFound=false;
00286         
00287    for(gen::ParameterCollector::const_iterator line = parameters.begin();
00288                                                line != parameters.end(); ++line) 
00289    {
00290           if(!strcmp((line->substr(0,5)).c_str(),"IPROC")) {
00291             if (!give(*line))
00292               throw edm::Exception(edm::errors::Configuration)
00293                 << "Herwig 6 did not accept the following: \""
00294                 << *line << "\"." << std::endl;
00295             else iprocFound=true;
00296           }
00297    }
00298         
00299    if (!iprocFound && !externalPartons)
00300           throw edm::Exception(edm::errors::Configuration)
00301             << "You have to define the process with IPROC."  << std::endl;
00302 
00303    // initialize other common blocks ...
00304    call(hwigin); // default init
00305    
00306    hwevnt.MAXER = 100000000;    // O(inf)
00307    hwpram.LWSUD = 0;            // don't write Sudakov form factors
00308    hwdspn.LWDEC = 0;            // don't write three/four body decays
00309                                 // (no fort.77 and fort.88 ...)
00310         // init LHAPDF glue
00311         std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
00312         for(unsigned int i = 0; i < 2; i++) {
00313                 hwpram.MODPDF[i] = -111;
00314                 std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
00315         }
00316 
00317         hwevnt.MAXPR = maxEventsToPrint;
00318         hwpram.IPRINT = herwigVerbosity;
00319 //      hwprop.RMASS[6] = 175.0;        //FIXME
00320 
00321         if (printCards) {
00322                 info << "\n";
00323                 info << "------------------------------------\n";
00324                 info << "Reading HERWIG parameters\n";
00325                 info << "------------------------------------\n";
00326     
00327         }
00328         for(gen::ParameterCollector::const_iterator line = parameters.begin();
00329             line != parameters.end(); ++line) {
00330                 if (printCards)
00331                         info << "   " << *line << "\n";
00332                 if (!give(*line))
00333                         throw edm::Exception(edm::errors::Configuration)
00334                                 << "Herwig 6 did not accept the following: \""
00335                                 << *line << "\"." << std::endl;
00336         }
00337 
00338         if (printCards)
00339                 info << "\n";
00340 
00341         if (externalPartons) {
00342                 std::vector<std::string> slha =
00343                                 lheRunInfo()->findHeader("slha");
00344                 if (!slha.empty())
00345                         setSLHAFromHeader(slha);
00346         }
00347 
00348         needClear = true;
00349 
00350         std::pair<int, int> pdfs(-1, -1);
00351         if (externalPartons)
00352                 pdfs = lheRunInfo()->pdfSetTranslation();
00353 
00354         if (hwpram.MODPDF[0] != -111 || hwpram.MODPDF[1] != -111) {
00355                 for(unsigned int i = 0; i < 2; i++)
00356                         if (hwpram.MODPDF[i] == -111)
00357                                 hwpram.MODPDF[i] = -1;
00358 
00359                 if (pdfs.first != -1 || pdfs.second != -1)
00360                         edm::LogError("Generator|Herwig6Hadronzier")
00361                                 << "Both external Les Houches event and "
00362                                    "config file specify a PDF set.  "
00363                                    "User PDF will override external one."
00364                                 << std::endl;
00365 
00366                 pdfs.first = hwpram.MODPDF[0] != -111 ? hwpram.MODPDF[0] : -1;
00367                 pdfs.second = hwpram.MODPDF[1] != -111 ? hwpram.MODPDF[1] : -1;
00368         }
00369 
00370         hwpram.MODPDF[0] = pdfs.first;
00371         hwpram.MODPDF[1] = pdfs.second;
00372 
00373         if (externalPartons)
00374                 hwproc.IPROC = -1;
00375 
00376 
00377         //Lars: lower EFFMIN threshold, to continue execution of IPROC=4000, lambda'_211=0.01 at LM7,10
00378         if( readParticleSpecFile ) {
00379           openParticleSpecFile(particleSpecFileName.c_str());
00380           hwpram.EFFMIN = 1e-5;
00381         }
00382         
00383 
00384    edm::LogInfo(info.str());
00385 
00386    return true;
00387 
00388 }
00389 
00390 bool Herwig6Hadronizer::initializeForInternalPartons()
00391 {
00392 
00393    if (useJimmy) call(jimmin);
00394    
00395    call(hwuinc);
00396 
00397    // initialize HERWIG event generation
00398    call(hweini);
00399 
00400    if (useJimmy) call(jminit);
00401 
00402    return true;
00403 
00404 }
00405 
00406 
00407 bool Herwig6Hadronizer::initialize(const lhef::HEPRUP *heprup)
00408 {
00409         clear();
00410 
00411         externalPartons = (heprup != 0);
00412 
00413         std::ostringstream info;
00414         info << "---------------------------------------------------\n"; 
00415         info << "Initializing Herwig6Hadronizer for "
00416              << (externalPartons ? "external" : "internal") << " partons\n";
00417         info << "---------------------------------------------------\n";
00418 
00419         info << "   Herwig verbosity level         = " << herwigVerbosity << "\n";
00420         info << "   HepMC verbosity                = " << hepmcVerbosity << "\n";
00421         info << "   Number of events to be printed = " << maxEventsToPrint << "\n";
00422 
00423         // Call hwudat to set up HERWIG block data
00424         hwudat();
00425 
00426         // Setting basic parameters
00427         if (externalPartons) {
00428                 hwproc.PBEAM1 = heprup->EBMUP.first;
00429                 hwproc.PBEAM2 = heprup->EBMUP.second;
00430                 pdgToHerwig(heprup->IDBMUP.first, hwbmch.PART1);
00431                 pdgToHerwig(heprup->IDBMUP.second, hwbmch.PART2);
00432         } else {
00433                 hwproc.PBEAM1 = 0.5 * comEnergy;
00434                 hwproc.PBEAM2 = 0.5 * comEnergy;
00435                 pdgToHerwig(2212, hwbmch.PART1);
00436                 pdgToHerwig(2212, hwbmch.PART2);
00437         }
00438 
00439         if (useJimmy) {
00440                 info << "   HERWIG will be using JIMMY for UE/MI.\n";
00441                 jmparm.MSFLAG = 1;
00442                 if (doMPInteraction)
00443                         info << "   JIMMY trying to generate multiple interactions.\n";
00444         }
00445 
00446         // set the IPROC already here... needed for VB pairs
00447 
00448         bool iprocFound=false;
00449         
00450         for(gen::ParameterCollector::const_iterator line = parameters.begin();
00451             line != parameters.end(); ++line) {
00452           if(!strcmp((line->substr(0,5)).c_str(),"IPROC")) {
00453             if (!give(*line))
00454               throw edm::Exception(edm::errors::Configuration)
00455                 << "Herwig 6 did not accept the following: \""
00456                 << *line << "\"." << std::endl;
00457             else iprocFound=true;
00458           }
00459         }
00460         
00461         if (!iprocFound && !externalPartons)
00462           throw edm::Exception(edm::errors::Configuration)
00463             << "You have to define the process with IPROC."  << std::endl;
00464 
00465         // initialize other common blocks ...
00466         call(hwigin);
00467         hwevnt.MAXER = 100000000;       // O(inf)
00468         hwpram.LWSUD = 0;               // don't write Sudakov form factors
00469         hwdspn.LWDEC = 0;               // don't write three/four body decays
00470                                         // (no fort.77 and fort.88 ...)
00471 
00472         // init LHAPDF glue
00473 
00474         std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
00475         for(unsigned int i = 0; i < 2; i++) {
00476                 hwpram.MODPDF[i] = -111;
00477                 std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
00478         }
00479 
00480         if (useJimmy)
00481                 call(jimmin);
00482 
00483         hwevnt.MAXPR = maxEventsToPrint;
00484         hwpram.IPRINT = herwigVerbosity;
00485 //      hwprop.RMASS[6] = 175.0;        //FIXME
00486 
00487         if (printCards) {
00488                 info << "\n";
00489                 info << "------------------------------------\n";
00490                 info << "Reading HERWIG parameters\n";
00491                 info << "------------------------------------\n";
00492     
00493         }
00494         for(gen::ParameterCollector::const_iterator line = parameters.begin();
00495             line != parameters.end(); ++line) {
00496                 if (printCards)
00497                         info << "   " << *line << "\n";
00498                 if (!give(*line))
00499                         throw edm::Exception(edm::errors::Configuration)
00500                                 << "Herwig 6 did not accept the following: \""
00501                                 << *line << "\"." << std::endl;
00502         }
00503 
00504         if (printCards)
00505                 info << "\n";
00506 
00507         if (externalPartons) {
00508                 std::vector<std::string> slha =
00509                                 lheRunInfo()->findHeader("slha");
00510                 if (!slha.empty())
00511                         setSLHAFromHeader(slha);
00512         }
00513 
00514         needClear = true;
00515 
00516         std::pair<int, int> pdfs(-1, -1);
00517         if (externalPartons)
00518                 pdfs = lheRunInfo()->pdfSetTranslation();
00519 
00520         if (hwpram.MODPDF[0] != -111 || hwpram.MODPDF[1] != -111) {
00521                 for(unsigned int i = 0; i < 2; i++)
00522                         if (hwpram.MODPDF[i] == -111)
00523                                 hwpram.MODPDF[i] = -1;
00524 
00525                 if (pdfs.first != -1 || pdfs.second != -1)
00526                         edm::LogError("Generator|Herwig6Hadronzier")
00527                                 << "Both external Les Houches event and "
00528                                    "config file specify a PDF set.  "
00529                                    "User PDF will override external one."
00530                                 << std::endl;
00531 
00532                 pdfs.first = hwpram.MODPDF[0] != -111 ? hwpram.MODPDF[0] : -1;
00533                 pdfs.second = hwpram.MODPDF[1] != -111 ? hwpram.MODPDF[1] : -1;
00534         }
00535 
00536         hwpram.MODPDF[0] = pdfs.first;
00537         hwpram.MODPDF[1] = pdfs.second;
00538 
00539         if (externalPartons)
00540                 hwproc.IPROC = -1;
00541 
00542         //Lars: lower EFFMIN threshold, to continue execution of IPROC=4000, lambda'_211=0.01 at LM7,10
00543         if( readParticleSpecFile ) {
00544           openParticleSpecFile(particleSpecFileName.c_str());
00545           hwpram.EFFMIN = 1e-5;
00546         }
00547         
00548         // HERWIG preparations ...
00549         call(hwuinc);
00550         markStable(13);          // MU+
00551         markStable(-13);         // MU-
00552         markStable(3112);        // SIGMA+
00553         markStable(-3112);       // SIGMABAR+
00554         markStable(3222);        // SIGMA-
00555         markStable(-3222);       // SIGMABAR-
00556         markStable(3122);        // LAMBDA0
00557         markStable(-3122);       // LAMBDABAR0
00558         markStable(3312);        // XI-
00559         markStable(-3312);       // XIBAR+
00560         markStable(3322);        // XI0
00561         markStable(-3322);       // XI0BAR
00562         markStable(3334);        // OMEGA-
00563         markStable(-3334);       // OMEGABAR+
00564         markStable(211);         // PI+
00565         markStable(-211);        // PI-
00566         markStable(321);         // K+
00567         markStable(-321);        // K-
00568         markStable(310);         // K_S0
00569         markStable(130);         // K_L0
00570 
00571         // better: merge with declareStableParticles
00572         // and get the list from configuration / Geant4 / Core somewhere
00573 
00574         // initialize HERWIG event generation
00575         call(hweini);
00576 
00577         if (useJimmy)
00578                 call(jminit);
00579 
00580         edm::LogInfo(info.str());
00581 
00582         return true;
00583 }
00584 
00585 bool Herwig6Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00586 {
00587         markStable(13);          // MU+
00588         markStable(-13);         // MU-
00589         markStable(3112);        // SIGMA+
00590         markStable(-3112);       // SIGMABAR+
00591         markStable(3222);        // SIGMA-
00592         markStable(-3222);       // SIGMABAR-
00593         markStable(3122);        // LAMBDA0
00594         markStable(-3122);       // LAMBDABAR0
00595         markStable(3312);        // XI-
00596         markStable(-3312);       // XIBAR+
00597         markStable(3322);        // XI0
00598         markStable(-3322);       // XI0BAR
00599         markStable(3334);        // OMEGA-
00600         markStable(-3334);       // OMEGABAR+
00601         markStable(211);         // PI+
00602         markStable(-211);        // PI-
00603         markStable(321);         // K+
00604         markStable(-321);        // K-
00605         markStable(310);         // K_S0
00606         markStable(130);         // K_L0
00607 
00608         for(std::vector<int>::const_iterator iter = pdgIds.begin();
00609             iter != pdgIds.end(); ++iter)
00610                 if (!markStable(*iter))
00611                         return false;
00612         return true;
00613 }
00614 
00615 bool Herwig6Hadronizer::declareSpecialSettings( const std::vector<std::string> )
00616 {
00617    return true;
00618 }
00619 
00620 void Herwig6Hadronizer::statistics()
00621 {
00622         if (!runInfo().internalXSec()) {
00623           // not set via LHE, so get it from HERWIG
00624           // the reason is that HERWIG doesn't compute the xsec
00625           // in all LHE modes
00626 
00627           double RNWGT = 1. / hwevnt.NWGTS;
00628           double AVWGT = hwevnt.WGTSUM * RNWGT;
00629 
00630           double xsec = 1.0e3 * AVWGT;
00631 
00632           runInfo().setInternalXSec(xsec);
00633         }
00634 }
00635 
00636 bool Herwig6Hadronizer::hadronize()
00637 {
00638         // hard process generation, parton shower, hadron formation
00639 
00640         InstanceWrapper wrapper(this);  // safe guard
00641 
00642         event().reset();
00643 
00644         // call herwig routines to create HEPEVT
00645 
00646         hwuine();       // initialize event
00647         
00648         if (callWithTimeout(10, hwepro)) { // process event and PS
00649           // We hung for more than 10 seconds
00650           int error = 199;
00651           hwwarn_("HWHGUP", &error);
00652         }
00653         
00654         hwbgen();       // parton cascades
00655 
00656         // call jimmy ... only if event is not killed yet by HERWIG
00657         if (useJimmy && doMPInteraction && !hwevnt.IERROR && call_hwmsct()) {
00658           if (lheEvent()) lheEvent()->count(lhef::LHERunInfo::kKilled);
00659           return false;
00660         }
00661         
00662         hwdhob();       // heavy quark decays
00663         hwcfor();       // cluster formation
00664         hwcdec();       // cluster decays
00665         
00666         // if event *not* killed by HERWIG, return true
00667         if (hwevnt.IERROR) {
00668           hwufne();     // finalize event, to keep system clean
00669           if (lheEvent()) lheEvent()->count(lhef::LHERunInfo::kKilled);
00670           return false;
00671         }
00672 
00673         if (lheEvent()) lheEvent()->count(lhef::LHERunInfo::kAccepted);
00674         return true;
00675 }
00676 
00677 void Herwig6Hadronizer::finalizeEvent()
00678 {
00679   lhef::LHEEvent::fixHepMCEventTimeOrdering(event().get());
00680   
00681   HepMC::PdfInfo pdfInfo;
00682   if(externalPartons) {
00683     lheEvent()->fillEventInfo( event().get() );
00684     lheEvent()->fillPdfInfo( &pdfInfo );
00685 
00686     // for MC@NLO: IDWRUP is not filled...
00687     if(event()->signal_process_id()==0)
00688       event()->set_signal_process_id( abs(hwproc.IPROC) );
00689 
00690   }
00691     
00692   HepMC::GenParticle* incomingParton = NULL;
00693   HepMC::GenParticle* targetParton = NULL;
00694   
00695   HepMC::GenParticle* incomingProton = NULL;
00696   HepMC::GenParticle* targetProton = NULL;
00697   
00698   // find incoming parton (first entry with IST=121)
00699   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00700       (it != event()->particles_end() && incomingParton==NULL); it++)
00701     if((*it)->status()==121) incomingParton = (*it);
00702   
00703   // find target parton (first entry with IST=122)
00704   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00705       (it != event()->particles_end() && targetParton==NULL); it++)
00706     if((*it)->status()==122) targetParton = (*it);
00707   
00708   // find incoming Proton (first entry ID=2212, IST=101)
00709   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00710       (it != event()->particles_end() && incomingProton==NULL); it++)
00711     if((*it)->status()==101 && (*it)->pdg_id()==2212) incomingProton = (*it);
00712   
00713   // find target Proton (first entry ID=2212, IST=102)
00714   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00715       (it != event()->particles_end() && targetProton==NULL); it++)
00716     if((*it)->status()==102 && (*it)->pdg_id()==2212) targetProton = (*it);
00717   
00718   // find hard scale Q (computed from colliding partons)
00719   if( incomingParton && targetParton ) {
00720     math::XYZTLorentzVector totMomentum(0,0,0,0);
00721     totMomentum+=incomingParton->momentum();
00722     totMomentum+=targetParton->momentum();
00723     double evScale = totMomentum.mass(); 
00724     double evScale2 = evScale*evScale; 
00725     
00726     // find alpha_QED & alpha_QCD         
00727     int one=1;
00728     double alphaQCD=hwualf_(&one,&evScale);
00729     double alphaQED=hwuaem_(&evScale2);
00730     
00731     if(!externalPartons || event()->event_scale() < 0) event()->set_event_scale(evScale);
00732     if(!externalPartons || event()->alphaQCD() < 0)    event()->set_alphaQCD(alphaQCD);
00733     if(!externalPartons || event()->alphaQED() < 0)    event()->set_alphaQED(alphaQED);
00734     
00735     if(!externalPartons || pdfInfo.x1() < 0) {
00736       // get the PDF information
00737       pdfInfo.set_id1( incomingParton->pdg_id()==21 ? 0 : incomingParton->pdg_id());
00738       pdfInfo.set_id2( targetParton->pdg_id()==21 ? 0 : targetParton->pdg_id());
00739       if( incomingProton && targetProton ) {
00740         double x1 = incomingParton->momentum().pz()/incomingProton->momentum().pz();
00741         double x2 = targetParton->momentum().pz()/targetProton->momentum().pz();        
00742         pdfInfo.set_x1(x1);
00743         pdfInfo.set_x2(x2);             
00744       }
00745       // we do not fill pdf1 & pdf2, since they are not easily available (what are they needed for anyways???)
00746       pdfInfo.set_scalePDF(evScale); // the same as Q above... does this make sense?
00747     } 
00748     
00749     if(!externalPartons || event()->signal_process_id() < 0) event()->set_signal_process_id( abs(hwproc.IPROC) );
00750     event()->set_pdf_info( pdfInfo );
00751   }
00752   
00753   // add event weight & PDF information
00754   if (lheRunInfo() != 0 && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
00755     // in LHE weighting mode 4 the weight is an xsec, so convert form HERWIG
00756     // to standard CMS unit "picobarn"
00757     event()->weights().push_back( 1.0e3 * hwevnt.EVWGT );
00758   else
00759     event()->weights().push_back( hwevnt.EVWGT );
00760 
00761     
00762   // find final parton (first entry with IST=123)
00763   HepMC::GenParticle* finalParton = NULL;
00764   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00765       (it != event()->particles_end() && finalParton==NULL); it++)
00766     if((*it)->status()==123) finalParton = (*it);
00767   
00768   
00769   // add GenEventInfo & binning Values
00770   eventInfo().reset(new GenEventInfoProduct(event().get()));    
00771   if(finalParton) {
00772     double thisPt=finalParton->momentum().perp();
00773     eventInfo()->setBinningValues(std::vector<double>(1, thisPt));
00774   }
00775   
00776   // emulate PY6 status codes, if switched on...
00777   if (emulatePythiaStatusCodes)
00778     pythiaStatusCodes();          
00779   
00780 }
00781 
00782 
00783 bool Herwig6Hadronizer::decay()
00784 {
00785         // hadron decays
00786 
00787         InstanceWrapper wrapper(this);  // safe guard
00788 
00789         hwdhad();       // unstable particle decays
00790         hwdhvy();       // heavy flavour decays
00791         hwmevt();       // soft underlying event
00792 
00793         hwufne();       // finalize event
00794 
00795         if (hwevnt.IERROR)
00796                 return false;
00797 
00798         event().reset(new HepMC::GenEvent);
00799         if (!conv.fill_next_event(event().get()))
00800                 throw cms::Exception("Herwig6Error")
00801                         << "HepMC Conversion problems in event." << std::endl;
00802 
00803         // do particle ID conversion Herwig->PDG, if requested
00804         if ( fConvertToPDG ) {
00805            for ( HepMC::GenEvent::particle_iterator part = event()->particles_begin(); part != event()->particles_end(); ++part) {
00806              if ((*part)->pdg_id() != HepPID::translateHerwigtoPDT((*part)->pdg_id()))
00807                (*part)->set_pdg_id(HepPID::translateHerwigtoPDT((*part)->pdg_id()));
00808            }
00809         }
00810 
00811         return true;
00812 }
00813 
00814 bool Herwig6Hadronizer::residualDecay()
00815 {
00816         return true;
00817 }
00818 
00819 void Herwig6Hadronizer::upInit()
00820 {
00821         lhef::CommonBlocks::fillHEPRUP(lheRunInfo()->getHEPRUP());
00822         heprup_.pdfgup[0] = heprup_.pdfgup[1] = -1;
00823         heprup_.pdfsup[0] = heprup_.pdfsup[1] = -1;
00824         // we set up the PDFs ourselves
00825         
00826         // pass HERWIG paramaters fomr header (if present)
00827         std::string mcnloHeader="herwig6header";
00828         std::vector<lhef::LHERunInfo::Header> headers=lheRunInfo()->getHeaders();
00829         for(std::vector<lhef::LHERunInfo::Header>::const_iterator hIter=headers.begin();hIter!=headers.end(); ++hIter) {
00830           if(hIter->tag()==mcnloHeader){
00831             readMCatNLOfile=true;
00832             for(lhef::LHERunInfo::Header::const_iterator lIter=hIter->begin(); lIter != hIter->end(); ++lIter) {
00833               if((lIter->c_str())[0]!='#' && (lIter->c_str())[0]!='\n') {   // it's not a comment) 
00834                 if (!give(*lIter))
00835                   throw edm::Exception(edm::errors::Configuration)
00836                     << "Herwig 6 did not accept the following: \""
00837                     << *lIter << "\"." << std::endl;
00838               }
00839             }
00840           }
00841         }
00842 }
00843 
00844 void Herwig6Hadronizer::upEvnt()
00845 {
00846         lhef::CommonBlocks::fillHEPEUP(lheEvent()->getHEPEUP());
00847 
00848         // if MCatNLO external file is read, read comment & pass IHPRO to HERWIG
00849         if(readMCatNLOfile) {
00850           for(std::vector<std::string>::const_iterator iter=lheEvent()->getComments().begin();
00851               iter!=lheEvent()->getComments().end(); ++iter) {
00852             std::string toParse(iter->substr(1));
00853             if (!give(toParse))
00854               throw edm::Exception(edm::errors::Configuration)
00855                     << "Herwig 6 did not accept the following: \""
00856                     << toParse << "\"." << std::endl;
00857           }
00858         }
00859         
00860 }
00861 
00862 int Herwig6Hadronizer::pythiaStatusCode(const HepMC::GenParticle *p) const
00863 {
00864         int status = p->status();
00865 
00866         // weird 9922212 particles...
00867         if (status == 3 && !p->end_vertex())
00868                 return 2;
00869 
00870         if (status >= 1 && status <= 3)
00871                 return status;
00872 
00873         if (!p->end_vertex())
00874                 return 1;
00875 
00876         // let's prevent particles having status 3, if the identical
00877         // particle downstream is a better status 3 candidate
00878         int currentId = p->pdg_id();
00879         int orig = status;
00880         if (status == 123 || status == 124 ||
00881             status == 155 || status == 156 || status == 160 ||
00882             (status >= 195 && status <= 197)) {
00883                 for(const HepMC::GenParticle *q = p;;) {
00884                         const HepMC::GenVertex *vtx = q->end_vertex();
00885                         if (!vtx)
00886                                 break;
00887 
00888                         HepMC::GenVertex::particles_out_const_iterator iter;
00889                         for(iter = vtx->particles_out_const_begin();
00890                             iter != vtx->particles_out_const_end(); ++iter)
00891                                 if ((*iter)->pdg_id() == currentId)
00892                                         break;
00893 
00894                         if (iter == vtx->particles_out_const_end())
00895                                 break;
00896 
00897                         q = *iter;
00898                         if (q->status() == 3 ||
00899                             ((status == 120 || status == 123 ||
00900                               status == 124) && orig > 124))
00901                                 return 4;
00902                 }
00903         }
00904 
00905         int nesting = 0;
00906         for(;;) {
00907                 if ((status >= 120 && status <= 122) || status == 3) {
00908                         // avoid flagging status 3 if there is a
00909                         // better status 3 candidate upstream
00910                         if (externalPartons)
00911                                 return ((orig >= 121 && orig <= 124) ||
00912                                         orig == 3) ? 3 : 4;
00913                         else
00914                                 return (nesting ||
00915                                         (status != 3 && orig <= 124)) ? 3 : 4;
00916                 }
00917 
00918                 // check whether we are leaving the hard process
00919                 // including heavy resonance decays
00920                 if (!(status == 4 || status == 123 || status == 124 ||
00921                       status == 155 || status == 156 || status == 160 ||
00922                       (status >= 195 && status <= 197)))
00923                         break;
00924 
00925                 const HepMC::GenVertex *vtx = p->production_vertex();
00926                 if (!vtx || !vtx->particles_in_size())
00927                         break;
00928 
00929                 p = *vtx->particles_in_const_begin();
00930                 status = p->status();
00931 
00932                 int newId = p->pdg_id();
00933 
00934                 if (!newId)
00935                         break;
00936 
00937                 // nesting increases if we move to the next-best mother
00938                 if (newId != currentId) {
00939                         if (++nesting > 1 && externalPartons)
00940                                 break;
00941                         currentId = newId;
00942                 }
00943         }
00944 
00945         return 2; 
00946 }
00947 
00948 void Herwig6Hadronizer::pythiaStatusCodes()
00949 {
00950         for(HepMC::GenEvent::particle_iterator iter =
00951                                                 event()->particles_begin();
00952             iter != event()->particles_end(); iter++)
00953                 (*iter)->set_status(pythiaStatusCode(*iter));
00954 
00955         for(HepMC::GenEvent::particle_iterator iter =
00956                                                 event()->particles_begin();
00957             iter != event()->particles_end(); iter++)
00958                 if ((*iter)->status() == 4)
00959                         (*iter)->set_status(2);
00960 }
00961 
00962 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00963 
00964 typedef edm::GeneratorFilter<Herwig6Hadronizer, gen::ExternalDecayDriver> Herwig6GeneratorFilter;
00965 DEFINE_FWK_MODULE(Herwig6GeneratorFilter);
00966 
00967 typedef edm::HadronizerFilter<Herwig6Hadronizer, gen::ExternalDecayDriver> Herwig6HadronizerFilter;
00968 DEFINE_FWK_MODULE(Herwig6HadronizerFilter);