CMS 3D CMS Logo

ThePEGInterface Class Reference

Id
ThePEGInterface.h,v 1.3 2008/06/27 16:59:51 stober Exp
More...

#include <GeneratorInterface/ThePEGInterface/interface/ThePEGInterface.h>

Inheritance diagram for ThePEGInterface:

lhef::ThePEGHadronisation ThePEGProducer ThePEGSource

List of all members.

Public Member Functions

 ThePEGInterface (const edm::ParameterSet &params)
virtual ~ThePEGInterface ()

Protected Member Functions

std::string dataFile (const edm::ParameterSet &pset, const std::string &paramName) const
std::string dataFile (const std::string &fileName) const
void initGenerator ()
void initRepository (const edm::ParameterSet &params) const

Static Protected Member Functions

static void clearAuxiliary (HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf)
static std::auto_ptr
< HepMC::GenEvent > 
convert (const ThePEG::EventPtr &event)
static void fillAuxiliary (HepMC::GenEvent *hepmc, HepMC::PdfInfo *pdf, const ThePEG::EventPtr &event)
static std::string resolveEnvVars (const std::string &s)

Protected Attributes

ThePEG::EGPtr eg_
std::auto_ptr
< HepMC::IO_BaseClass > 
iobc_

Private Member Functions

void readParameterSet (const edm::ParameterSet &base, const std::string &paramSet) const

Private Attributes

const std::string dataLocation_
const std::string generator_
const std::string run_
const int skipEvents_


Detailed Description

Id
ThePEGInterface.h,v 1.3 2008/06/27 16:59:51 stober Exp

Id
ThePEGInterface.cc,v 1.3 2008/06/27 16:59:51 stober Exp

Oliver Oberst <oberst@ekp.uni-karlsruhe.de> Fred-Markus Stober <stober@ekp.uni-karlsruhe.de>

Definition at line 23 of file ThePEGInterface.h.


Constructor & Destructor Documentation

ThePEGInterface::ThePEGInterface ( const edm::ParameterSet params  ) 

Definition at line 39 of file ThePEGInterface.cc.

References edm::ParameterSet::getUntrackedParameter(), iobc_, and out.

00039                                                             :
00040         dataLocation_(resolveEnvVars(pset.getParameter<string>("dataLocation"))),
00041         generator_(pset.getParameter<string>("generatorModule")),
00042         run_(pset.getParameter<string>("run")),
00043         skipEvents_(pset.getUntrackedParameter<int>("skipEvents", 0))
00044 {
00045         // Write events in hepmc ascii format for debugging purposes
00046         string eventLog = pset.getUntrackedParameter<string>("printEvents", "");
00047         if (!eventLog.empty()) {
00048                 iobc_.reset(new HepMC::IO_ExtendedAscii(eventLog.c_str(), ios::out));
00049                 edm::LogInfo("ThePEGSource") << "Event logging switched on (=> " << eventLog << ")";
00050         }
00051 }

ThePEGInterface::~ThePEGInterface (  )  [virtual]

Definition at line 53 of file ThePEGInterface.cc.

References eg_.

00054 {
00055         if (eg_)
00056                 eg_->finalize();
00057         edm::LogInfo("ThePEGInterface") << "Event generator finalized";
00058 }


Member Function Documentation

void ThePEGInterface::clearAuxiliary ( HepMC::GenEvent *  hepmc,
HepMC::PdfInfo *  pdf 
) [static, protected]

Definition at line 213 of file ThePEGInterface.cc.

Referenced by lhef::ThePEGHadronisation::doHadronisation(), ThePEGSource::produce(), and ThePEGProducer::produce().

00215 {
00216         if (hepmc) {
00217                 hepmc->set_event_scale(-1.0);
00218                 hepmc->set_alphaQCD(-1.0);
00219                 hepmc->set_alphaQED(-1.0);
00220         }
00221         if (pdf) {
00222                 pdf->set_id1(-100);
00223                 pdf->set_id2(-100);
00224                 pdf->set_x1(-1.0);
00225                 pdf->set_x2(-1.0);
00226                 pdf->set_scalePDF(-1.0);
00227                 pdf->set_pdf1(-1.0);
00228                 pdf->set_pdf2(-1.0);
00229         }
00230 }

auto_ptr< HepMC::GenEvent > ThePEGInterface::convert ( const ThePEG::EventPtr &  event  )  [static, protected]

Definition at line 206 of file ThePEGInterface.cc.

References ThePEG::HepMCConverter< HepMCEventT, Traits >::convert().

Referenced by lhef::ThePEGHadronisation::doHadronisation(), ThePEGProducer::produce(), and ThePEGSource::produce().

00208 {
00209         return std::auto_ptr<HepMC::GenEvent>(
00210                 ThePEG::HepMCConverter<HepMC::GenEvent>::convert(*event));
00211 }

std::string ThePEGInterface::dataFile ( const edm::ParameterSet pset,
const std::string &  paramName 
) const [protected]

std::string ThePEGInterface::dataFile ( const std::string &  fileName  )  const [protected]

Referenced by initRepository().

void ThePEGInterface::fillAuxiliary ( HepMC::GenEvent *  hepmc,
HepMC::PdfInfo *  pdf,
const ThePEG::EventPtr &  event 
) [static, protected]

Definition at line 232 of file ThePEGInterface.cc.

References scale, sqr(), funct::sqrt(), and v.

Referenced by lhef::ThePEGHadronisation::doHadronisation(), ThePEGSource::produce(), and ThePEGProducer::produce().

00235 {
00236         if (!event->primaryCollision())
00237                 return;
00238 
00239         ThePEG::tcEHPtr eh = ThePEG::dynamic_ptr_cast<ThePEG::tcEHPtr>(
00240                                 event->primaryCollision()->handler());
00241         double scale = eh->lastScale();
00242 
00243         if (hepmc) {
00244                 if (hepmc->event_scale() < 0 && scale > 0)
00245                         hepmc->set_event_scale(std::sqrt(scale) / ThePEG::GeV);
00246 
00247                 if (hepmc->alphaQCD() < 0)
00248                         hepmc->set_alphaQCD(eh->lastAlphaS());
00249                 if (hepmc->alphaQED() < 0)
00250                         hepmc->set_alphaQED(eh->lastAlphaEM());
00251         }
00252 
00253         if (pdf) {
00254                 const ThePEG::PPair &beams = eh->lastParticles();
00255                 const ThePEG::PPair &partons = eh->lastPartons();
00256                 ThePEG::tcPDFPtr pdf1 = eh->lastExtractor()->getPDF(
00257                                                 beams.first->dataPtr());
00258                 ThePEG::tcPDFPtr pdf2 = eh->lastExtractor()->getPDF(
00259                                                 beams.second->dataPtr());
00260                 double x1, x2;
00261 
00262                 if (pdf->id1() == -100) {
00263                         int id = partons.first->id();
00264                         pdf->set_id1(id == 21 ? 0 : id);
00265                 }
00266                 if (pdf->id2() == -100) {
00267                         int id = partons.second->id();
00268                         pdf->set_id2(id == 21 ? 0 : id);
00269                 }
00270 
00271                 if (pdf->scalePDF() < 0)
00272                         pdf->set_scalePDF(std::sqrt(scale) / ThePEG::GeV);
00273                 else
00274                         scale = ThePEG::sqr(pdf->scalePDF()) * ThePEG::GeV;
00275 
00276                 if (pdf->x1() < 0) {
00277                         x1 = eh->lastX1();
00278                         pdf->set_x1(x1);
00279                 } else
00280                         x1 = pdf->x1();
00281 
00282                 if (pdf->x2() < 0) {
00283                         x2 = eh->lastX2();
00284                         pdf->set_x2(x2);
00285                 } else
00286                         x2 = pdf->x2();
00287 
00288                 if (pdf1 && pdf->pdf1() < 0) {
00289                         double v = pdf1->xfx(beams.first->dataPtr(),
00290                                              partons.first->dataPtr(),
00291                                              scale, x1);
00292                         if (v > 0)
00293                                 pdf->set_pdf1(v);
00294                         else
00295                                 pdf->set_pdf2(-1.0);
00296                 }
00297                 if (pdf2 && pdf->pdf2() < 0) {
00298                         double v = pdf2->xfx(beams.first->dataPtr(),
00299                                              partons.first->dataPtr(),
00300                                              scale, x2);
00301                         if (v > 0)
00302                                 pdf->set_pdf2(v);
00303                         else
00304                                 pdf->set_pdf2(-1.0);
00305                 }
00306         }
00307 
00308 }

void ThePEGInterface::initGenerator (  )  [protected]

Definition at line 186 of file ThePEGInterface.cc.

References eg_, lat::endl(), Exception, generator_, i, run_, skipEvents_, and tmp.

Referenced by ThePEGProducer::beginRun(), ThePEGSource::beginRun(), and lhef::ThePEGHadronisation::newRunInfo().

00187 {
00188         // Get generator from the repository and initialize it
00189         ThePEG::BaseRepository::CheckObjectDirectory(generator_);
00190         ThePEG::EGPtr tmp = ThePEG::BaseRepository::GetObject<ThePEG::EGPtr>(generator_);
00191         if (tmp) {
00192                 eg_ = ThePEG::Repository::makeRun(tmp, run_);
00193                 eg_->initialize();
00194                 edm::LogInfo("ThePEGInterface") << "EventGenerator initialized";
00195         } else
00196                 throw cms::Exception("ThePEGInterface")
00197                         << "EventGenerator could not be initialized!" << endl;
00198 
00199         // Skip events
00200         for (int i = 0; i < skipEvents_; i++) {
00201                 eg_->shoot();
00202                 edm::LogInfo("ThePEGInterface") << "Event discarded";
00203         }
00204 }

void ThePEGInterface::initRepository ( const edm::ParameterSet params  )  const [protected]

Definition at line 130 of file ThePEGInterface.cc.

References dataFile(), lat::endl(), edm::ParameterSet::getParameter(), iter, SiStripLorentzAngle_cfi::read, and readParameterSet().

Referenced by lhef::ThePEGHadronisation::ThePEGHadronisation(), ThePEGProducer::ThePEGProducer(), and ThePEGSource::ThePEGSource().

00131 {
00132         /* Initialize the repository from
00133          * 1. the repository file (default values)
00134          * 2. the ThePEG config files
00135          * 3. the CMSSW config blocks
00136          */
00137         stringstream logstream;
00138 
00139         // Read the repository of serialized default values
00140         string repository = dataFile(pset, "repository");
00141         if (!repository.empty()) {
00142                 edm::LogInfo("ThePEGInterface") << "Loading repository (" << repository << ")";
00143                 ThePEG::Repository::load(repository);
00144         }
00145 
00146         if (!getenv("ThePEG_INSTALL_PATH")) {
00147                 vector<string> libdirlist = ThePEG::DynamicLoader::allPaths();
00148                 for(vector<string>::const_iterator libdir = libdirlist.begin();
00149                     libdir < libdirlist.end(); ++libdir) {
00150                         if (libdir->empty() || (*libdir)[0] != '/')
00151                                 continue;
00152                         if (boost::filesystem::exists(*libdir +
00153                                         "/../../share/ThePEG/PDFsets.index")) {
00154                                 setenv("ThePEG_INSTALL_PATH",
00155                                        libdir->c_str(), 0);
00156                                 break;
00157                         }
00158                 }
00159         }
00160 
00161         // Read ThePEG config files to read
00162         vector<string> configFiles = pset.getParameter<vector<string> >("configFiles");
00163 
00164         // Loop over the config files
00165         for(vector<string>::const_iterator iter = configFiles.begin();
00166             iter != configFiles.end(); ++iter) {
00167                 edm::LogInfo("ThePEGInterface") << "Reading config file (" << *iter << ")";
00168                 ThePEG::Repository::read(dataFile(*iter), logstream);
00169                 edm::LogInfo("ThePEGSource") << logstream.str();
00170         }
00171 
00172         // Read CMSSW config file parameter sets starting from "parameterSets"
00173         readParameterSet(pset, "parameterSets");
00174 
00175         // Print the directories where ThePEG looks for libs
00176         vector<string> libdirlist = ThePEG::DynamicLoader::allPaths();
00177         for(vector<string>::const_iterator libdir = libdirlist.begin();
00178             libdir < libdirlist.end(); ++libdir)
00179                 edm::LogInfo("ThePEGInterface") << "DynamicLoader path = " << *libdir << endl;
00180 
00181         // Output status information about the repository
00182         ThePEG::Repository::stats(logstream);
00183         edm::LogInfo("ThePEGInterface") << logstream.str();
00184 }

void ThePEGInterface::readParameterSet ( const edm::ParameterSet base,
const std::string &  paramSet 
) const [private]

Referenced by initRepository().

static std::string ThePEGInterface::resolveEnvVars ( const std::string &  s  )  [static, protected]


Member Data Documentation

const std::string ThePEGInterface::dataLocation_ [private]

Definition at line 50 of file ThePEGInterface.h.

ThePEG::EGPtr ThePEGInterface::eg_ [protected]

Definition at line 46 of file ThePEGInterface.h.

Referenced by ThePEGProducer::beginRun(), ThePEGSource::beginRun(), lhef::ThePEGHadronisation::doHadronisation(), ThePEGSource::endRun(), ThePEGProducer::endRun(), initGenerator(), lhef::ThePEGHadronisation::newRunInfo(), ThePEGSource::produce(), ThePEGProducer::produce(), and ~ThePEGInterface().

const std::string ThePEGInterface::generator_ [private]

Definition at line 51 of file ThePEGInterface.h.

Referenced by initGenerator().

std::auto_ptr<HepMC::IO_BaseClass> ThePEGInterface::iobc_ [protected]

Definition at line 47 of file ThePEGInterface.h.

Referenced by ThePEGSource::produce(), ThePEGProducer::produce(), and ThePEGInterface().

const std::string ThePEGInterface::run_ [private]

Definition at line 52 of file ThePEGInterface.h.

Referenced by initGenerator().

const int ThePEGInterface::skipEvents_ [private]

Definition at line 53 of file ThePEGInterface.h.

Referenced by initGenerator().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:33:24 2009 for CMSSW by  doxygen 1.5.4