CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/GeneratorInterface/PomwigInterface/src/PomwigHadronizer.cc

Go to the documentation of this file.
00001 #include "GeneratorInterface/PomwigInterface/interface/PomwigHadronizer.h"
00002 
00003 #include <cstring>
00004 #include <sstream>
00005 #include <string>
00006 #include <vector>
00007 #include <memory>
00008 #include <map>
00009 #include <set>
00010 
00011 #include <boost/shared_ptr.hpp>
00012 #include <boost/algorithm/string/classification.hpp>
00013 #include <boost/algorithm/string/split.hpp>
00014 
00015 #include <HepMC/GenEvent.h>
00016 #include <HepMC/GenParticle.h>
00017 #include <HepMC/GenVertex.h>
00018 #include <HepMC/PdfInfo.h>
00019 #include <HepMC/HerwigWrapper6_4.h>
00020 #include <HepMC/HEPEVT_Wrapper.h>
00021 #include <HepMC/IO_HERWIG.h>
00022 #include "HepPID/ParticleIDTranslations.hh"
00023 
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 
00026 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00027 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00028 
00029 #include "SimDataFormats/GeneratorProducts/interface/LesHouches.h"
00030 #include "SimDataFormats/GeneratorProducts/interface/LHECommonBlocks.h"
00031 
00032 #include "GeneratorInterface/Core/interface/ParameterCollector.h"
00033 #include "GeneratorInterface/Core/interface/BaseHadronizer.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 namespace gen
00041 {
00042 extern "C" {
00043         void hwuidt_(int *iopt, int *ipdg, int *iwig, char nwig[8]);
00044 }
00045 
00046 // helpers
00047 namespace {
00048         static inline bool call_hwmsct()
00049         {
00050                 int result;
00051                 hwmsct(&result);
00052                 return result;
00053         }
00054 
00055         static int pdgToHerwig(int ipdg, char *nwig)
00056         {
00057                 int iopt = 1;
00058                 int iwig = 0;
00059                 hwuidt_(&iopt, &ipdg, &iwig, nwig);
00060                 return ipdg ? iwig : 0;
00061         }
00062 
00063         static bool markStable(int pdgId)
00064         {
00065                 char nwig[9] = "        ";
00066                 if (!pdgToHerwig(pdgId, nwig))
00067                         return false;
00068                 hwusta(nwig, 1);
00069                 return true;
00070         }
00071 }
00072 
00073 #define qcd_1994 qcd_1994_
00074 extern "C" {
00075     void qcd_1994(double&,double&,double*,int&);
00076 }
00077 // For H1 2006 fits
00078 #define qcd_2006 qcd_2006_
00079 extern "C" {
00080     void qcd_2006(double&,double&,int&,double*,double*,double*,double*,double*);
00081 }
00082 
00083 extern "C" {
00084         void hwwarn_(const char *method, int *id);
00085         void setherwpdf_(void);
00086         void mysetpdfpath_(const char *path);
00087 }
00088 
00089 PomwigHadronizer::PomwigHadronizer(const edm::ParameterSet &params) :
00090         BaseHadronizer(params),
00091         needClear(false),
00092         parameters(params.getParameter<edm::ParameterSet>("HerwigParameters")),
00093         herwigVerbosity(params.getUntrackedParameter<int>("herwigVerbosity", 0)),
00094         hepmcVerbosity(params.getUntrackedParameter<int>("hepmcVerbosity", 0)),
00095         maxEventsToPrint(params.getUntrackedParameter<int>("maxEventsToPrint", 0)),
00096         printCards(params.getUntrackedParameter<bool>("printCards", false)),
00097         comEnergy(params.getParameter<double>("comEnergy")),
00098         survivalProbability(params.getParameter<double>("survivalProbability")),
00099         diffTopology(params.getParameter<int>("diffTopology")),
00100         h1fit(params.getParameter<int>("h1fit")),
00101         useJimmy(params.getParameter<bool>("useJimmy")),
00102         doMPInteraction(params.getParameter<bool>("doMPInteraction")),
00103         numTrials(params.getUntrackedParameter<int>("numTrialsMPI", 100)),
00104         doPDGConvert(false)
00105 {
00106         if ( params.exists( "doPDGConvert" ) )
00107            doPDGConvert = params.getParameter<bool>("doPDGConvert");
00108 }
00109 
00110 PomwigHadronizer::~PomwigHadronizer()
00111 {
00112         clear();
00113 }
00114 
00115 void PomwigHadronizer::clear()
00116 {
00117         if (!needClear)
00118                 return;
00119 
00120         // teminate elementary process
00121         call(hwefin);
00122         if (useJimmy)
00123                 call(jmefin);
00124 
00125         needClear = false;
00126 }
00127 
00128 bool PomwigHadronizer::initializeForExternalPartons()
00129 {
00130    return false;
00131 }
00132 
00133 bool PomwigHadronizer::initializeForInternalPartons()
00134 {
00135    clear();
00136 
00137    edm::LogVerbatim("") << "----------------------------------------------\n"
00138                         << "Initializing PomwigHadronizer\n"
00139                         << "----------------------------------------------\n";
00140 
00141    // Call hwudat to set up HERWIG block data
00142    //hwudat();
00143  
00144    // Setting basic parameters ...
00145    hwproc.PBEAM1 = comEnergy/2.;
00146    hwproc.PBEAM2 = comEnergy/2.;
00147    // Choose beam particles for POMWIG depending on topology
00148    switch (diffTopology){
00149         case 0: //DPE
00150                 hwbmch.PART1[0]  = 'E';
00151                 hwbmch.PART1[1]  = '-';
00152                 hwbmch.PART2[0]  = 'E';
00153                 hwbmch.PART2[1]  = '-';
00154                 break;
00155         case 1: //SD survive PART1
00156                 hwbmch.PART1[0]  = 'E';
00157                 hwbmch.PART1[1]  = '-';
00158                 hwbmch.PART2[0]  = 'P';
00159                 hwbmch.PART2[1]  = ' ';
00160                 break;
00161         case 2: //SD survive PART2
00162                 hwbmch.PART1[0]  = 'P';
00163                 hwbmch.PART1[1]  = ' ';
00164                 hwbmch.PART2[0]  = 'E';
00165                 hwbmch.PART2[1]  = '-';
00166                 break;
00167         case 3: //Non diffractive
00168                 hwbmch.PART1[0]  = 'P';
00169                 hwbmch.PART1[1]  = ' ';
00170                 hwbmch.PART2[0]  = 'P';
00171                 hwbmch.PART2[1]  = ' ';
00172                 break;
00173         default:
00174                 throw edm::Exception(edm::errors::Configuration,"PomwigError")
00175           <<" Invalid Diff. Topology. Must be DPE(diffTopology = 0), SD particle 1 (diffTopology = 1), SD particle 2 (diffTopology = 2) and Non diffractive (diffTopology = 3)";
00176                 break;
00177    }
00178    for(int i=2;i<8;++i){
00179     hwbmch.PART1[i]  = ' ';
00180     hwbmch.PART2[i]  = ' ';}
00181 
00182    // initialize other common blocks ...
00183    call(hwigin);
00184 
00185    hwevnt.MAXER = 100000000;    // O(inf)
00186    hwpram.LWSUD = 0;            // don't write Sudakov form factors
00187    hwdspn.LWDEC = 0;            // don't write three/four body decays
00188                                 // (no fort.77 and fort.88 ...)a
00189 
00190    std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
00191    for(unsigned int i = 0; i < 2; i++) {
00192         hwpram.MODPDF[i] = -111;
00193         std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
00194    }
00195 
00196    hwevnt.MAXPR = maxEventsToPrint;
00197    hwpram.IPRINT = herwigVerbosity;
00198 
00199    edm::LogVerbatim("") << "------------------------------------\n"
00200                         << "Reading HERWIG parameters\n"
00201                         << "------------------------------------\n";
00202     
00203    for(gen::ParameterCollector::const_iterator line = parameters.begin();
00204                                                line != parameters.end(); ++line) {
00205       edm::LogVerbatim("") << "   " << *line;
00206       if (!give(*line))
00207           throw edm::Exception(edm::errors::Configuration)
00208                 << "Herwig 6 did not accept the following: \""
00209                 << *line << "\"." << std::endl;
00210    }
00211 
00212    needClear = true;
00213  
00214    call(hwuinc);
00215 
00216    hwusta("PI0     ",1);
00217  
00218    if(!initializeDPDF()) return false;
00219 
00220    call(hweini);
00221 
00222    return true;
00223 }
00224 
00225 bool PomwigHadronizer::initializeDPDF()
00226 {
00227    // Initialize H1 pomeron/reggeon
00228 
00229    if(diffTopology == 3) return true;
00230  
00231    if((diffTopology != 0)&&(diffTopology != 1)&&(diffTopology != 2)) return false;
00232 
00233    int nstru = hwpram.NSTRU;
00234    int ifit = h1fit;
00235    if((nstru == 9)||(nstru == 10)){
00236          if((ifit <= 0)||(ifit >= 7)){
00237               throw edm::Exception(edm::errors::Configuration,"PomwigError")
00238               <<" Attempted to set non existant H1 1997 fit index. Has to be 1...6";
00239          }
00240          std::string aux((nstru == 9)?"Pomeron":"Reggeon"); 
00241          edm::LogVerbatim("") << "   H1 1997 pdf's: " << aux << "\n"
00242                               << "   IFIT = " << ifit;
00243          double xp = 0.1;
00244          double Q2 = 75.0;
00245          double xpq[13];
00246          qcd_1994(xp,Q2,xpq,ifit);
00247     } else if((nstru >= 12)&&(nstru <= 15)){
00248          bool isPom = (nstru == 12)||(nstru == 14);
00249          bool isFitA = (nstru == 12)||(nstru == 13);
00250          ifit = (isFitA)?1:2;
00251          std::string aux_0((isFitA)?"A":"B");
00252          std::string aux_1((isPom)?"Pomeron":"Reggeon"); 
00253          edm::LogVerbatim("") << "   H1 2006 Fit " << aux_0 << " " << aux_1 << "\n"
00254                               << "   IFIT = "<< ifit;
00255          double xp = 0.1;
00256          double Q2 = 75.0;
00257          double xpq[13];
00258          double f2[2];
00259          double fl[2];
00260          double c2[2];
00261          double cl[2];
00262          qcd_2006(xp,Q2,ifit,xpq,f2,fl,c2,cl);
00263     } else{
00264          throw edm::Exception(edm::errors::Configuration,"PomwigError")
00265                <<" Only running Pomeron H1 1997 (NSTRU=9), H1 2006 fit A (NSTRU=12) and H1 2006 fit B (NSTRU=14) or Reggeon H1 1997 (NSTRU=10), H1 2006 fit A (NSTRU=13) and H1 2006 fit B (NSTRU=15)";
00266     }
00267  
00268     return true;
00269 }
00270 
00271 bool PomwigHadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00272 {
00273         for(std::vector<int>::const_iterator iter = pdgIds.begin();
00274             iter != pdgIds.end(); ++iter)
00275                 if (!markStable(*iter))
00276                         return false;
00277         return true;
00278 }
00279 
00280 void PomwigHadronizer::statistics()
00281 {
00282         double RNWGT = 1. / hwevnt.NWGTS;
00283         double AVWGT = hwevnt.WGTSUM * RNWGT;
00284 
00285         double xsec = 1.0e3 * AVWGT;
00286         xsec = survivalProbability*xsec;
00287 
00288         runInfo().setInternalXSec(xsec);
00289 }
00290 
00291 bool PomwigHadronizer::hadronize()
00292 {
00293    return false;
00294 }
00295 
00296 bool PomwigHadronizer::generatePartonsAndHadronize()
00297 {
00298         // hard process generation, parton shower, hadron formation
00299 
00300         InstanceWrapper wrapper(this);  // safe guard
00301 
00302         event().reset();
00303 
00304         // call herwig routines to create HEPEVT
00305 
00306         hwuine();       // initialize event
00307         
00308         if (callWithTimeout(10, hwepro)) { // process event and PS
00309           // We hung for more than 10 seconds
00310           int error = 199;
00311           hwwarn_("HWHGUP", &error);
00312         }
00313         
00314         hwbgen();       // parton cascades
00315 
00316         // call jimmy ... only if event is not killed yet by HERWIG
00317         if (useJimmy && doMPInteraction && !hwevnt.IERROR && call_hwmsct()) 
00318           return false;
00319         
00320         hwdhob();       // heavy quark decays
00321         hwcfor();       // cluster formation
00322         hwcdec();       // cluster decays
00323         
00324         // if event *not* killed by HERWIG, return true
00325         if (!hwevnt.IERROR) return true;
00326 
00327         hwufne();       // finalize event       
00328         return false;
00329 }
00330 
00331 void PomwigHadronizer::finalizeEvent()
00332 {
00333         lhef::LHEEvent::fixHepMCEventTimeOrdering(event().get());
00334 
00335         event()->set_signal_process_id(hwproc.IPROC);
00336 
00337         event()->weights().push_back(hwevnt.EVWGT);
00338 }
00339 
00340 bool PomwigHadronizer::decay()
00341 {
00342         // hadron decays
00343 
00344         InstanceWrapper wrapper(this);  // safe guard
00345 
00346         hwdhad();       // unstable particle decays
00347         hwdhvy();       // heavy flavour decays
00348         hwmevt();       // soft underlying event
00349 
00350         hwufne();       // finalize event
00351 
00352         if (hwevnt.IERROR)
00353                 return false;
00354 
00355         event().reset(new HepMC::GenEvent);
00356         if (!conv.fill_next_event(event().get()))
00357                 throw cms::Exception("PomwigError")
00358                         << "HepMC Conversion problems in event." << std::endl;
00359 
00360         // do particle ID conversion Herwig->PDG, if requested
00361         if(doPDGConvert){
00362            for ( HepMC::GenEvent::particle_iterator part = event()->particles_begin(); part != event()->particles_end(); ++part) {
00363              if ((*part)->pdg_id() != HepPID::translateHerwigtoPDT((*part)->pdg_id()))
00364                (*part)->set_pdg_id(HepPID::translateHerwigtoPDT((*part)->pdg_id()));
00365            }
00366         }
00367 
00368         return true;
00369 }
00370 
00371 bool PomwigHadronizer::residualDecay()
00372 {
00373         return true;
00374 }
00375 
00376 }//namespace gen