CMS 3D CMS Logo

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