CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
AlpgenSource Class Reference
Inheritance diagram for AlpgenSource:
edm::ExternalInputSource edm::ConfigurableInputSource edm::InputSource edm::ProductRegistryHelper

Public Member Functions

 AlpgenSource (const edm::ParameterSet &params, const edm::InputSourceDescription &desc)
 Constructor. More...
 
virtual ~AlpgenSource ()
 Destructor. More...
 
- Public Member Functions inherited from edm::ExternalInputSource
InputFileCatalogcatalog ()
 
 ExternalInputSource (ParameterSet const &pset, InputSourceDescription const &desc, bool realData=true)
 
std::vector< std::string > const & fileNames () const
 
std::vector< std::string > const & logicalFileNames () const
 
virtual ~ExternalInputSource ()
 
- Public Member Functions inherited from edm::ConfigurableInputSource
 ConfigurableInputSource (ParameterSet const &pset, InputSourceDescription const &desc, bool realData=true)
 
EventNumber_t event () const
 
unsigned int eventCreationDelay () const
 
LuminosityBlockNumber_t luminosityBlock () const
 
unsigned int numberEventsInLumi () const
 
unsigned int numberEventsInRun () const
 
unsigned int numberEventsInThisLumi () const
 
unsigned int numberEventsInThisRun () const
 
TimeValue_t presentTime () const
 
RunNumber_t run () const
 
unsigned int timeBetweenEvents () const
 
virtual ~ConfigurableInputSource ()
 
- Public Member Functions inherited from edm::InputSource
boost::shared_ptr
< ActivityRegistry
actReg () const
 Accessor for Activity Registry. More...
 
void closeFile (boost::shared_ptr< FileBlock >, bool cleaningUpAfterException)
 close current file More...
 
void doBeginJob ()
 Called by framework at beginning of job. More...
 
void doBeginLumi (LuminosityBlockPrincipal &lbp)
 Called by framework at beginning of lumi block. More...
 
void doBeginRun (RunPrincipal &rp)
 Called by framework at beginning of run. More...
 
void doEndJob ()
 Called by framework at end of job. More...
 
void doEndLumi (LuminosityBlockPrincipal &lbp, bool cleaningUpAfterException)
 Called by framework at end of lumi block. More...
 
void doEndRun (RunPrincipal &rp, bool cleaningUpAfterException)
 Called by framework at end of run. More...
 
void doPostForkReacquireResources (boost::shared_ptr< multicore::MessageReceiverForSource >)
 
void doPreForkReleaseResources ()
 Called by the framework before forking the process. More...
 
ProcessingController::ForwardState forwardState () const
 
bool goToEvent (EventID const &eventID)
 
 InputSource (ParameterSet const &, InputSourceDescription const &)
 Constructor. More...
 
void issueReports (EventID const &eventID)
 issue an event report More...
 
LuminosityBlockNumber_t luminosityBlock () const
 Accessor for current luminosity block number. More...
 
boost::shared_ptr
< LuminosityBlockAuxiliary
luminosityBlockAuxiliary () const
 Called by the framework to merge or insert lumi in principal cache. More...
 
int markLumi ()
 Mark lumi as read. More...
 
int markRun ()
 Mark run as read. More...
 
int maxEvents () const
 
int maxLuminosityBlocks () const
 
ModuleDescription const & moduleDescription () const
 Accessor for 'module' description. More...
 
ItemType nextItemType ()
 
bool primary () const
 Accessor for primary input source flag. More...
 
ProcessConfiguration const & processConfiguration () const
 Accessor for Process Configuration. More...
 
std::string const & processGUID () const
 Accessor for global process identifier. More...
 
ProcessingMode processingMode () const
 RunsLumisAndEvents (default), RunsAndLumis, or Runs. More...
 
boost::shared_ptr
< ProductRegistry const > 
productRegistry () const
 Accessor for product registry. More...
 
bool randomAccess () const
 
void readAndCacheLumi (bool merge, HistoryAppender &historyAppender)
 Read next luminosity block. More...
 
void readAndCacheRun (bool merge, HistoryAppender &historyAppender)
 Read next run. More...
 
EventPrincipalreadEvent (boost::shared_ptr< LuminosityBlockPrincipal > lbCache)
 
EventPrincipalreadEvent (EventID const &)
 Read a specific event. More...
 
boost::shared_ptr< FileBlockreadFile ()
 Read next file. More...
 
boost::shared_ptr
< LuminosityBlockAuxiliary
readLuminosityBlockAuxiliary ()
 Read next luminosity block Auxilary. More...
 
boost::shared_ptr< RunAuxiliaryreadRunAuxiliary ()
 Read next run Auxiliary. More...
 
ProcessHistoryID const & reducedProcessHistoryID () const
 
void registerProducts ()
 Register any produced products. More...
 
int remainingEvents () const
 
int remainingLuminosityBlocks () const
 
void repeat ()
 Reset the remaining number of events/lumis to the maximum number. More...
 
ProcessingController::ReverseState reverseState () const
 
void rewind ()
 Begin again at the first event. More...
 
RunNumber_t run () const
 Accessor for current run number. More...
 
boost::shared_ptr< RunAuxiliaryrunAuxiliary () const
 Called by the framework to merge or insert run in principal cache. More...
 
void setLuminosityBlockNumber_t (LuminosityBlockNumber_t lb)
 Set the luminosity block ID. More...
 
void setRunNumber (RunNumber_t r)
 Set the run number. More...
 
void skipEvents (int offset)
 
Timestamp const & timestamp () const
 Accessor for the current time, as seen by the input source. More...
 
virtual ~InputSource ()
 Destructor. More...
 

Private Member Functions

virtual void beginRun (edm::Run &run)
 
template<typename T >
T getParameter (AlpgenHeader::Parameter index) const
 Function to get parameter by name from AlpgenHeader. More...
 
template<typename T >
T getParameter (AlpgenHeader::Parameter index, const T &defValue) const
 Function to get parameter by name from AlpgenHeader, w/ default. More...
 
unsigned int processID () const
 
virtual bool produce (edm::Event &event)
 
bool readAlpgenEvent (lhef::HEPEUP &hepeup)
 Read an event and put it into the HEPEUP. More...
 
std::string slhaMassLine (int pdgId, AlpgenHeader::Masses mass, const std::string &comment) const
 

Private Attributes

std::string extraHeaderFileName_
 Name of the extra header file. More...
 
std::string extraHeaderName_
 Name given to the extra header. More...
 
std::string fileName_
 Name of the input file. More...
 
AlpgenHeader header
 Alpgen _unw.par file as an AlpgenHeader. More...
 
std::auto_ptr< std::ifstream > inputFile_
 Pointer to the input file. More...
 
LHERunInfoProduct::Header lheAlpgenUnwParHeader
 Alpgen _unw.par file as a LHE header. More...
 
unsigned long nEvent_
 Number of events. More...
 
unsigned long skipEvents_
 Number of events to skip. More...
 
bool writeAlpgenParFile
 
bool writeAlpgenWgtFile
 configuration flags More...
 
bool writeExtraHeader
 

Additional Inherited Members

- Public Types inherited from edm::InputSource
enum  ItemType {
  IsInvalid, IsStop, IsFile, IsRun,
  IsLumi, IsEvent, IsRepeat
}
 
enum  ProcessingMode { Runs, RunsAndLumis, RunsLumisAndEvents }
 
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::ExternalInputSource
static void fillDescription (ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from edm::ConfigurableInputSource
static void fillDescription (ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from edm::InputSource
static const std::string & baseType ()
 
static void fillDescription (ParameterSetDescription &desc)
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::ExternalInputSource
void incrementFileIndex ()
 
- Protected Member Functions inherited from edm::ConfigurableInputSource
void reallyReadEvent ()
 
void setEventNumber (EventNumber_t e)
 
void setTime (TimeValue_t t)
 
- Protected Member Functions inherited from edm::InputSource
void decreaseRemainingEventsBy (int iSkipped)
 
EventPrincipaleventPrincipalCache ()
 
boost::shared_ptr
< LuminosityBlockPrincipal >
const 
luminosityBlockPrincipal () const
 
PrincipalCache const & principalCache () const
 
PrincipalCacheprincipalCache ()
 
ProductRegistryproductRegistryUpdate () const
 
void reset () const
 
void resetLuminosityBlockAuxiliary () const
 
void resetRunAuxiliary () const
 
boost::shared_ptr
< RunPrincipal > const 
runPrincipal () const
 
void setLuminosityBlockAuxiliary (LuminosityBlockAuxiliary *lbp)
 
void setRunAuxiliary (RunAuxiliary *rp)
 
void setTimestamp (Timestamp const &theTime)
 To set the current time, as seen by the input source. More...
 
ItemType state () const
 

Detailed Description

Definition at line 29 of file AlpgenSource.cc.

Constructor & Destructor Documentation

AlpgenSource::AlpgenSource ( const edm::ParameterSet params,
const edm::InputSourceDescription desc 
)

Constructor.

Definition at line 91 of file AlpgenSource.cc.

References LHERunInfoProduct::Header::addLine(), LHERunInfoProduct::Header::begin(), LHERunInfoProduct::Header::end(), edm::hlt::Exception, fileName_, edm::ExternalInputSource::fileNames(), header, lheAlpgenUnwParHeader, AlpgenHeader::parse(), and matplotRender::reader.

92  :
93  edm::ExternalInputSource(params, desc, false),
94  skipEvents_(params.getUntrackedParameter<unsigned int>("skipEvents", 0)),
95  nEvent_(0), lheAlpgenUnwParHeader("AlpgenUnwParFile"),
96  extraHeaderFileName_(params.getUntrackedParameter<std::string>("extraHeaderFileName","")),
97  extraHeaderName_(params.getUntrackedParameter<std::string>("extraHeaderName","")),
98  writeAlpgenWgtFile(params.getUntrackedParameter<bool>("writeAlpgenWgtFile", true)),
99  writeAlpgenParFile(params.getUntrackedParameter<bool>("writeAlpgenParFile", true)),
100  writeExtraHeader(params.getUntrackedParameter<bool>("writeExtraHeader", false))
101 {
102  std::vector<std::string> allFileNames = fileNames();
103 
104  // Only one filename
105  if (allFileNames.size() != 1)
106  throw cms::Exception("Generator|AlpgenInterface")
107  << "AlpgenSource needs exactly one file specified "
108  "for now." << std::endl;
109 
110  fileName_ = allFileNames[0];
111 
112  // Strip the "file:" prefix
113  if (fileName_.find("file:") != 0)
114  throw cms::Exception("Generator|AlpgenInterface") << "AlpgenSource only supports the file: scheme "
115  "for now." << std::endl;
116  fileName_.erase(0, 5);
117 
118  // Open the _unw.par file to store additional
119  // informations in the LHERunInfoProduct
120  std::ifstream reader((fileName_ + "_unw.par").c_str());
121  if (!reader.good())
122  throw cms::Exception("Generator|AlpgenInterface")
123  << "AlpgenSource was unable to open the file \""
124  << fileName_ << "_unw.par\"." << std::endl;
125 
126  // A full copy of the _unw.par file in the LHE header.
127  char buffer[256];
128  while(reader.getline(buffer, sizeof buffer))
129  lheAlpgenUnwParHeader.addLine(std::string(buffer) + "\n");
130 
131  // Parse that header to setup an Alpgen header,
132  // which will be used in the production itself.
135  throw cms::Exception("Generator|AlpgenInterface")
136  << "AlpgenSource was unable to parse the Alpgen "
137  << "unweighted parameter file." << std::endl;
138 
139  // Declare the products.
140  produces<LHERunInfoProduct, edm::InRun>();
141  produces<LHEEventProduct>();
142 }
T getUntrackedParameter(std::string const &, T const &) const
bool parse(const std::vector< std::string >::const_iterator &begin, const std::vector< std::string >::const_iterator &end)
Definition: AlpgenHeader.cc:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
const_iterator end() const
bool writeAlpgenParFile
Definition: AlpgenSource.cc:87
bool writeExtraHeader
Definition: AlpgenSource.cc:88
unsigned long nEvent_
Number of events.
Definition: AlpgenSource.cc:72
unsigned long skipEvents_
Number of events to skip.
Definition: AlpgenSource.cc:69
bool writeAlpgenWgtFile
configuration flags
Definition: AlpgenSource.cc:86
LHERunInfoProduct::Header lheAlpgenUnwParHeader
Alpgen _unw.par file as a LHE header.
Definition: AlpgenSource.cc:75
std::string extraHeaderFileName_
Name of the extra header file.
Definition: AlpgenSource.cc:80
void addLine(const std::string &line)
std::string fileName_
Name of the input file.
Definition: AlpgenSource.cc:63
std::string extraHeaderName_
Name given to the extra header.
Definition: AlpgenSource.cc:83
std::vector< std::string > const & fileNames() const
const_iterator begin() const
AlpgenSource::~AlpgenSource ( )
virtual

Destructor.

Definition at line 144 of file AlpgenSource.cc.

145 {
146 }

Member Function Documentation

void AlpgenSource::beginRun ( edm::Run run)
privatevirtual

Reimplemented from edm::ConfigurableInputSource.

Definition at line 158 of file AlpgenSource.cc.

References LHERunInfoProduct::Header::addLine(), indexGen::comments, AlpgenHeader::ebeam, lhef::HEPRUP::EBMUP, edm::hlt::Exception, extraHeaderFileName_, extraHeaderName_, fileName_, header, lhef::HEPRUP::IDBMUP, lhef::HEPRUP::IDWTUP, AlpgenHeader::ih2, inputFile_, lheAlpgenUnwParHeader, lhef::HEPRUP::LPRUP, AlpgenHeader::mb, AlpgenHeader::mc, AlpgenHeader::mh, AlpgenHeader::mt, AlpgenHeader::mw, AlpgenHeader::mz, lhef::HEPRUP::PDFGUP, lhef::HEPRUP::PDFSUP, processID(), edm::Run::put(), lhef::HEPRUP::resize(), slhaMassLine(), writeAlpgenParFile, writeAlpgenWgtFile, writeExtraHeader, lhef::HEPRUP::XERRUP, lhef::HEPRUP::XMAXUP, AlpgenHeader::xsec, AlpgenHeader::xsecErr, and lhef::HEPRUP::XSECUP.

159 {
160  // At this point, the lheUnwParHeader has the full contents of the _unw.par
161  // file. So we can get the HEPRUP information from the LHE header itself.
162  lhef::HEPRUP heprup;
163 
164  // Get basic run information.
165  // Beam identity.
166  heprup.IDBMUP.first = 2212;
167  switch(getParameter<int>(AlpgenHeader::ih2)) {
168  case 1:
169  heprup.IDBMUP.second = 2212;
170  break;
171  case -1:
172  heprup.IDBMUP.second = -2212;
173  break;
174  default:
175  throw cms::Exception("Generator|AlpgenInterface")
176  << "AlpgenSource was unable to understand the ih2 "
177  << "parameter." << std::endl;
178  }
179 
180  // Beam energy.
181  heprup.EBMUP.second = heprup.EBMUP.first =
182  getParameter<double>(AlpgenHeader::ebeam);
183 
184  // PDF info. Initially, Alpgen doesn't fill it.
185  heprup.PDFGUP.first = -1;
186  heprup.PDFGUP.second = -1;
187  heprup.PDFSUP.first = -1;
188  heprup.PDFSUP.second = -1;
189 
190  // Unweighted events.
191  heprup.IDWTUP = 3;
192 
193  // Only one process.
194  heprup.resize(1);
195 
196  // Cross section and error.
197  heprup.XSECUP[0] = header.xsec;
198  heprup.XERRUP[0] = header.xsecErr;
199 
200  // Maximum weight.
201  heprup.XMAXUP[0] = header.xsec;
202 
203  // Process code for Pythia.
204  heprup.LPRUP[0] = processID();
205 
206  // Comments on top.
208  comments.addLine("\n");
209  comments.addLine("\tExtracted by AlpgenInterface\n");
210 
211  // Add SLHA header containing particle masses from Alpgen.
212  // Pythia6Hadronisation will feed the masses to Pythia automatically.
213  LHERunInfoProduct::Header slha("slha");
214  slha.addLine("\n# SLHA header containing masses from Alpgen\n");
215  slha.addLine("Block MASS # Mass spectrum (kinematic masses)\n");
216  slha.addLine("# PDG Mass\n");
217  slha.addLine(slhaMassLine(4, AlpgenHeader::mc, "charm pole mass"));
218  slha.addLine(slhaMassLine(5, AlpgenHeader::mb, "bottom pole mass"));
219  slha.addLine(slhaMassLine(6, AlpgenHeader::mt, "top pole mass"));
220  slha.addLine(slhaMassLine(23, AlpgenHeader::mz, "Z mass"));
221  slha.addLine(slhaMassLine(24, AlpgenHeader::mw, "W mass"));
222  slha.addLine(slhaMassLine(25, AlpgenHeader::mh, "H mass"));
223 
224  char buffer[512];
225 
226  // We also add the information on weighted events.
227  LHERunInfoProduct::Header lheAlpgenWgtHeader("AlpgenWgtFile");
228  if (writeAlpgenWgtFile) {
229  std::ifstream wgtascii((fileName_+".wgt").c_str());
230  while(wgtascii.getline(buffer,512)) {
231  lheAlpgenWgtHeader.addLine(std::string(buffer) + "\n");
232  }
233  }
234 
235  LHERunInfoProduct::Header lheAlpgenParHeader("AlpgenParFile");
236  if (writeAlpgenParFile) {
237  std::ifstream parascii((fileName_+".par").c_str());
238  while(parascii.getline(buffer,512)) {
239  lheAlpgenParHeader.addLine(std::string(buffer) + "\n");
240  }
241  }
242 
243  // If requested by the user, we also add any specific header provided.
244  // Nota bene: the header is put in the LHERunInfoProduct AS IT WAS GIVEN.
245  // That means NO CROSS-CHECKS WHATSOEVER. Use with care.
246  LHERunInfoProduct::Header extraHeader(extraHeaderName_.c_str());
247  if(writeExtraHeader) {
248  std::ifstream extraascii(extraHeaderFileName_.c_str());
249  while(extraascii.getline(buffer,512)) {
250  extraHeader.addLine(std::string(buffer) + "\n");
251  }
252  }
253 
254  // Build the final Run info object. Backwards-compatible order.
255  std::auto_ptr<LHERunInfoProduct> runInfo(new LHERunInfoProduct(heprup));
256  runInfo->addHeader(comments);
257  runInfo->addHeader(lheAlpgenUnwParHeader);
258  if (writeAlpgenWgtFile)
259  runInfo->addHeader(lheAlpgenWgtHeader);
260  if (writeAlpgenParFile)
261  runInfo->addHeader(lheAlpgenParHeader);
262  runInfo->addHeader(slha);
263  if(writeExtraHeader)
264  runInfo->addHeader(extraHeader);
265  run.put(runInfo);
266 
267  // Open the .unw file in the heap, and set the global pointer to it.
268  inputFile_.reset(new std::ifstream((fileName_ + ".unw").c_str()));
269  if (!inputFile_->good())
270  throw cms::Exception("Generator|AlpgenInterface")
271  << "AlpgenSource was unable to open the file \""
272  << fileName_ << ".unw\"." << std::endl;
273 
274 }
void resize(int nrup)
Definition: LesHouches.h:52
std::pair< double, double > EBMUP
Definition: LesHouches.h:78
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
std::pair< int, int > IDBMUP
Definition: LesHouches.h:73
bool writeAlpgenParFile
Definition: AlpgenSource.cc:87
bool writeExtraHeader
Definition: AlpgenSource.cc:88
std::pair< int, int > PDFGUP
Definition: LesHouches.h:84
bool writeAlpgenWgtFile
configuration flags
Definition: AlpgenSource.cc:86
LHERunInfoProduct::Header lheAlpgenUnwParHeader
Alpgen _unw.par file as a LHE header.
Definition: AlpgenSource.cc:75
std::string extraHeaderFileName_
Name of the extra header file.
Definition: AlpgenSource.cc:80
void addLine(const std::string &line)
std::vector< double > XERRUP
Definition: LesHouches.h:114
std::string fileName_
Name of the input file.
Definition: AlpgenSource.cc:63
std::vector< double > XMAXUP
Definition: LesHouches.h:119
std::string slhaMassLine(int pdgId, AlpgenHeader::Masses mass, const std::string &comment) const
std::string extraHeaderName_
Name given to the extra header.
Definition: AlpgenSource.cc:83
std::pair< int, int > PDFSUP
Definition: LesHouches.h:90
std::auto_ptr< std::ifstream > inputFile_
Pointer to the input file.
Definition: AlpgenSource.cc:66
double xsecErr
Definition: AlpgenHeader.h:70
void put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Run.h:81
std::vector< double > XSECUP
Definition: LesHouches.h:108
std::vector< int > LPRUP
Definition: LesHouches.h:124
unsigned int processID() const
template<typename T >
T AlpgenSource::getParameter ( AlpgenHeader::Parameter  index) const
private

Function to get parameter by name from AlpgenHeader.

Definition at line 277 of file AlpgenSource.cc.

References header, AlpgenHeader::parameterName(), AlpgenHeader::params, and pos.

278 {
279  std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
280  header.params.find(index);
281  if (pos == header.params.end())
282  throw cms::Exception("Generator|AlpgenInterface")
283  << "Requested Alpgen parameter \""
285  "not found in Alpgen parameter file." << std::endl;
286 
287  return T(pos->second);
288 }
std::map< Parameter, double > params
Definition: AlpgenHeader.h:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
static std::string parameterName(Parameter index)
long double T
template<typename T >
T AlpgenSource::getParameter ( AlpgenHeader::Parameter  index,
const T defValue 
) const
private

Function to get parameter by name from AlpgenHeader, w/ default.

Definition at line 291 of file AlpgenSource.cc.

References header, AlpgenHeader::params, and pos.

293 {
294  std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
295  header.params.find(index);
296  if (pos == header.params.end())
297  return defValue;
298  else
299  return T(pos->second);
300 }
std::map< Parameter, double > params
Definition: AlpgenHeader.h:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
long double T
unsigned int AlpgenSource::processID ( ) const
private

The Alpgen process ID. This is defined as processID() = 100*X + 10*Y + Z, where = ihrd, Y = ihvy, Z = njets

Definition at line 302 of file AlpgenSource.cc.

References header, AlpgenHeader::ihrd, AlpgenHeader::ihvy, and AlpgenHeader::njets.

Referenced by beginRun(), and readAlpgenEvent().

303 {
304  // return 661; // The original, old thing.
305  // digits #XYZ: X = ihrd, Y = ihvy, Z = njets
306  return header.ihrd * 100 +
307  getParameter<unsigned int>(AlpgenHeader::ihvy, 0) * 10 +
308  getParameter<unsigned int>(AlpgenHeader::njets, 0);
309 }
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
unsigned int ihrd
Definition: AlpgenHeader.h:68
bool AlpgenSource::produce ( edm::Event event)
privatevirtual

Implements edm::ConfigurableInputSource.

Definition at line 403 of file AlpgenSource.cc.

References edm::hlt::Exception, alpgen::fixEventHiggsTTbar(), alpgen::fixEventMultiBoson(), alpgen::fixEventSingleTop(), alpgen::fixEventTTbar(), alpgen::fixEventWZ(), header, AlpgenHeader::ihrd, inputFile_, AlpgenHeader::itopprc, AlpgenHeader::masses, AlpgenHeader::mb, nEvent_, AlpgenHeader::params, readAlpgenEvent(), and skipEvents_.

404 {
405  // The LHE Event Record
406  lhef::HEPEUP hepeup;
407 
408  // Read the .unw file until it is over.
409  for(;;) {
410  if (!readAlpgenEvent(hepeup)) {
411  if (inputFile_->eof())
412  return false;
413 
414  throw cms::Exception("Generator|AlpgenInterface")
415  << "AlpgenSource is not able read event no. "
416  << nEvent_ << std::endl;
417  }
418 
419  nEvent_++;
420  if (skipEvents_ > 0)
421  skipEvents_--;
422  else
423  break;
424  }
425 
426  // Here are the Alpgen routines for filling up the rest
427  // of the LHE Event Record. The .unw file has the information
428  // in a compressed way, e.g. it doesn't list the W boson -
429  // one has to reconstruct it from the e nu pair.
430 
431  switch(header.ihrd) {
432  case 1:
433  case 2:
434  case 3:
435  case 4:
436  case 10:
437  case 14:
438  case 15:
439  alpgen::fixEventWZ(hepeup);
440  break;
441  case 5:
443  break;
444  case 6:
445  alpgen::fixEventTTbar(hepeup);
446  break;
447  case 8:
449  break;
450  case 13:
453  break;
454  case 7:
455  case 9:
456  case 11:
457  case 12:
458  case 16:
459  // No fixes needed.
460  break;
461 
462  default:
463  throw cms::Exception("Generator|AlpgenInterface")
464  << "Unrecognized IHRD process code" << std::endl;
465  }
466 
467  // Create the LHEEventProduct and put it into the Event.
468  std::auto_ptr<LHEEventProduct> lheEvent(new LHEEventProduct(hepeup));
469  event.put(lheEvent);
470 
471  return true;
472 }
std::map< Parameter, double > params
Definition: AlpgenHeader.h:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
unsigned long nEvent_
Number of events.
Definition: AlpgenSource.cc:72
unsigned long skipEvents_
Number of events to skip.
Definition: AlpgenSource.cc:69
void fixEventWZ(lhef::HEPEUP &hepeup)
Fixes Event Record for ihrd = 1,2,3,4,10,14,15.
void fixEventMultiBoson(lhef::HEPEUP &hepeup)
Fixes Event Record for ihrd = 5.
void fixEventHiggsTTbar(lhef::HEPEUP &hepeup)
Fixes Event Record for ihrd = 8.
void fixEventSingleTop(lhef::HEPEUP &hepeup, double mb, int itopprc)
Fixes Event Record for ihrd = 13.
unsigned int ihrd
Definition: AlpgenHeader.h:68
std::auto_ptr< std::ifstream > inputFile_
Pointer to the input file.
Definition: AlpgenSource.cc:66
void fixEventTTbar(lhef::HEPEUP &hepeup)
Fixes Event Record for ihrd = 6.
double masses[MASS_MAX]
Definition: AlpgenHeader.h:73
bool readAlpgenEvent(lhef::HEPEUP &hepeup)
Read an event and put it into the HEPEUP.
bool AlpgenSource::readAlpgenEvent ( lhef::HEPEUP hepeup)
private

Read an event and put it into the HEPEUP.

Definition at line 311 of file AlpgenSource.cc.

References relval_parameters_module::energy, reco::flavour(), i, lhef::HEPEUP::ICOLUP, lhef::HEPEUP::IDPRUP, lhef::HEPEUP::IDUP, inputFile_, lhef::HEPEUP::ISTUP, python.rootplot.utilities::ls(), scaleCards::mass, lhef::HEPEUP::MOTHUP, processID(), lhef::HEPEUP::PUP, lhef::HEPEUP::resize(), lhef::HEPEUP::SCALUP, mathSSE::sqrt(), and lhef::HEPEUP::XWGTUP.

Referenced by produce().

312 {
313  char buffer[512];
314  double dummy;
315  int nPart;
316  double sWgtRes;
317  double sQ;
318 
319  inputFile_->getline(buffer, sizeof buffer);
320  if (!inputFile_->good())
321  return false;
322 
323  std::istringstream ls(buffer);
324 
325  // Event number and process don't matter (or do they?)
326  ls >> dummy >> dummy;
327 
328  // Number of particles in the record, sample's average weight and Q scale
329  ls >> nPart >> sWgtRes >> sQ;
330 
331  if (ls.bad() || nPart < 1 || nPart > 1000)
332  return false;
333 
334  // Make room for the particles listed in the Alpgen file
335  hepeup.resize(nPart);
336 
337  // Scales, weights and process ID.
338  hepeup.SCALUP = sQ;
339  hepeup.XWGTUP = sWgtRes;
340  hepeup.IDPRUP = processID();
341 
342  // Incoming lines
343  for(int i = 0; i != 2; ++i) {
344  inputFile_->getline(buffer, sizeof buffer);
345  std::istringstream ls(buffer);
346  int flavour; ls >> flavour;
347  int colour1; ls >> colour1;
348  int colour2; ls >> colour2;
349  double zMomentum; ls >> zMomentum;
350 
351  if (inputFile_->bad())
352  return false;
353 
354  // Setting the HEPEUP of the incoming lines.
355  hepeup.IDUP[i] = flavour;
356  hepeup.ISTUP[i] = -1;
357  hepeup.MOTHUP[i].first = 0;
358  hepeup.MOTHUP[i].second = 0;
359  hepeup.PUP[i][0] = 0.;
360  hepeup.PUP[i][1] = 0.;
361  hepeup.PUP[i][2] = zMomentum;
362  hepeup.PUP[i][3] = std::fabs(zMomentum);
363  hepeup.PUP[i][4] = 0.;
364  if (colour1) colour1 += 500;
365  if (colour2) colour2 += 500;
366  hepeup.ICOLUP[i].first = colour1;
367  hepeup.ICOLUP[i].second = colour2;
368  }
369 
370  // Outgoing lines
371  for(int i = 2; i != nPart; ++i) {
372  inputFile_->getline(buffer, sizeof buffer);
373  std::istringstream ls(buffer);
374  int flavour; ls >> flavour;
375  int colour1; ls >> colour1;
376  int colour2; ls >> colour2;
377  double px, py, pz, mass;
378  ls >> px >> py >> pz >> mass;
379  double energy = std::sqrt(px*px + py*py + pz*pz + mass*mass);
380 
381  if (inputFile_->bad())
382  return false;
383 
384  // Setting the HEPEUP of the outgoing lines.
385  hepeup.IDUP[i] = flavour;
386  hepeup.ISTUP[i] = 1;
387  hepeup.MOTHUP[i].first = 1;
388  hepeup.MOTHUP[i].second = 2;
389  hepeup.PUP[i][0] = px;
390  hepeup.PUP[i][1] = py;
391  hepeup.PUP[i][2] = pz;
392  hepeup.PUP[i][3] = energy;
393  hepeup.PUP[i][4] = mass;
394  if (colour1) colour1 += 500;
395  if (colour2) colour2 += 500;
396  hepeup.ICOLUP[i].first = colour1;
397  hepeup.ICOLUP[i].second = colour2;
398  }
399 
400  return true;
401 }
int i
Definition: DBlmapReader.cc:9
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
void resize(int nup)
Definition: LesHouches.h:161
T sqrt(T t)
Definition: SSEVec.h:46
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< int > ISTUP
Definition: LesHouches.h:230
std::vector< int > IDUP
Definition: LesHouches.h:225
tuple mass
Definition: scaleCards.py:27
std::auto_ptr< std::ifstream > inputFile_
Pointer to the input file.
Definition: AlpgenSource.cc:66
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31
double XWGTUP
Definition: LesHouches.h:196
std::vector< std::pair< int, int > > ICOLUP
Definition: LesHouches.h:242
double SCALUP
Definition: LesHouches.h:210
unsigned int processID() const
std::string AlpgenSource::slhaMassLine ( int  pdgId,
AlpgenHeader::Masses  mass,
const std::string &  comment 
) const
private

Converts the AlpgenHeader::Masses to a std::string formatted as a slhaMassLine to facilitate passing them to Alpgen.

Definition at line 148 of file AlpgenSource.cc.

References header, scaleCards::mass, and AlpgenHeader::masses.

Referenced by beginRun().

150 {
151  std::ostringstream ss;
152  ss << std::setw(9) << pdgId << " " << std::scientific
153  << std::setprecision(9) << header.masses[mass] << " # "
154  << comment << std::endl;
155  return ss.str();
156 }
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:77
tuple mass
Definition: scaleCards.py:27
double masses[MASS_MAX]
Definition: AlpgenHeader.h:73
#define comment(par)
Definition: vmac.h:162

Member Data Documentation

std::string AlpgenSource::extraHeaderFileName_
private

Name of the extra header file.

Definition at line 80 of file AlpgenSource.cc.

Referenced by beginRun().

std::string AlpgenSource::extraHeaderName_
private

Name given to the extra header.

Definition at line 83 of file AlpgenSource.cc.

Referenced by beginRun().

std::string AlpgenSource::fileName_
private

Name of the input file.

Definition at line 63 of file AlpgenSource.cc.

Referenced by AlpgenSource(), and beginRun().

AlpgenHeader AlpgenSource::header
private

Alpgen _unw.par file as an AlpgenHeader.

Definition at line 77 of file AlpgenSource.cc.

Referenced by AlpgenSource(), beginRun(), getParameter(), processID(), produce(), and slhaMassLine().

std::auto_ptr<std::ifstream> AlpgenSource::inputFile_
private

Pointer to the input file.

Definition at line 66 of file AlpgenSource.cc.

Referenced by beginRun(), produce(), and readAlpgenEvent().

LHERunInfoProduct::Header AlpgenSource::lheAlpgenUnwParHeader
private

Alpgen _unw.par file as a LHE header.

Definition at line 75 of file AlpgenSource.cc.

Referenced by AlpgenSource(), and beginRun().

unsigned long AlpgenSource::nEvent_
private

Number of events.

Definition at line 72 of file AlpgenSource.cc.

Referenced by produce().

unsigned long AlpgenSource::skipEvents_
private

Number of events to skip.

Definition at line 69 of file AlpgenSource.cc.

Referenced by produce().

bool AlpgenSource::writeAlpgenParFile
private

Definition at line 87 of file AlpgenSource.cc.

Referenced by beginRun().

bool AlpgenSource::writeAlpgenWgtFile
private

configuration flags

Definition at line 86 of file AlpgenSource.cc.

Referenced by beginRun().

bool AlpgenSource::writeExtraHeader
private

Definition at line 88 of file AlpgenSource.cc.

Referenced by beginRun().