CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/GeneratorInterface/PartonShowerVeto/src/JetMatchingAlpgen.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <cstring>
00003 #include <vector>
00004 #include <memory>
00005 #include <string>
00006 
00007 #include <HepMC/GenEvent.h>
00008 
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 
00012 #include "GeneratorInterface/LHEInterface/interface/LHERunInfo.h"
00013 #include "GeneratorInterface/LHEInterface/interface/LHEEvent.h"
00014 #include "GeneratorInterface/PartonShowerVeto/interface/JetMatchingAlpgen.h"
00015 
00016 namespace gen {
00017 
00018 extern "C" {
00019     // Put here all the functions for interfacing Fortran code.
00020     
00021 //      extern void alinit_();
00022         extern void alsetp_();
00023 //      extern void alevnt_();
00024         extern void alveto_(int *ipveto);
00025 
00026       extern void dbpart_();
00027       extern void pyupre_();
00028 
00029         void alshcd_(char csho[3]);
00030         void alshen_();
00031 } // extern "C"
00032 
00033 // Constructor
00034 JetMatchingAlpgen::JetMatchingAlpgen(const edm::ParameterSet &params) :
00035   JetMatching(params),
00036   applyMatching(params.getParameter<bool>("applyMatching")),
00037   runInitialized(false)
00038 {
00039   ahopts_.etclus = params.getParameter<double>("etMin");
00040   ahopts_.rclus = params.getParameter<double>("drMin");
00041   ahopts_.iexc = params.getParameter<bool>("exclusive");
00042 }
00043 
00044 // Destructor   
00045 JetMatchingAlpgen::~JetMatchingAlpgen()
00046 {
00047 }
00048 
00049 std::set<std::string> JetMatchingAlpgen::capabilities() const
00050 {
00051         std::set<std::string> result;
00052         result.insert("psFinalState");
00053         result.insert("hepevt");
00054         result.insert("pythia6");
00055         return result;
00056 }
00057 
00058 // Implements the Alpgen MLM method - use alpg_match.F
00059 
00060 void JetMatchingAlpgen::init(const lhef::LHERunInfo* runInfo)
00061 {
00062 
00063   // Read Alpgen run card stored in the LHERunInfo object.
00064   std::vector<std::string> headerLines = runInfo->findHeader("AlpgenUnwParFile");
00065   if (headerLines.empty())
00066     throw cms::Exception("Generator|PartonShowerVeto")
00067       << "In order to use Alpgen jet matching, "
00068       "the input file has to contain the corresponding "
00069       "Alpgen headers." << std::endl;
00070   
00071   // Parse the header using its bultin function.
00072   header.parse(headerLines.begin(), headerLines.end());
00073 
00074   // I don't want to print this right now.
00075 //      std::cout << "Alpgen header" << std::endl;
00076 //      std::cout << "========================" << std::endl;
00077 //      std::cout << "\tihrd = " << header.ihrd << std::endl;
00078 //      std::cout << "\tmc = " << header.masses[AlpgenHeader::mc]
00079 //                << ", mb = " << header.masses[AlpgenHeader::mb]
00080 //                << ", mt = " << header.masses[AlpgenHeader::mt]
00081 //                << ", mw = " << header.masses[AlpgenHeader::mw]
00082 //                << ", mz = " << header.masses[AlpgenHeader::mz]
00083 //                << ", mh = " << header.masses[AlpgenHeader::mh]
00084 //                << std::endl;
00085 //      for(std::map<AlpgenHeader::Parameter, double>::const_iterator iter =
00086 //              header.params.begin(); iter != header.params.end(); ++iter)
00087 //              std::cout << "\t" << AlpgenHeader::parameterName(iter->first)
00088 //                        << " = " << iter->second << std::endl;
00089 //      std::cout << "\txsec = " << header.xsec
00090 //                << " +-" << header.xsecErr << std::endl;
00091 //      std::cout << "========================" << std::endl;
00092         
00093   // Here we pass a few header variables to common block and 
00094   // call Alpgen init routine to do the rest.
00095   // The variables passed first are the ones directly that
00096   // need to be set up "manually": IHRD and the masses.
00097   // (ebeam is set just to mimic the original code)
00098   // Afterwards, we pass the full spectrum of Alpgen
00099   // parameters directly into the AHPARS structure, to be
00100   // treated by AHSPAR which is called inside alsetp_().
00101 
00102   std::copy(header.masses, header.masses + AlpgenHeader::MASS_MAX, ahppara_.masses);
00103   ahppara_.ihrd = header.ihrd;
00104   ahppara_.ebeam = header.ebeam;
00105   
00106   for(std::map<AlpgenHeader::Parameter, double>::const_iterator iter =
00107         header.params.begin(); iter != header.params.end(); ++iter) {
00108     if (iter->first <= 0 || iter->first >= (int)AHPARS::nparam - 1)
00109       continue;
00110     ahpars_.parval[(int)iter->first - 1] = iter->second;
00111   }
00112   
00113   // Run the rest of the setup.
00114   alsetp_();
00115   
00116   // When we reach this point, the run is fully initialized.
00117   runInitialized = true;
00118 }
00119 
00120 void JetMatchingAlpgen::beforeHadronisation(const lhef::LHEEvent* event)
00121 {
00122   // We can't continue if the run has not been initialized.
00123   if (!runInitialized)
00124     throw cms::Exception("Generator|PartonShowerVeto")
00125       << "Run not initialized in JetMatchingAlpgen"
00126       << std::endl;
00127 
00128   // We are called just after LHEInterface has filled in
00129   // the Fortran common block (and Pythia6 called UPEVNT).
00130   
00131   // Possibly not interesting for us.
00132   // (except perhaps for debugging?)
00133   //  pyupre_();
00134   //  dbpart_();
00135   eventInitialized = true;
00136 }
00137 
00138 int JetMatchingAlpgen::match(const HepMC::GenEvent *partonLevel,
00139                              const HepMC::GenEvent *finalState,
00140                              bool showeredFinalState)
00141 {
00142   if (!showeredFinalState)
00143     throw cms::Exception("Generator|PartonShowerVeto")
00144       << "Alpgen matching expected parton shower "
00145       "final state." << std::endl;
00146 
00147   if (!runInitialized)
00148     throw cms::Exception("Generator|PartonShowerVeto")
00149       << "Run not initialized in JetMatchingAlpgen"
00150       << std::endl;
00151 
00152   if (!eventInitialized)
00153     throw cms::Exception("Generator|PartonShowerVeto")
00154       << "Event not initialized in JetMatchingAlpgen"
00155       << std::endl;
00156 
00157   // If matching not required (e.g., icckw = 0), don't run the
00158   // FORTRAN veto code.
00159   if(!applyMatching) return 0;
00160   
00161   // Call the Fortran veto code. 
00162   int veto = 0;
00163   alveto_(&veto);
00164   
00165   eventInitialized = false;
00166 
00167   // If event was vetoed, the variable veto will contain the number 1. 
00168   // In this case, we must return 1 - that will be used as the return value from UPVETO.
00169   // If event was accepted, the variable veto will contain the number 0.
00170   // In this case, we must return 0 - that will be used as the return value from UPVETO.
00171   return veto ? 1 : 0;
00172 }
00173 
00174 void alshcd_(char csho[3])
00175 {
00176         std::strncpy(csho, "PYT", 3);   // or "HER"
00177 }
00178 
00179 void alshen_()
00180 {
00181 }
00182 
00183 } // end namespace gen