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::ProducerSourceFromFiles edm::ProducerSourceBase edm::FromFiles 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::ProducerSourceFromFiles
virtual bool noFiles () const override
 
 ProducerSourceFromFiles (ParameterSet const &pset, InputSourceDescription const &desc, bool realData)
 
virtual ~ProducerSourceFromFiles ()
 
- Public Member Functions inherited from edm::ProducerSourceBase
EventNumber_t event () const
 
unsigned int eventCreationDelay () const
 
EventID const & eventID () 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
 
 ProducerSourceBase (ParameterSet const &pset, InputSourceDescription const &desc, bool realData)
 
RunNumber_t run () const
 
unsigned int timeBetweenEvents () const
 
virtual ~ProducerSourceBase ()
 
- Public Member Functions inherited from edm::InputSource
std::shared_ptr< ActivityRegistryactReg () const
 Accessor for Activity Registry. More...
 
std::shared_ptr
< BranchIDListHelper
branchIDListHelper () const
 Accessor for branchIDListHelper. More...
 
void closeFile (FileBlock *, bool cleaningUpAfterException)
 close current file More...
 
void doBeginJob ()
 Called by framework at beginning of job. More...
 
void doBeginLumi (LuminosityBlockPrincipal &lbp, ProcessContext const *)
 Called by framework at beginning of lumi block. More...
 
void doBeginRun (RunPrincipal &rp, ProcessContext const *)
 Called by framework at beginning of run. More...
 
void doEndJob ()
 Called by framework at end of job. More...
 
void doEndLumi (LuminosityBlockPrincipal &lbp, bool cleaningUpAfterException, ProcessContext const *)
 Called by framework at end of lumi block. More...
 
void doEndRun (RunPrincipal &rp, bool cleaningUpAfterException, ProcessContext const *)
 Called by framework at end of run. More...
 
void doPostForkReacquireResources (std::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...
 
 InputSource (InputSource const &)=delete
 
void issueReports (EventID const &eventID)
 issue an event report More...
 
LuminosityBlockNumber_t luminosityBlock () const
 Accessor for current luminosity block number. More...
 
std::shared_ptr
< LuminosityBlockAuxiliary
luminosityBlockAuxiliary () const
 Called by the framework to merge or insert lumi in principal cache. More...
 
int maxEvents () const
 
int maxLuminosityBlocks () const
 
ModuleDescription const & moduleDescription () const
 Accessor for 'module' description. More...
 
ItemType nextItemType ()
 Advances the source to the next item. More...
 
InputSourceoperator= (InputSource const &)=delete
 
ProcessConfiguration const & processConfiguration () const
 Accessor for Process Configuration. More...
 
std::string const & processGUID () const
 Accessor for global process identifier. More...
 
ProcessHistoryRegistry const & processHistoryRegistry () const
 Const accessor for process history registry. More...
 
ProcessingMode processingMode () const
 RunsLumisAndEvents (default), RunsAndLumis, or Runs. More...
 
std::shared_ptr
< ProductRegistry const > 
productRegistry () const
 Accessor for product registry. More...
 
bool randomAccess () const
 
void readAndMergeLumi (LuminosityBlockPrincipal &lbp)
 Read next luminosity block (same as a prior lumi) More...
 
void readAndMergeRun (RunPrincipal &rp)
 Read next run (same as a prior run) More...
 
void readEvent (EventPrincipal &ep, StreamContext &)
 Read next event. More...
 
bool readEvent (EventPrincipal &ep, EventID const &, StreamContext &)
 Read a specific event. More...
 
std::unique_ptr< FileBlockreadFile ()
 Read next file. More...
 
void readLuminosityBlock (LuminosityBlockPrincipal &lumiPrincipal, HistoryAppender &historyAppender)
 Read next luminosity block (new lumi) More...
 
std::shared_ptr
< LuminosityBlockAuxiliary
readLuminosityBlockAuxiliary ()
 Read next luminosity block Auxilary. More...
 
void readRun (RunPrincipal &runPrincipal, HistoryAppender &historyAppender)
 Read next run (new run) More...
 
std::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...
 
SharedResourcesAcquirerresourceSharedWithDelayedReader () const
 Returns nullptr if no resource shared between the Source and a DelayedReader. More...
 
ProcessingController::ReverseState reverseState () const
 
void rewind ()
 Begin again at the first event. More...
 
RunNumber_t run () const
 Accessor for current run number. More...
 
std::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)
 
bool skipForForking ()
 
std::shared_ptr
< ThinnedAssociationsHelper
thinnedAssociationsHelper () const
 Accessor for thinnedAssociationsHelper. More...
 
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) override
 
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 void produce (edm::Event &event) override
 
bool readAlpgenEvent (lhef::HEPEUP &hepeup)
 Read an event and put it into the HEPEUP. More...
 
virtual bool setRunAndEventInfo (edm::EventID &, edm::TimeValue_t &, edm::EventAuxiliary::ExperimentType &) override
 
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::unique_ptr< lhef::HEPEUPhepeup_
 
std::unique_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, IsSynchronize
}
 
enum  ProcessingMode { Runs, RunsAndLumis, RunsLumisAndEvents }
 
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::ProducerSourceFromFiles
static void fillDescription (ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from edm::ProducerSourceBase
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::InputSource
void decreaseRemainingEventsBy (int iSkipped)
 
bool eventCached () const
 
std::shared_ptr
< LuminosityBlockPrincipal >
const 
luminosityBlockPrincipal () const
 
bool newLumi () const
 
bool newRun () const
 
ProcessHistoryRegistryprocessHistoryRegistryForUpdate () const
 
ProductRegistryproductRegistryUpdate () const
 
void reset () const
 
void resetEventCached ()
 
void resetLuminosityBlockAuxiliary (bool isNewLumi=true) const
 
void resetNewLumi ()
 
void resetNewRun ()
 
void resetRunAuxiliary (bool isNewRun=true) const
 
std::shared_ptr< RunPrincipal >
const 
runPrincipal () const
 
void setEventCached ()
 Called by the framework to merge or ached() const {return eventCached_;}. More...
 
void setLuminosityBlockAuxiliary (LuminosityBlockAuxiliary *lbp)
 
void setNewLumi ()
 
void setNewRun ()
 
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 94 of file AlpgenSource.cc.

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

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

Destructor.

Definition at line 147 of file AlpgenSource.cc.

148 {
149 }

Member Function Documentation

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

Reimplemented from edm::ProducerSourceBase.

Definition at line 161 of file AlpgenSource.cc.

References LHERunInfoProduct::Header::addLine(), indexGen::comments, AlpgenHeader::ebeam, lhef::HEPRUP::EBMUP, 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(), AlCaHLTBitMon_QueryRunRegistry::string, writeAlpgenParFile, writeAlpgenWgtFile, writeExtraHeader, lhef::HEPRUP::XERRUP, lhef::HEPRUP::XMAXUP, AlpgenHeader::xsec, AlpgenHeader::xsecErr, and lhef::HEPRUP::XSECUP.

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

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

281 {
282  std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
283  header.params.find(index);
284  if (pos == header.params.end())
285  throw cms::Exception("Generator|AlpgenInterface")
286  << "Requested Alpgen parameter \""
288  "not found in Alpgen parameter file." << std::endl;
289 
290  return T(pos->second);
291 }
std::map< Parameter, double > params
Definition: AlpgenHeader.h:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:78
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 294 of file AlpgenSource.cc.

References header, and AlpgenHeader::params.

296 {
297  std::map<AlpgenHeader::Parameter, double>::const_iterator pos =
298  header.params.find(index);
299  if (pos == header.params.end())
300  return defValue;
301  else
302  return T(pos->second);
303 }
std::map< Parameter, double > params
Definition: AlpgenHeader.h:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:78
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 305 of file AlpgenSource.cc.

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

Referenced by beginRun(), and readAlpgenEvent().

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

Implements edm::ProducerSourceBase.

Definition at line 433 of file AlpgenSource.cc.

References Exception, alpgen::fixEventHiggsTTbar(), alpgen::fixEventMultiBoson(), alpgen::fixEventSingleTop(), alpgen::fixEventTTbar(), alpgen::fixEventWZ(), header, hepeup_, AlpgenHeader::ihrd, AlpgenHeader::itopprc, AlpgenHeader::masses, AlpgenHeader::mb, and AlpgenHeader::params.

434 {
435 
436  // Here are the Alpgen routines for filling up the rest
437  // of the LHE Event Record. The .unw file has the information
438  // in a compressed way, e.g. it doesn't list the W boson -
439  // one has to reconstruct it from the e nu pair.
440  lhef::HEPEUP& hepeup = *hepeup_;
441 
442  switch(header.ihrd) {
443  case 1:
444  case 2:
445  case 3:
446  case 4:
447  case 10:
448  case 14:
449  case 15:
450  alpgen::fixEventWZ(hepeup);
451  break;
452  case 5:
454  break;
455  case 6:
456  alpgen::fixEventTTbar(hepeup);
457  break;
458  case 8:
460  break;
461  case 13:
464  break;
465  case 7:
466  case 9:
467  case 11:
468  case 12:
469  case 16:
470  // No fixes needed.
471  break;
472 
473  default:
474  throw cms::Exception("Generator|AlpgenInterface")
475  << "Unrecognized IHRD process code" << std::endl;
476  }
477 
478  // Create the LHEEventProduct and put it into the Event.
479  std::auto_ptr<LHEEventProduct> lheEvent(new LHEEventProduct(hepeup));
480  event.put(lheEvent);
481 
482  hepeup_.reset();
483 }
std::map< Parameter, double > params
Definition: AlpgenHeader.h:67
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:78
std::unique_ptr< lhef::HEPEUP > hepeup_
Definition: AlpgenSource.cc:80
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
void fixEventTTbar(lhef::HEPEUP &hepeup)
Fixes Event Record for ihrd = 6.
double masses[MASS_MAX]
Definition: AlpgenHeader.h:73
bool AlpgenSource::readAlpgenEvent ( lhef::HEPEUP hepeup)
private

Read an event and put it into the HEPEUP.

Definition at line 314 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(), lhef::HEPEUP::MOTHUP, nPart(), processID(), lhef::HEPEUP::PUP, lhef::HEPEUP::resize(), lhef::HEPEUP::SCALUP, mathSSE::sqrt(), and lhef::HEPEUP::XWGTUP.

Referenced by setRunAndEventInfo().

315 {
316  char buffer[512];
317  double dummy;
318  int nPart;
319  double sWgtRes;
320  double sQ;
321 
322  inputFile_->getline(buffer, sizeof buffer);
323  if (!inputFile_->good())
324  return false;
325 
326  std::istringstream ls(buffer);
327 
328  // Event number and process don't matter (or do they?)
329  ls >> dummy >> dummy;
330 
331  // Number of particles in the record, sample's average weight and Q scale
332  ls >> nPart >> sWgtRes >> sQ;
333 
334  if (ls.bad() || nPart < 1 || nPart > 1000)
335  return false;
336 
337  // Make room for the particles listed in the Alpgen file
338  hepeup.resize(nPart);
339 
340  // Scales, weights and process ID.
341  hepeup.SCALUP = sQ;
342  hepeup.XWGTUP = sWgtRes;
343  hepeup.IDPRUP = processID();
344 
345  // Incoming lines
346  for(int i = 0; i != 2; ++i) {
347  inputFile_->getline(buffer, sizeof buffer);
348  std::istringstream ls(buffer);
349  int flavour; ls >> flavour;
350  int colour1; ls >> colour1;
351  int colour2; ls >> colour2;
352  double zMomentum; ls >> zMomentum;
353 
354  if (inputFile_->bad())
355  return false;
356 
357  // Setting the HEPEUP of the incoming lines.
358  hepeup.IDUP[i] = flavour;
359  hepeup.ISTUP[i] = -1;
360  hepeup.MOTHUP[i].first = 0;
361  hepeup.MOTHUP[i].second = 0;
362  hepeup.PUP[i][0] = 0.;
363  hepeup.PUP[i][1] = 0.;
364  hepeup.PUP[i][2] = zMomentum;
365  hepeup.PUP[i][3] = std::fabs(zMomentum);
366  hepeup.PUP[i][4] = 0.;
367  if (colour1) colour1 += 500;
368  if (colour2) colour2 += 500;
369  hepeup.ICOLUP[i].first = colour1;
370  hepeup.ICOLUP[i].second = colour2;
371  }
372 
373  // Outgoing lines
374  for(int i = 2; i != nPart; ++i) {
375  inputFile_->getline(buffer, sizeof buffer);
376  std::istringstream ls(buffer);
377  int flavour; ls >> flavour;
378  int colour1; ls >> colour1;
379  int colour2; ls >> colour2;
380  double px, py, pz, mass;
381  ls >> px >> py >> pz >> mass;
382  double energy = std::sqrt(px*px + py*py + pz*pz + mass*mass);
383 
384  if (inputFile_->bad())
385  return false;
386 
387  // Setting the HEPEUP of the outgoing lines.
388  hepeup.IDUP[i] = flavour;
389  hepeup.ISTUP[i] = 1;
390  hepeup.MOTHUP[i].first = 1;
391  hepeup.MOTHUP[i].second = 2;
392  hepeup.PUP[i][0] = px;
393  hepeup.PUP[i][1] = py;
394  hepeup.PUP[i][2] = pz;
395  hepeup.PUP[i][3] = energy;
396  hepeup.PUP[i][4] = mass;
397  if (colour1) colour1 += 500;
398  if (colour2) colour2 += 500;
399  hepeup.ICOLUP[i].first = colour1;
400  hepeup.ICOLUP[i].second = colour2;
401  }
402 
403  return true;
404 }
int i
Definition: DBlmapReader.cc:9
std::unique_ptr< std::ifstream > inputFile_
Pointer to the input file.
Definition: AlpgenSource.cc:67
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:48
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< int > ISTUP
Definition: LesHouches.h:230
std::vector< int > IDUP
Definition: LesHouches.h:225
TString nPart(Int_t part, TString string, TString delimit=";", Bool_t removerest=true)
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
bool AlpgenSource::setRunAndEventInfo ( edm::EventID ,
edm::TimeValue_t ,
edm::EventAuxiliary::ExperimentType  
)
overrideprivatevirtual

Implements edm::ProducerSourceBase.

Definition at line 406 of file AlpgenSource.cc.

References Exception, hepeup_, inputFile_, nEvent_, readAlpgenEvent(), and skipEvents_.

407 {
408  // The LHE Event Record
409  hepeup_.reset(new lhef::HEPEUP);
410 
411  lhef::HEPEUP& hepeup = *hepeup_;
412 
413  // Read the .unw file until it is over.
414  for(;;) {
415  if (!readAlpgenEvent(hepeup)) {
416  if (inputFile_->eof())
417  return false;
418 
419  throw cms::Exception("Generator|AlpgenInterface")
420  << "AlpgenSource is not able read event no. "
421  << nEvent_ << std::endl;
422  }
423 
424  nEvent_++;
425  if (skipEvents_ > 0)
426  skipEvents_--;
427  else
428  break;
429  }
430  return true;
431 }
std::unique_ptr< std::ifstream > inputFile_
Pointer to the input file.
Definition: AlpgenSource.cc:67
unsigned long nEvent_
Number of events.
Definition: AlpgenSource.cc:73
unsigned long skipEvents_
Number of events to skip.
Definition: AlpgenSource.cc:70
std::unique_ptr< lhef::HEPEUP > hepeup_
Definition: AlpgenSource.cc:80
bool readAlpgenEvent(lhef::HEPEUP &hepeup)
Read an event and put it into the HEPEUP.
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 151 of file AlpgenSource.cc.

References header, AlpgenHeader::masses, and contentValuesCheck::ss.

Referenced by beginRun().

153 {
154  std::ostringstream ss;
155  ss << std::setw(9) << pdgId << " " << std::scientific
156  << std::setprecision(9) << header.masses[mass] << " # "
157  << comment << std::endl;
158  return ss.str();
159 }
AlpgenHeader header
Alpgen _unw.par file as an AlpgenHeader.
Definition: AlpgenSource.cc:78
double masses[MASS_MAX]
Definition: AlpgenHeader.h:73
#define comment(par)
Definition: vmac.h:161

Member Data Documentation

std::string AlpgenSource::extraHeaderFileName_
private

Name of the extra header file.

Definition at line 83 of file AlpgenSource.cc.

Referenced by beginRun().

std::string AlpgenSource::extraHeaderName_
private

Name given to the extra header.

Definition at line 86 of file AlpgenSource.cc.

Referenced by beginRun().

std::string AlpgenSource::fileName_
private

Name of the input file.

Definition at line 64 of file AlpgenSource.cc.

Referenced by AlpgenSource(), and beginRun().

AlpgenHeader AlpgenSource::header
private

Alpgen _unw.par file as an AlpgenHeader.

Definition at line 78 of file AlpgenSource.cc.

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

std::unique_ptr<lhef::HEPEUP> AlpgenSource::hepeup_
private

Definition at line 80 of file AlpgenSource.cc.

Referenced by produce(), and setRunAndEventInfo().

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

Pointer to the input file.

Definition at line 67 of file AlpgenSource.cc.

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

LHERunInfoProduct::Header AlpgenSource::lheAlpgenUnwParHeader
private

Alpgen _unw.par file as a LHE header.

Definition at line 76 of file AlpgenSource.cc.

Referenced by AlpgenSource(), and beginRun().

unsigned long AlpgenSource::nEvent_
private

Number of events.

Definition at line 73 of file AlpgenSource.cc.

Referenced by setRunAndEventInfo().

unsigned long AlpgenSource::skipEvents_
private

Number of events to skip.

Definition at line 70 of file AlpgenSource.cc.

Referenced by setRunAndEventInfo().

bool AlpgenSource::writeAlpgenParFile
private

Definition at line 90 of file AlpgenSource.cc.

Referenced by beginRun().

bool AlpgenSource::writeAlpgenWgtFile
private

configuration flags

Definition at line 89 of file AlpgenSource.cc.

Referenced by beginRun().

bool AlpgenSource::writeExtraHeader
private

Definition at line 91 of file AlpgenSource.cc.

Referenced by beginRun().