CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/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/HerwigWrapper.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::readSettings( int )
00134 {
00135 
00136    clear();
00137 
00138    edm::LogVerbatim("") << "----------------------------------------------\n"
00139                         << "Initializing PomwigHadronizer\n"
00140                         << "----------------------------------------------\n";
00141 
00142    // Call hwudat to set up HERWIG block data
00143    hwudat();
00144  
00145    // Setting basic parameters ...
00146    hwproc.PBEAM1 = comEnergy/2.;
00147    hwproc.PBEAM2 = comEnergy/2.;
00148    // Choose beam particles for POMWIG depending on topology
00149    switch (diffTopology){
00150         case 0: //DPE
00151                 hwbmch.PART1[0]  = 'E';
00152                 hwbmch.PART1[1]  = '-';
00153                 hwbmch.PART2[0]  = 'E';
00154                 hwbmch.PART2[1]  = '-';
00155                 break;
00156         case 1: //SD survive PART1
00157                 hwbmch.PART1[0]  = 'E';
00158                 hwbmch.PART1[1]  = '-';
00159                 hwbmch.PART2[0]  = 'P';
00160                 hwbmch.PART2[1]  = ' ';
00161                 break;
00162         case 2: //SD survive PART2
00163                 hwbmch.PART1[0]  = 'P';
00164                 hwbmch.PART1[1]  = ' ';
00165                 hwbmch.PART2[0]  = 'E';
00166                 hwbmch.PART2[1]  = '-';
00167                 break;
00168         case 3: //Non diffractive
00169                 hwbmch.PART1[0]  = 'P';
00170                 hwbmch.PART1[1]  = ' ';
00171                 hwbmch.PART2[0]  = 'P';
00172                 hwbmch.PART2[1]  = ' ';
00173                 break;
00174         default:
00175                 throw edm::Exception(edm::errors::Configuration,"PomwigError")
00176           <<" Invalid Diff. Topology. Must be DPE(diffTopology = 0), SD particle 1 (diffTopology = 1), SD particle 2 (diffTopology = 2) and Non diffractive (diffTopology = 3)";
00177                 break;
00178    }
00179    for(int i=2;i<8;++i){
00180     hwbmch.PART1[i]  = ' ';
00181     hwbmch.PART2[i]  = ' ';}
00182 
00183    // initialize other common blocks ...
00184    call(hwigin);
00185 
00186    hwevnt.MAXER = 100000000;    // O(inf)
00187    hwpram.LWSUD = 0;            // don't write Sudakov form factors
00188    hwdspn.LWDEC = 0;            // don't write three/four body decays
00189                                 // (no fort.77 and fort.88 ...)a
00190 
00191    std::memset(hwprch.AUTPDF, ' ', sizeof hwprch.AUTPDF);
00192    for(unsigned int i = 0; i < 2; i++) {
00193         hwpram.MODPDF[i] = -111;
00194         std::memcpy(hwprch.AUTPDF[i], "HWLHAPDF", 8);
00195    }
00196 
00197    hwevnt.MAXPR = maxEventsToPrint;
00198    hwpram.IPRINT = herwigVerbosity;
00199 
00200    edm::LogVerbatim("") << "------------------------------------\n"
00201                         << "Reading HERWIG parameters\n"
00202                         << "------------------------------------\n";
00203     
00204    for(gen::ParameterCollector::const_iterator line = parameters.begin();
00205                                                line != parameters.end(); ++line) {
00206       edm::LogVerbatim("") << "   " << *line;
00207       if (!give(*line))
00208           throw edm::Exception(edm::errors::Configuration)
00209                 << "Herwig 6 did not accept the following: \""
00210                 << *line << "\"." << std::endl;
00211    }
00212 
00213    needClear = true;
00214 
00215    return true;
00216 
00217 }
00218 
00219 bool PomwigHadronizer::initializeForInternalPartons()
00220 {
00221  
00222    call(hwuinc);
00223 
00224    hwusta("PI0     ",1);
00225  
00226    if(!initializeDPDF()) return false;
00227 
00228    call(hweini);
00229 
00230    return true;
00231 }
00232 
00233 bool PomwigHadronizer::initializeDPDF()
00234 {
00235    // Initialize H1 pomeron/reggeon
00236 
00237    if(diffTopology == 3) return true;
00238  
00239    if((diffTopology != 0)&&(diffTopology != 1)&&(diffTopology != 2)) return false;
00240 
00241    int nstru = hwpram.NSTRU;
00242    int ifit = h1fit;
00243    if((nstru == 9)||(nstru == 10)){
00244          if((ifit <= 0)||(ifit >= 7)){
00245               throw edm::Exception(edm::errors::Configuration,"PomwigError")
00246               <<" Attempted to set non existant H1 1997 fit index. Has to be 1...6";
00247          }
00248          std::string aux((nstru == 9)?"Pomeron":"Reggeon"); 
00249          edm::LogVerbatim("") << "   H1 1997 pdf's: " << aux << "\n"
00250                               << "   IFIT = " << ifit;
00251          double xp = 0.1;
00252          double Q2 = 75.0;
00253          double xpq[13];
00254          qcd_1994(xp,Q2,xpq,ifit);
00255     } else if((nstru >= 12)&&(nstru <= 15)){
00256          bool isPom = (nstru == 12)||(nstru == 14);
00257          bool isFitA = (nstru == 12)||(nstru == 13);
00258          ifit = (isFitA)?1:2;
00259          std::string aux_0((isFitA)?"A":"B");
00260          std::string aux_1((isPom)?"Pomeron":"Reggeon"); 
00261          edm::LogVerbatim("") << "   H1 2006 Fit " << aux_0 << " " << aux_1 << "\n"
00262                               << "   IFIT = "<< ifit;
00263          double xp = 0.1;
00264          double Q2 = 75.0;
00265          double xpq[13];
00266          double f2[2];
00267          double fl[2];
00268          double c2[2];
00269          double cl[2];
00270          qcd_2006(xp,Q2,ifit,xpq,f2,fl,c2,cl);
00271     } else{
00272          throw edm::Exception(edm::errors::Configuration,"PomwigError")
00273                <<" 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)";
00274     }
00275  
00276     return true;
00277 }
00278 
00279 bool PomwigHadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00280 {
00281         for(std::vector<int>::const_iterator iter = pdgIds.begin();
00282             iter != pdgIds.end(); ++iter)
00283                 if (!markStable(*iter))
00284                         return false;
00285         return true;
00286 }
00287 
00288 void PomwigHadronizer::statistics()
00289 {
00290         double RNWGT = 1. / hwevnt.NWGTS;
00291         double AVWGT = hwevnt.WGTSUM * RNWGT;
00292 
00293         double xsec = 1.0e3 * AVWGT;
00294         xsec = survivalProbability*xsec;
00295 
00296         runInfo().setInternalXSec(xsec);
00297 }
00298 
00299 bool PomwigHadronizer::hadronize()
00300 {
00301    return false;
00302 }
00303 
00304 bool PomwigHadronizer::generatePartonsAndHadronize()
00305 {
00306         // hard process generation, parton shower, hadron formation
00307 
00308         InstanceWrapper wrapper(this);  // safe guard
00309 
00310         event().reset();
00311 
00312         // call herwig routines to create HEPEVT
00313 
00314         hwuine();       // initialize event
00315         
00316         if (callWithTimeout(10, hwepro)) { // process event and PS
00317           // We hung for more than 10 seconds
00318           int error = 199;
00319           hwwarn_("HWHGUP", &error);
00320         }
00321         
00322         hwbgen();       // parton cascades
00323 
00324         // call jimmy ... only if event is not killed yet by HERWIG
00325         if (useJimmy && doMPInteraction && !hwevnt.IERROR && call_hwmsct()) 
00326           return false;
00327         
00328         hwdhob();       // heavy quark decays
00329         hwcfor();       // cluster formation
00330         hwcdec();       // cluster decays
00331         
00332         // if event *not* killed by HERWIG, return true
00333         if (!hwevnt.IERROR) return true;
00334 
00335         hwufne();       // finalize event       
00336         return false;
00337 }
00338 
00339 void PomwigHadronizer::finalizeEvent()
00340 {
00341         lhef::LHEEvent::fixHepMCEventTimeOrdering(event().get());
00342 
00343         event()->set_signal_process_id(hwproc.IPROC);
00344 
00345         event()->weights().push_back(hwevnt.EVWGT);
00346 }
00347 
00348 bool PomwigHadronizer::decay()
00349 {
00350         // hadron decays
00351 
00352         InstanceWrapper wrapper(this);  // safe guard
00353 
00354         hwdhad();       // unstable particle decays
00355         hwdhvy();       // heavy flavour decays
00356         hwmevt();       // soft underlying event
00357 
00358         hwufne();       // finalize event
00359 
00360         if (hwevnt.IERROR)
00361                 return false;
00362 
00363         event().reset(new HepMC::GenEvent);
00364         if (!conv.fill_next_event(event().get()))
00365                 throw cms::Exception("PomwigError")
00366                         << "HepMC Conversion problems in event." << std::endl;
00367 
00368         // do particle ID conversion Herwig->PDG, if requested
00369         if(doPDGConvert){
00370            for ( HepMC::GenEvent::particle_iterator part = event()->particles_begin(); part != event()->particles_end(); ++part) {
00371              if ((*part)->pdg_id() != HepPID::translateHerwigtoPDT((*part)->pdg_id()))
00372                (*part)->set_pdg_id(HepPID::translateHerwigtoPDT((*part)->pdg_id()));
00373            }
00374         }
00375 
00376         return true;
00377 }
00378 
00379 bool PomwigHadronizer::residualDecay()
00380 {
00381         return true;
00382 }
00383 
00384 }//namespace gen