CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/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/HerwigWrapper6_4.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 initializeForInternalPartons() { return initialize(0); }
00085         bool initializeForExternalPartons() { return initialize(lheRunInfo()->getHEPRUP()); }
00086 
00087         bool declareStableParticles( const std::vector<int> &pdgIds);
00088         bool declareSpecialSettings( const std::vector<std::string> );
00089 
00090         void statistics();
00091 
00092 
00093         bool generatePartonsAndHadronize() { return hadronize(); }
00094         bool hadronize();
00095         bool decay();
00096         bool residualDecay();
00097         void finalizeEvent();
00098 
00099         const char *classname() const { return "Herwig6Hadronizer"; }
00100 
00101     private:
00102         void clear();
00103 
00104         int pythiaStatusCode(const HepMC::GenParticle *p) const;
00105         void pythiaStatusCodes();
00106 
00107         virtual void upInit();
00108         virtual void upEvnt();
00109 
00110         HepMC::IO_HERWIG                conv;
00111         bool                            needClear;
00112         bool                            externalPartons;
00113 
00114         gen::ParameterCollector         parameters;
00115         int                             herwigVerbosity;
00116         int                             hepmcVerbosity;
00117         int                             maxEventsToPrint;
00118         bool                            printCards;
00119         bool                            emulatePythiaStatusCodes;
00120         double                          comEnergy;
00121         bool                            useJimmy;
00122         bool                            doMPInteraction;
00123         int                             numTrials;
00124         bool                            fConvertToPDG;  // convert PIDs 
00125 
00126         bool                            readMCatNLOfile;
00127 
00128 };
00129 
00130 extern "C" {
00131         void hwwarn_(const char *method, int *id);
00132         void setherwpdf_(void);
00133         void mysetpdfpath_(const char *path);
00134 } // extern "C"
00135 
00136 Herwig6Hadronizer::Herwig6Hadronizer(const edm::ParameterSet &params) :
00137         BaseHadronizer(params),
00138         needClear(false),
00139         parameters(params.getParameter<edm::ParameterSet>("HerwigParameters")),
00140         herwigVerbosity(params.getUntrackedParameter<int>("herwigVerbosity", 0)),
00141         hepmcVerbosity(params.getUntrackedParameter<int>("hepmcVerbosity", 0)),
00142         maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
00143         printCards(params.getUntrackedParameter<bool>("printCards", false)),
00144         emulatePythiaStatusCodes(params.getUntrackedParameter<bool>("emulatePythiaStatusCodes", false)),
00145         comEnergy(params.getParameter<double>("comEnergy")),
00146         useJimmy(params.getParameter<bool>("useJimmy")),
00147         doMPInteraction(params.getParameter<bool>("doMPInteraction")),
00148         numTrials(params.getUntrackedParameter<int>("numTrialsMPI", 100)),
00149         readMCatNLOfile(false)
00150 {
00151   
00152   fConvertToPDG = false;
00153   if ( params.exists( "doPDGConvert" ) )
00154     fConvertToPDG = params.getParameter<bool>("doPDGConvert");
00155 }
00156 
00157 Herwig6Hadronizer::~Herwig6Hadronizer()
00158 {
00159         clear();
00160 }
00161 
00162 void Herwig6Hadronizer::clear()
00163 {
00164         if (!needClear)
00165                 return;
00166 
00167         // teminate elementary process
00168         call(hwefin);
00169         if (useJimmy)
00170                 call(jmefin);
00171 
00172         needClear = false;
00173 }
00174 
00175 void Herwig6Hadronizer::setSLHAFromHeader(
00176                                 const std::vector<std::string> &lines)
00177 {
00178         std::set<std::string> blocks;
00179         std::string block;
00180         for(std::vector<std::string>::const_iterator iter = lines.begin();
00181             iter != lines.end(); ++iter) {
00182                 std::string line = *iter;
00183                 std::transform(line.begin(), line.end(),
00184                                line.begin(), (int(*)(int))std::toupper);
00185                 std::string::size_type pos = line.find('#');
00186                 if (pos != std::string::npos)
00187                         line.resize(pos);
00188 
00189                 if (line.empty())
00190                         continue;
00191 
00192                 if (!boost::algorithm::is_space()(line[0])) {
00193                         std::vector<std::string> tokens;
00194                         boost::split(tokens, line,
00195                                      boost::algorithm::is_space(),
00196                                      boost::token_compress_on);
00197                         if (!tokens.size())
00198                                 continue;
00199                         block.clear();
00200                         if (tokens.size() < 2)
00201                                 continue;
00202                         if (tokens[0] == "BLOCK")
00203                                 block = tokens[1];
00204                         else if (tokens[0] == "DECAY")
00205                                 block = "DECAY";
00206 
00207                         if (block.empty())
00208                                 continue;
00209 
00210                         if (!blocks.count(block)) {
00211                                 blocks.insert(block);
00212                                 edm::LogWarning("Generator|Herwig6Hadronzier")
00213                                         << "Unsupported SLHA block \"" << block
00214                                         << "\".  It will be ignored."
00215                                         << std::endl;
00216                         }
00217                 }
00218         }
00219 }
00220 
00221 bool Herwig6Hadronizer::initialize(const lhef::HEPRUP *heprup)
00222 {
00223         clear();
00224 
00225         externalPartons = (heprup != 0);
00226 
00227         std::ostringstream info;
00228         info << "---------------------------------------------------\n"; 
00229         info << "Initializing Herwig6Hadronizer for "
00230              << (externalPartons ? "external" : "internal") << " partons\n";
00231         info << "---------------------------------------------------\n";
00232 
00233         info << "   Herwig verbosity level         = " << herwigVerbosity << "\n";
00234         info << "   HepMC verbosity                = " << hepmcVerbosity << "\n";
00235         info << "   Number of events to be printed = " << maxEventsToPrint << "\n";
00236 
00237         // Call hwudat to set up HERWIG block data
00238         hwudat();
00239 
00240         // Setting basic parameters
00241         if (externalPartons) {
00242                 hwproc.PBEAM1 = heprup->EBMUP.first;
00243                 hwproc.PBEAM2 = heprup->EBMUP.second;
00244                 pdgToHerwig(heprup->IDBMUP.first, hwbmch.PART1);
00245                 pdgToHerwig(heprup->IDBMUP.second, hwbmch.PART2);
00246         } else {
00247                 hwproc.PBEAM1 = 0.5 * comEnergy;
00248                 hwproc.PBEAM2 = 0.5 * comEnergy;
00249                 pdgToHerwig(2212, hwbmch.PART1);
00250                 pdgToHerwig(2212, hwbmch.PART2);
00251         }
00252 
00253         if (useJimmy) {
00254                 info << "   HERWIG will be using JIMMY for UE/MI.\n";
00255                 jmparm.MSFLAG = 1;
00256                 if (doMPInteraction)
00257                         info << "   JIMMY trying to generate multiple interactions.\n";
00258         }
00259 
00260 
00261         // set the IPROC already here... needed for VB pairs
00262 
00263         bool iprocFound=false;
00264         
00265         for(gen::ParameterCollector::const_iterator line = parameters.begin();
00266             line != parameters.end(); ++line) {
00267           if(!strcmp((line->substr(0,5)).c_str(),"IPROC")) {
00268             if (!give(*line))
00269               throw edm::Exception(edm::errors::Configuration)
00270                 << "Herwig 6 did not accept the following: \""
00271                 << *line << "\"." << std::endl;
00272             else iprocFound=true;
00273           }
00274         }
00275         
00276         if (!iprocFound && !externalPartons)
00277           throw edm::Exception(edm::errors::Configuration)
00278             << "You have to define the process with IPROC."  << std::endl;
00279 
00280         // initialize other common blocks ...
00281         call(hwigin);
00282         hwevnt.MAXER = 100000000;       // O(inf)
00283         hwpram.LWSUD = 0;               // don't write Sudakov form factors
00284         hwdspn.LWDEC = 0;               // don't write three/four body decays
00285                                         // (no fort.77 and fort.88 ...)
00286 
00287         // init LHAPDF glue
00288 
00289         std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
00290         for(unsigned int i = 0; i < 2; i++) {
00291                 hwpram.MODPDF[i] = -111;
00292                 std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
00293         }
00294 
00295         if (useJimmy)
00296                 call(jimmin);
00297 
00298         hwevnt.MAXPR = maxEventsToPrint;
00299         hwpram.IPRINT = herwigVerbosity;
00300 //      hwprop.RMASS[6] = 175.0;        //FIXME
00301 
00302         if (printCards) {
00303                 info << "\n";
00304                 info << "------------------------------------\n";
00305                 info << "Reading HERWIG parameters\n";
00306                 info << "------------------------------------\n";
00307     
00308         }
00309         for(gen::ParameterCollector::const_iterator line = parameters.begin();
00310             line != parameters.end(); ++line) {
00311                 if (printCards)
00312                         info << "   " << *line << "\n";
00313                 if (!give(*line))
00314                         throw edm::Exception(edm::errors::Configuration)
00315                                 << "Herwig 6 did not accept the following: \""
00316                                 << *line << "\"." << std::endl;
00317         }
00318 
00319         if (printCards)
00320                 info << "\n";
00321 
00322         if (externalPartons) {
00323                 std::vector<std::string> slha =
00324                                 lheRunInfo()->findHeader("slha");
00325                 if (!slha.empty())
00326                         setSLHAFromHeader(slha);
00327         }
00328 
00329         needClear = true;
00330 
00331         std::pair<int, int> pdfs(-1, -1);
00332         if (externalPartons)
00333                 pdfs = lheRunInfo()->pdfSetTranslation();
00334 
00335         if (hwpram.MODPDF[0] != -111 || hwpram.MODPDF[1] != -111) {
00336                 for(unsigned int i = 0; i < 2; i++)
00337                         if (hwpram.MODPDF[i] == -111)
00338                                 hwpram.MODPDF[i] = -1;
00339 
00340                 if (pdfs.first != -1 || pdfs.second != -1)
00341                         edm::LogError("Generator|Herwig6Hadronzier")
00342                                 << "Both external Les Houches event and "
00343                                    "config file specify a PDF set.  "
00344                                    "User PDF will override external one."
00345                                 << std::endl;
00346 
00347                 pdfs.first = hwpram.MODPDF[0] != -111 ? hwpram.MODPDF[0] : -1;
00348                 pdfs.second = hwpram.MODPDF[1] != -111 ? hwpram.MODPDF[1] : -1;
00349         }
00350 
00351         hwpram.MODPDF[0] = pdfs.first;
00352         hwpram.MODPDF[1] = pdfs.second;
00353 
00354         if (externalPartons)
00355                 hwproc.IPROC = -1;
00356 
00357         // HERWIG preparations ...
00358         call(hwuinc);
00359         markStable(13);          // MU+
00360         markStable(-13);         // MU-
00361         markStable(3112);        // SIGMA+
00362         markStable(-3112);       // SIGMABAR+
00363         markStable(3222);        // SIGMA-
00364         markStable(-3222);       // SIGMABAR-
00365         markStable(3122);        // LAMBDA0
00366         markStable(-3122);       // LAMBDABAR0
00367         markStable(3312);        // XI-
00368         markStable(-3312);       // XIBAR+
00369         markStable(3322);        // XI0
00370         markStable(-3322);       // XI0BAR
00371         markStable(3334);        // OMEGA-
00372         markStable(-3334);       // OMEGABAR+
00373         markStable(211);         // PI+
00374         markStable(-211);        // PI-
00375         markStable(321);         // K+
00376         markStable(-321);        // K-
00377         markStable(310);         // K_S0
00378         markStable(130);         // K_L0
00379 
00380         // better: merge with declareStableParticles
00381         // and get the list from configuration / Geant4 / Core somewhere
00382 
00383         // initialize HERWIG event generation
00384         call(hweini);
00385 
00386         if (useJimmy)
00387                 call(jminit);
00388 
00389         edm::LogInfo(info.str());
00390 
00391         return true;
00392 }
00393 
00394 bool Herwig6Hadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00395 {
00396         for(std::vector<int>::const_iterator iter = pdgIds.begin();
00397             iter != pdgIds.end(); ++iter)
00398                 if (!markStable(*iter))
00399                         return false;
00400         return true;
00401 }
00402 
00403 bool Herwig6Hadronizer::declareSpecialSettings( const std::vector<std::string> )
00404 {
00405    return true;
00406 }
00407 
00408 void Herwig6Hadronizer::statistics()
00409 {
00410         if (!runInfo().internalXSec()) {
00411           // not set via LHE, so get it from HERWIG
00412           // the reason is that HERWIG doesn't compute the xsec
00413           // in all LHE modes
00414 
00415           double RNWGT = 1. / hwevnt.NWGTS;
00416           double AVWGT = hwevnt.WGTSUM * RNWGT;
00417 
00418           double xsec = 1.0e3 * AVWGT;
00419 
00420           runInfo().setInternalXSec(xsec);
00421         }
00422 }
00423 
00424 bool Herwig6Hadronizer::hadronize()
00425 {
00426         // hard process generation, parton shower, hadron formation
00427 
00428         InstanceWrapper wrapper(this);  // safe guard
00429 
00430         event().reset();
00431 
00432         // call herwig routines to create HEPEVT
00433 
00434         hwuine();       // initialize event
00435         
00436         if (callWithTimeout(10, hwepro)) { // process event and PS
00437           // We hung for more than 10 seconds
00438           int error = 199;
00439           hwwarn_("HWHGUP", &error);
00440         }
00441         
00442         hwbgen();       // parton cascades
00443 
00444         // call jimmy ... only if event is not killed yet by HERWIG
00445         if (useJimmy && doMPInteraction && !hwevnt.IERROR && call_hwmsct()) {
00446           if (lheEvent()) lheEvent()->count(lhef::LHERunInfo::kKilled);
00447           return false;
00448         }
00449         
00450         hwdhob();       // heavy quark decays
00451         hwcfor();       // cluster formation
00452         hwcdec();       // cluster decays
00453         
00454         // if event *not* killed by HERWIG, return true
00455         if (hwevnt.IERROR) {
00456           hwufne();     // finalize event, to keep system clean
00457           if (lheEvent()) lheEvent()->count(lhef::LHERunInfo::kKilled);
00458           return false;
00459         }
00460 
00461         if (lheEvent()) lheEvent()->count(lhef::LHERunInfo::kAccepted);
00462         return true;
00463 }
00464 
00465 void Herwig6Hadronizer::finalizeEvent()
00466 {
00467   lhef::LHEEvent::fixHepMCEventTimeOrdering(event().get());
00468   
00469   HepMC::PdfInfo pdfInfo;
00470   if(externalPartons) {
00471     lheEvent()->fillEventInfo( event().get() );
00472     lheEvent()->fillPdfInfo( &pdfInfo );
00473 
00474     // for MC@NLO: IDWRUP is not filled...
00475     if(event()->signal_process_id()==0)
00476       event()->set_signal_process_id( abs(hwproc.IPROC) );
00477 
00478   }
00479     
00480   HepMC::GenParticle* incomingParton = NULL;
00481   HepMC::GenParticle* targetParton = NULL;
00482   
00483   HepMC::GenParticle* incomingProton = NULL;
00484   HepMC::GenParticle* targetProton = NULL;
00485   
00486   // find incoming parton (first entry with IST=121)
00487   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00488       (it != event()->particles_end() && incomingParton==NULL); it++)
00489     if((*it)->status()==121) incomingParton = (*it);
00490   
00491   // find target parton (first entry with IST=122)
00492   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00493       (it != event()->particles_end() && targetParton==NULL); it++)
00494     if((*it)->status()==122) targetParton = (*it);
00495   
00496   // find incoming Proton (first entry ID=2212, IST=101)
00497   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00498       (it != event()->particles_end() && incomingProton==NULL); it++)
00499     if((*it)->status()==101 && (*it)->pdg_id()==2212) incomingProton = (*it);
00500   
00501   // find target Proton (first entry ID=2212, IST=102)
00502   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00503       (it != event()->particles_end() && targetProton==NULL); it++)
00504     if((*it)->status()==102 && (*it)->pdg_id()==2212) targetProton = (*it);
00505   
00506   // find hard scale Q (computed from colliding partons)
00507   if( incomingParton && targetParton ) {
00508     math::XYZTLorentzVector totMomentum(0,0,0,0);
00509     totMomentum+=incomingParton->momentum();
00510     totMomentum+=targetParton->momentum();
00511     double evScale = totMomentum.mass(); 
00512     double evScale2 = evScale*evScale; 
00513     
00514     // find alpha_QED & alpha_QCD         
00515     int one=1;
00516     double alphaQCD=hwualf_(&one,&evScale);
00517     double alphaQED=hwuaem_(&evScale2);
00518     
00519     if(!externalPartons || event()->event_scale() < 0) event()->set_event_scale(evScale);
00520     if(!externalPartons || event()->alphaQCD() < 0)    event()->set_alphaQCD(alphaQCD);
00521     if(!externalPartons || event()->alphaQED() < 0)    event()->set_alphaQED(alphaQED);
00522     
00523     if(!externalPartons || pdfInfo.x1() < 0) {
00524       // get the PDF information
00525       pdfInfo.set_id1( incomingParton->pdg_id()==21 ? 0 : incomingParton->pdg_id());
00526       pdfInfo.set_id2( targetParton->pdg_id()==21 ? 0 : targetParton->pdg_id());
00527       if( incomingProton && targetProton ) {
00528         double x1 = incomingParton->momentum().pz()/incomingProton->momentum().pz();
00529         double x2 = targetParton->momentum().pz()/targetProton->momentum().pz();        
00530         pdfInfo.set_x1(x1);
00531         pdfInfo.set_x2(x2);             
00532       }
00533       // we do not fill pdf1 & pdf2, since they are not easily available (what are they needed for anyways???)
00534       pdfInfo.set_scalePDF(evScale); // the same as Q above... does this make sense?
00535     } 
00536     
00537     if(!externalPartons || event()->signal_process_id() < 0) event()->set_signal_process_id( abs(hwproc.IPROC) );
00538     event()->set_pdf_info( pdfInfo );
00539   }
00540   
00541   // add event weight & PDF information
00542   if (lheRunInfo() != 0 && std::abs(lheRunInfo()->getHEPRUP()->IDWTUP) == 4)
00543     // in LHE weighting mode 4 the weight is an xsec, so convert form HERWIG
00544     // to standard CMS unit "picobarn"
00545     event()->weights().push_back( 1.0e3 * hwevnt.EVWGT );
00546   else
00547     event()->weights().push_back( hwevnt.EVWGT );
00548 
00549     
00550   // find final parton (first entry with IST=123)
00551   HepMC::GenParticle* finalParton = NULL;
00552   for(HepMC::GenEvent::particle_const_iterator it = event()->particles_begin(); 
00553       (it != event()->particles_end() && finalParton==NULL); it++)
00554     if((*it)->status()==123) finalParton = (*it);
00555   
00556   
00557   // add GenEventInfo & binning Values
00558   eventInfo().reset(new GenEventInfoProduct(event().get()));    
00559   if(finalParton) {
00560     double thisPt=finalParton->momentum().perp();
00561     eventInfo()->setBinningValues(std::vector<double>(1, thisPt));
00562   }
00563   
00564   // emulate PY6 status codes, if switched on...
00565   if (emulatePythiaStatusCodes)
00566     pythiaStatusCodes();          
00567   
00568 }
00569 
00570 
00571 bool Herwig6Hadronizer::decay()
00572 {
00573         // hadron decays
00574 
00575         InstanceWrapper wrapper(this);  // safe guard
00576 
00577         hwdhad();       // unstable particle decays
00578         hwdhvy();       // heavy flavour decays
00579         hwmevt();       // soft underlying event
00580 
00581         hwufne();       // finalize event
00582 
00583         if (hwevnt.IERROR)
00584                 return false;
00585 
00586         event().reset(new HepMC::GenEvent);
00587         if (!conv.fill_next_event(event().get()))
00588                 throw cms::Exception("Herwig6Error")
00589                         << "HepMC Conversion problems in event." << std::endl;
00590 
00591         // do particle ID conversion Herwig->PDG, if requested
00592         if ( fConvertToPDG ) {
00593            for ( HepMC::GenEvent::particle_iterator part = event()->particles_begin(); part != event()->particles_end(); ++part) {
00594              if ((*part)->pdg_id() != HepPID::translateHerwigtoPDT((*part)->pdg_id()))
00595                (*part)->set_pdg_id(HepPID::translateHerwigtoPDT((*part)->pdg_id()));
00596            }
00597         }
00598 
00599         return true;
00600 }
00601 
00602 bool Herwig6Hadronizer::residualDecay()
00603 {
00604         return true;
00605 }
00606 
00607 void Herwig6Hadronizer::upInit()
00608 {
00609         lhef::CommonBlocks::fillHEPRUP(lheRunInfo()->getHEPRUP());
00610         heprup_.pdfgup[0] = heprup_.pdfgup[1] = -1;
00611         heprup_.pdfsup[0] = heprup_.pdfsup[1] = -1;
00612         // we set up the PDFs ourselves
00613         
00614         // pass HERWIG paramaters fomr header (if present)
00615         std::string mcnloHeader="herwig6header";
00616         std::vector<lhef::LHERunInfo::Header> headers=lheRunInfo()->getHeaders();
00617         for(std::vector<lhef::LHERunInfo::Header>::const_iterator hIter=headers.begin();hIter!=headers.end(); ++hIter) {
00618           if(hIter->tag()==mcnloHeader){
00619             readMCatNLOfile=true;
00620             for(lhef::LHERunInfo::Header::const_iterator lIter=hIter->begin(); lIter != hIter->end(); ++lIter) {
00621               if((lIter->c_str())[0]!='#' && (lIter->c_str())[0]!='\n') {   // it's not a comment) 
00622                 if (!give(*lIter))
00623                   throw edm::Exception(edm::errors::Configuration)
00624                     << "Herwig 6 did not accept the following: \""
00625                     << *lIter << "\"." << std::endl;
00626               }
00627             }
00628           }
00629         }
00630 }
00631 
00632 void Herwig6Hadronizer::upEvnt()
00633 {
00634         lhef::CommonBlocks::fillHEPEUP(lheEvent()->getHEPEUP());
00635 
00636         // if MCatNLO external file is read, read comment & pass IHPRO to HERWIG
00637         if(readMCatNLOfile) {
00638           for(std::vector<std::string>::const_iterator iter=lheEvent()->getComments().begin();
00639               iter!=lheEvent()->getComments().end(); ++iter) {
00640             std::string toParse(iter->substr(1));
00641             if (!give(toParse))
00642               throw edm::Exception(edm::errors::Configuration)
00643                     << "Herwig 6 did not accept the following: \""
00644                     << toParse << "\"." << std::endl;
00645           }
00646         }
00647         
00648 }
00649 
00650 int Herwig6Hadronizer::pythiaStatusCode(const HepMC::GenParticle *p) const
00651 {
00652         int status = p->status();
00653 
00654         // weird 9922212 particles...
00655         if (status == 3 && !p->end_vertex())
00656                 return 2;
00657 
00658         if (status >= 1 && status <= 3)
00659                 return status;
00660 
00661         if (!p->end_vertex())
00662                 return 1;
00663 
00664         // let's prevent particles having status 3, if the identical
00665         // particle downstream is a better status 3 candidate
00666         int currentId = p->pdg_id();
00667         int orig = status;
00668         if (status == 123 || status == 124 ||
00669             status == 155 || status == 156 || status == 160 ||
00670             (status >= 195 && status <= 197)) {
00671                 for(const HepMC::GenParticle *q = p;;) {
00672                         const HepMC::GenVertex *vtx = q->end_vertex();
00673                         if (!vtx)
00674                                 break;
00675 
00676                         HepMC::GenVertex::particles_out_const_iterator iter;
00677                         for(iter = vtx->particles_out_const_begin();
00678                             iter != vtx->particles_out_const_end(); ++iter)
00679                                 if ((*iter)->pdg_id() == currentId)
00680                                         break;
00681 
00682                         if (iter == vtx->particles_out_const_end())
00683                                 break;
00684 
00685                         q = *iter;
00686                         if (q->status() == 3 ||
00687                             ((status == 120 || status == 123 ||
00688                               status == 124) && orig > 124))
00689                                 return 4;
00690                 }
00691         }
00692 
00693         int nesting = 0;
00694         for(;;) {
00695                 if ((status >= 120 && status <= 122) || status == 3) {
00696                         // avoid flagging status 3 if there is a
00697                         // better status 3 candidate upstream
00698                         if (externalPartons)
00699                                 return ((orig >= 121 && orig <= 124) ||
00700                                         orig == 3) ? 3 : 4;
00701                         else
00702                                 return (nesting ||
00703                                         (status != 3 && orig <= 124)) ? 3 : 4;
00704                 }
00705 
00706                 // check whether we are leaving the hard process
00707                 // including heavy resonance decays
00708                 if (!(status == 4 || status == 123 || status == 124 ||
00709                       status == 155 || status == 156 || status == 160 ||
00710                       (status >= 195 && status <= 197)))
00711                         break;
00712 
00713                 const HepMC::GenVertex *vtx = p->production_vertex();
00714                 if (!vtx || !vtx->particles_in_size())
00715                         break;
00716 
00717                 p = *vtx->particles_in_const_begin();
00718                 status = p->status();
00719 
00720                 int newId = p->pdg_id();
00721 
00722                 if (!newId)
00723                         break;
00724 
00725                 // nesting increases if we move to the next-best mother
00726                 if (newId != currentId) {
00727                         if (++nesting > 1 && externalPartons)
00728                                 break;
00729                         currentId = newId;
00730                 }
00731         }
00732 
00733         return 2; 
00734 }
00735 
00736 void Herwig6Hadronizer::pythiaStatusCodes()
00737 {
00738         for(HepMC::GenEvent::particle_iterator iter =
00739                                                 event()->particles_begin();
00740             iter != event()->particles_end(); iter++)
00741                 (*iter)->set_status(pythiaStatusCode(*iter));
00742 
00743         for(HepMC::GenEvent::particle_iterator iter =
00744                                                 event()->particles_begin();
00745             iter != event()->particles_end(); iter++)
00746                 if ((*iter)->status() == 4)
00747                         (*iter)->set_status(2);
00748 }
00749 
00750 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00751 
00752 typedef edm::GeneratorFilter<Herwig6Hadronizer, gen::ExternalDecayDriver> Herwig6GeneratorFilter;
00753 DEFINE_FWK_MODULE(Herwig6GeneratorFilter);
00754 
00755 typedef edm::HadronizerFilter<Herwig6Hadronizer, gen::ExternalDecayDriver> Herwig6HadronizerFilter;
00756 DEFINE_FWK_MODULE(Herwig6HadronizerFilter);