Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:01:12

0001 // -*- C++ -*-
0002 //
0003 // This file is part of HepMC
0004 // Copyright (C) 2014-2023 The HepMC collaboration (see AUTHORS for details)
0005 //
0006 #ifndef HEPMC3_CROSS_SECTION_H
0007 #define HEPMC3_CROSS_SECTION_H
0008 /**
0009  *  @file GenCrossSection.h
0010  *  @brief Definition of attribute \b class GenCrossSection
0011  *
0012  *  @class HepMC3::GenCrossSection
0013  *  @brief Stores additional information about cross-section
0014  *
0015  *  This is an example of event attribute used to store cross-section information
0016  *
0017  *  This class is meant to be used to pass, on an event by event basis,
0018  *  the current best guess of the total cross section.
0019  *  It is expected that the final cross section will be stored elsewhere.
0020  *
0021  *    - double cross_section;       // cross section in pb.
0022  *    - double cross_section_error; // error associated with this cross section.
0023  *    - long accepted_events;       ///< The number of events generated so far.
0024  *    - long attempted_events;      ///< The number of events attempted so far.
0025  *
0026  *  In addition, several cross sections and related info can be
0027  *  included in case of runs with mulltiple weights.
0028  *
0029  *  The units of cross_section and cross_section_error are expected to be pb.
0030  *
0031  *  @ingroup attributes
0032  *
0033  */
0034 #include <iostream>
0035 #include <algorithm>
0036 #include "HepMC3/Attribute.h"
0037 
0038 namespace HepMC3 {
0039 /** Deprecated */
0040 using namespace std;
0041 
0042 class GenCrossSection : public Attribute {
0043 
0044 //
0045 // Fields
0046 //
0047 private:
0048 
0049     long accepted_events;       ///< The number of events generated so far.
0050     long attempted_events;      ///< The number of events attempted so far.
0051 
0052     std::vector<double> cross_sections;       ///< Per-weight cross-section.
0053     std::vector<double> cross_section_errors; ///< Per-weight errors.
0054 //
0055 // Functions
0056 //
0057 public:
0058     /** @brief Implementation of Attribute::from_string */
0059     bool from_string(const std::string &att)  override;
0060 
0061     /** @brief Implementation of Attribute::to_string */
0062     bool to_string(std::string &att) const override;
0063     /// @name Deprecated functionality
0064     /// @{
0065     /** @brief Set all fields
0066         @deprecated Use set_cross_section(const std::vector<double>& xs, const std::vector<double>& xs_err instead
0067 
0068         @param xs Cross section
0069         @param xs_err Uncertainty on cross-section
0070         @param n_acc Number of accepted events
0071         @param n_att Number of attempted events
0072     */
0073     void set_cross_section(const double& xs, const double& xs_err,const long& n_acc = -1, const long& n_att = -1);
0074 
0075     /// @}
0076     /** @brief Set all fields */
0077     void set_cross_section(const std::vector<double>& xs, const std::vector<double>& xs_err,const long& n_acc = -1, const long& n_att = -1);
0078 
0079     /** @brief Get the cross-sections
0080      */
0081     const std::vector<double>& xsecs() const { return cross_sections; }
0082 
0083     /** @brief Get the cross-section errors
0084      */
0085     const std::vector<double>& xsec_errs() const { return cross_section_errors; }
0086 
0087 
0088     /** @brief Set the number of accepted events
0089      */
0090     void set_accepted_events(const long& n_acc ) {
0091         accepted_events=n_acc;
0092     }
0093 
0094     /** @brief Set the number of attempted events
0095      */
0096     void set_attempted_events(const long& n_att ) {
0097         attempted_events=n_att;
0098     }
0099 
0100     /** @brief Get the number of accepted events
0101      */
0102     long get_accepted_events() const {
0103         return accepted_events;
0104     }
0105 
0106     /** @brief Get the number of attempted events
0107      */
0108     long get_attempted_events() const {
0109         return  attempted_events;
0110     }
0111 
0112     /** @brief Set the cross section  corresponding to the weight
0113         named \a wName.
0114      */
0115     void set_xsec(const std::string& wName,const double& xs) {
0116         int pos = windx(wName);
0117         if ( pos < 0 ) throw std::runtime_error("GenCrossSection::set_xsec(const std::string&,const double&): no weight with given name in this run");
0118         set_xsec(pos, xs);
0119     }
0120 
0121     /** @brief Set the cross section corresponding to the weight with
0122         index \a indx.
0123      */
0124     void set_xsec(const unsigned long& index, const double& xs) {
0125         if ( index >= cross_sections.size() ) {throw std::runtime_error("GenCrossSection::set_xsec(const unsigned long&): index outside of range");}
0126         cross_sections[index] = xs;
0127     }
0128 
0129     /** @brief Set the cross section error corresponding to the weight
0130         named \a wName.
0131      */
0132     void set_xsec_err(const std::string& wName, const double& xs_err) {
0133         int pos = windx(wName);
0134         if ( pos < 0 ) throw std::runtime_error("GenCrossSection::set_xsec_err(const std::string&,const double&): no weight with given name in this run");
0135         set_xsec_err(pos, xs_err);
0136     }
0137 
0138     /** @brief Set the cross section error corresponding to the weight
0139         with index \a indx.
0140      */
0141     void set_xsec_err(const unsigned long& index, const double& xs_err) {
0142         if ( index >= cross_section_errors.size() ) {throw std::runtime_error("GenCrossSection::set_xsec_err(const unsigned long&): index outside of range");}
0143         cross_section_errors[index] = xs_err;
0144     }
0145 
0146     /** @brief Get the cross section corresponding to the weight named
0147         \a wName.
0148      */
0149     double xsec(const std::string& wName) const {
0150         int pos = windx(wName);
0151         if ( pos < 0 ) throw std::runtime_error("GenCrossSection::xsec(const std::string&): no weight with given name in this run");
0152         return xsec(pos);
0153     }
0154 
0155     /** @brief Get the cross section corresponding to the weight with index
0156         \a indx.
0157      */
0158     double xsec(const unsigned long& index = 0) const {
0159         if ( index < cross_sections.size() ) { return cross_sections.at(index); }
0160         else  { throw std::runtime_error("GenCrossSection::xsec(const unsigned long&): index outside of range");}
0161         return 0.0;
0162     }
0163 
0164     /** @brief Get the cross section error corresponding to the weight
0165         named \a wName.
0166      */
0167     double xsec_err(const std::string& wName) const {
0168         int pos = windx(wName);
0169         if ( pos < 0 ) throw std::runtime_error("GenCrossSection::xsec_err(const std::string&): no weight with given name in this run");
0170         return xsec_err(pos);
0171     }
0172 
0173     /** @brief Get the cross section error corresponding to the weight
0174         with index \a indx.
0175      */
0176     double xsec_err(const unsigned long& index = 0) const {
0177         if ( index < cross_section_errors.size() ) {return cross_section_errors.at(index);}
0178         else  { throw std::runtime_error("GenCrossSection::xsec_err(const unsigned long&): index outside of range");}
0179         return 0.0;
0180     }
0181 
0182     bool operator==( const GenCrossSection& ) const; ///< Operator ==
0183     bool operator!=( const GenCrossSection& ) const; ///< Operator !=
0184     bool is_valid()                           const; ///< Verify that the instance contains non-zero information
0185 
0186 private:
0187 
0188     /** @brief get the weight index given a weight name. */
0189     int windx(const std::string& wName) const;
0190 
0191 };
0192 
0193 } // namespace HepMC3
0194 
0195 #endif