CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

gen::JetMatchingAlpgen Class Reference

#include <JetMatchingAlpgen.h>

Inheritance diagram for gen::JetMatchingAlpgen:
gen::JetMatching

List of all members.

Public Member Functions

 JetMatchingAlpgen (const edm::ParameterSet &params)
 ~JetMatchingAlpgen ()

Private Member Functions

void beforeHadronisation (const lhef::LHEEvent *event)
std::set< std::string > capabilities () const
void init (const lhef::LHERunInfo *runInfo)
int match (const HepMC::GenEvent *partonLevel, const HepMC::GenEvent *finalState, bool showeredFinalState)

Private Attributes

bool applyMatching
bool eventInitialized
AlpgenHeader header
bool runInitialized

Detailed Description

Definition at line 12 of file JetMatchingAlpgen.h.


Constructor & Destructor Documentation

gen::JetMatchingAlpgen::JetMatchingAlpgen ( const edm::ParameterSet params)

Definition at line 34 of file JetMatchingAlpgen.cc.

References ahopts_, AHOPTS::etclus, edm::ParameterSet::getParameter(), AHOPTS::iexc, and AHOPTS::rclus.

                                                                  :
  JetMatching(params),
  applyMatching(params.getParameter<bool>("applyMatching")),
  runInitialized(false)
{
  ahopts_.etclus = params.getParameter<double>("etMin");
  ahopts_.rclus = params.getParameter<double>("drMin");
  ahopts_.iexc = params.getParameter<bool>("exclusive");
}
gen::JetMatchingAlpgen::~JetMatchingAlpgen ( )

Definition at line 45 of file JetMatchingAlpgen.cc.

{
}

Member Function Documentation

void gen::JetMatchingAlpgen::beforeHadronisation ( const lhef::LHEEvent event) [private, virtual]

Reimplemented from gen::JetMatching.

Definition at line 120 of file JetMatchingAlpgen.cc.

References eventInitialized, Exception, and runInitialized.

{
  // We can't continue if the run has not been initialized.
  if (!runInitialized)
    throw cms::Exception("Generator|PartonShowerVeto")
      << "Run not initialized in JetMatchingAlpgen"
      << std::endl;

  // We are called just after LHEInterface has filled in
  // the Fortran common block (and Pythia6 called UPEVNT).
  
  // Possibly not interesting for us.
  // (except perhaps for debugging?)
  //  pyupre_();
  //  dbpart_();
  eventInitialized = true;
}
std::set< std::string > gen::JetMatchingAlpgen::capabilities ( ) const [private, virtual]

Reimplemented from gen::JetMatching.

Definition at line 49 of file JetMatchingAlpgen.cc.

References query::result.

{
        std::set<std::string> result;
        result.insert("psFinalState");
        result.insert("hepevt");
        result.insert("pythia6");
        return result;
}
void gen::JetMatchingAlpgen::init ( const lhef::LHERunInfo runInfo) [private, virtual]

Reimplemented from gen::JetMatching.

Definition at line 60 of file JetMatchingAlpgen.cc.

References ahpars_, ahppara_, gen::alsetp_(), filterCSVwithJSON::copy, AHPPARA::ebeam, AlpgenHeader::ebeam, lhef::LHERunInfo::findHeader(), header, AHPPARA::ihrd, AlpgenHeader::ihrd, AlpgenHeader::MASS_MAX, AlpgenHeader::masses, AHPPARA::masses, AHPARS::nparam, AlpgenHeader::params, AlpgenHeader::parse(), AHPARS::parval, and runInitialized.

{

  // Read Alpgen run card stored in the LHERunInfo object.
  std::vector<std::string> headerLines = runInfo->findHeader("AlpgenUnwParFile");
  if (headerLines.empty())
    throw cms::Exception("Generator|PartonShowerVeto")
      << "In order to use Alpgen jet matching, "
      "the input file has to contain the corresponding "
      "Alpgen headers." << std::endl;
  
  // Parse the header using its bultin function.
  header.parse(headerLines.begin(), headerLines.end());

  // I don't want to print this right now.
//      std::cout << "Alpgen header" << std::endl;
//      std::cout << "========================" << std::endl;
//      std::cout << "\tihrd = " << header.ihrd << std::endl;
//      std::cout << "\tmc = " << header.masses[AlpgenHeader::mc]
//                << ", mb = " << header.masses[AlpgenHeader::mb]
//                << ", mt = " << header.masses[AlpgenHeader::mt]
//                << ", mw = " << header.masses[AlpgenHeader::mw]
//                << ", mz = " << header.masses[AlpgenHeader::mz]
//                << ", mh = " << header.masses[AlpgenHeader::mh]
//                << std::endl;
//      for(std::map<AlpgenHeader::Parameter, double>::const_iterator iter =
//              header.params.begin(); iter != header.params.end(); ++iter)
//              std::cout << "\t" << AlpgenHeader::parameterName(iter->first)
//                        << " = " << iter->second << std::endl;
//      std::cout << "\txsec = " << header.xsec
//                << " +-" << header.xsecErr << std::endl;
//      std::cout << "========================" << std::endl;
        
  // Here we pass a few header variables to common block and 
  // call Alpgen init routine to do the rest.
  // The variables passed first are the ones directly that
  // need to be set up "manually": IHRD and the masses.
  // (ebeam is set just to mimic the original code)
  // Afterwards, we pass the full spectrum of Alpgen
  // parameters directly into the AHPARS structure, to be
  // treated by AHSPAR which is called inside alsetp_().

  std::copy(header.masses, header.masses + AlpgenHeader::MASS_MAX, ahppara_.masses);
  ahppara_.ihrd = header.ihrd;
  ahppara_.ebeam = header.ebeam;
  
  for(std::map<AlpgenHeader::Parameter, double>::const_iterator iter =
        header.params.begin(); iter != header.params.end(); ++iter) {
    if (iter->first <= 0 || iter->first >= (int)AHPARS::nparam - 1)
      continue;
    ahpars_.parval[(int)iter->first - 1] = iter->second;
  }
  
  // Run the rest of the setup.
  alsetp_();
  
  // When we reach this point, the run is fully initialized.
  runInitialized = true;
}
int gen::JetMatchingAlpgen::match ( const HepMC::GenEvent *  partonLevel,
const HepMC::GenEvent *  finalState,
bool  showeredFinalState 
) [private, virtual]

Implements gen::JetMatching.

Definition at line 138 of file JetMatchingAlpgen.cc.

References gen::alveto_(), applyMatching, eventInitialized, Exception, and runInitialized.

{
  if (!showeredFinalState)
    throw cms::Exception("Generator|PartonShowerVeto")
      << "Alpgen matching expected parton shower "
      "final state." << std::endl;

  if (!runInitialized)
    throw cms::Exception("Generator|PartonShowerVeto")
      << "Run not initialized in JetMatchingAlpgen"
      << std::endl;

  if (!eventInitialized)
    throw cms::Exception("Generator|PartonShowerVeto")
      << "Event not initialized in JetMatchingAlpgen"
      << std::endl;

  // If matching not required (e.g., icckw = 0), don't run the
  // FORTRAN veto code.
  if(!applyMatching) return 0;
  
  // Call the Fortran veto code. 
  int veto = 0;
  alveto_(&veto);
  
  eventInitialized = false;

  // If event was vetoed, the variable veto will contain the number 1. 
  // In this case, we must return 1 - that will be used as the return value from UPVETO.
  // If event was accepted, the variable veto will contain the number 0.
  // In this case, we must return 0 - that will be used as the return value from UPVETO.
  return veto ? 1 : 0;
}

Member Data Documentation

Definition at line 27 of file JetMatchingAlpgen.h.

Referenced by match().

Definition at line 29 of file JetMatchingAlpgen.h.

Referenced by beforeHadronisation(), and match().

Definition at line 31 of file JetMatchingAlpgen.h.

Referenced by init().

Definition at line 28 of file JetMatchingAlpgen.h.

Referenced by beforeHadronisation(), init(), and match().