Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-03 04:18:17

0001 #ifndef SimDataFormats_GeneratorProducts_GenEventInfoProduct3_h
0002 #define SimDataFormats_GeneratorProducts_GenEventInfoProduct3_h
0003 
0004 #include <vector>
0005 #include <memory>
0006 
0007 #include "SimDataFormats/GeneratorProducts/interface/PdfInfo.h"
0008 
0009 namespace HepMC3 {
0010   class GenEvent;
0011 }  // namespace HepMC3
0012 
0013 /** \class GenEventInfoProduct3
0014  *
0015  */
0016 
0017 class GenEventInfoProduct3 {
0018 public:
0019   GenEventInfoProduct3();
0020   GenEventInfoProduct3(const HepMC3::GenEvent *evt);
0021   GenEventInfoProduct3(const GenEventInfoProduct3 &other);
0022   GenEventInfoProduct3(GenEventInfoProduct3 &&other) = default;
0023   ~GenEventInfoProduct3();
0024 
0025   GenEventInfoProduct3 &operator=(const GenEventInfoProduct3 &other);
0026   GenEventInfoProduct3 &operator=(GenEventInfoProduct3 &&other) = default;
0027 
0028   typedef gen::PdfInfo PDF;
0029 
0030   // getters
0031 
0032   std::vector<double> &weights() { return weights_; }
0033   const std::vector<double> &weights() const { return weights_; }
0034 
0035   double weight() const { return weights_.empty() ? 1.0 : weights_[0]; }
0036 
0037   double weightProduct() const;
0038 
0039   unsigned int signalProcessID() const { return signalProcessID_; }
0040 
0041   double qScale() const { return qScale_; }
0042   double alphaQCD() const { return alphaQCD_; }
0043   double alphaQED() const { return alphaQED_; }
0044 
0045   const PDF *pdf() const { return pdf_.get(); }
0046   bool hasPDF() const { return pdf() != nullptr; }
0047 
0048   const std::vector<double> &binningValues() const { return binningValues_; }
0049   bool hasBinningValues() const { return !binningValues_.empty(); }
0050 
0051   const std::vector<float> &DJRValues() const { return DJRValues_; }
0052   bool hasDJRValues() const { return !DJRValues_.empty(); }
0053 
0054   int nMEPartons() const { return nMEPartons_; }
0055 
0056   int nMEPartonsFiltered() const { return nMEPartonsFiltered_; }
0057 
0058   // setters
0059 
0060   void setWeights(const std::vector<double> &weights) { weights_ = weights; }
0061 
0062   void setSignalProcessID(unsigned int procID) { signalProcessID_ = procID; }
0063 
0064   void setScales(double q = -1., double qcd = -1., double qed = -1.) { qScale_ = q, alphaQCD_ = qcd, alphaQED_ = qed; }
0065 
0066   void setPDF(const PDF *pdf) { pdf_.reset(pdf ? new PDF(*pdf) : nullptr); }
0067 
0068   void setBinningValues(const std::vector<double> &values) { binningValues_ = values; }
0069 
0070   void setDJR(const std::vector<float> &values) { DJRValues_ = values; }
0071 
0072   void setNMEPartons(int n) { nMEPartons_ = n; }
0073 
0074   void setNMEPartonsFiltered(int n) { nMEPartonsFiltered_ = n; }
0075 
0076 private:
0077   // HepMC3::GenEvent provides a list of weights
0078   std::vector<double> weights_;
0079 
0080   // generator-dependent process ID
0081   unsigned int signalProcessID_;
0082 
0083   // information about scales
0084   double qScale_;
0085   double alphaQCD_, alphaQED_;
0086 
0087   // optional PDF info
0088   std::unique_ptr<PDF> pdf_;
0089 
0090   // If event was produced in bis, this contains
0091   // the values that were used to define which
0092   // bin the event belongs in
0093   // This replaces the genEventScale, which only
0094   // corresponds to Pythia pthat.  The RunInfo
0095   // will contain the information what physical
0096   // quantity these values actually belong to
0097   std::vector<double> binningValues_;
0098   std::vector<float> DJRValues_;
0099   int nMEPartons_;
0100   int nMEPartonsFiltered_;
0101 };
0102 
0103 #endif  // SimDataFormats_GeneratorProducts_GenEventInfoProduct3_h