Fork me on GitHub

source: git/external/HepMC3/GenCrossSection.h@ 95a917c

Last change on this file since 95a917c was 95a917c, checked in by Pavel Demin <pavel.demin@…>, 3 years ago

add HepMC3 library

  • Property mode set to 100644
File size: 4.6 KB
RevLine 
[95a917c]1// -*- C++ -*-
2//
3// This file is part of HepMC
4// Copyright (C) 2014-2020 The HepMC collaboration (see AUTHORS for details)
5//
6#ifndef HEPMC3_CROSS_SECTION_H
7#define HEPMC3_CROSS_SECTION_H
8/**
9 * @file GenCrossSection.h
10 * @brief Definition of attribute \b class GenCrossSection
11 *
12 * @class HepMC3::GenCrossSection
13 * @brief Stores additional information about cross-section
14 *
15 * This is an example of event attribute used to store cross-section information
16 *
17 * This class is meant to be used to pass, on an event by event basis,
18 * the current best guess of the total cross section.
19 * It is expected that the final cross section will be stored elsewhere.
20 *
21 * - double cross_section; // cross section in pb.
22 * - double cross_section_error; // error associated with this cross section.
23 * - long accepted_events; ///< The number of events generated so far.
24 * - long attempted_events; ///< The number of events attempted so far.
25 *
26 * In addition, several cross sections and related info can be
27 * included in case of runs with mulltiple weights.
28 *
29 * The units of cross_section and cross_section_error are expected to be pb.
30 *
31 * @ingroup attributes
32 *
33 */
34#include <iostream>
35#include <algorithm>
36#include "HepMC3/Attribute.h"
37
38namespace HepMC3 {
39/** Deprecated */
40using namespace std;
41
42class GenCrossSection : public Attribute {
43
44//
45// Fields
46//
47private:
48
49 long accepted_events; ///< The number of events generated so far.
50 long attempted_events; ///< The number of events attempted so far.
51
52 std::vector<double> cross_sections; ///< Per-weight cross-section.
53 std::vector<double> cross_section_errors; ///< Per-weight errors.
54//
55// Functions
56//
57public:
58 /** @brief Implementation of Attribute::from_string */
59 bool from_string(const std::string &att) override;
60
61 /** @brief Implementation of Attribute::to_string */
62 bool to_string(std::string &att) const override;
63
64 /** @brief Set all fields */
65 void set_cross_section(const double& xs, const double& xs_err,const long& n_acc = -1, const long& n_att = -1);
66
67 /** @brief Set the number of accepted events
68 */
69 void set_accepted_events(const long& n_acc ) {
70 accepted_events=n_acc;
71 }
72
73 /** @brief Set the number of attempted events
74 */
75 void set_attempted_events(const long& n_att ) {
76 attempted_events=n_att;
77 }
78
79 /** @brief Get the number of accepted events
80 */
81 long get_accepted_events() const {
82 return accepted_events;
83 }
84
85 /** @brief Get the number of attempted events
86 */
87 long get_attempted_events() const {
88 return attempted_events;
89 }
90
91 /** @brief Set the cross section corresponding to the weight
92 named \a wName.
93 */
94 void set_xsec(const std::string& wName,const double& xs) {
95 set_xsec(windx(wName), xs);
96 }
97
98 /** @brief Set the cross section corresponding to the weight with
99 index \a indx.
100 */
101 void set_xsec(const int& indx, const double& xs) {
102 cross_sections[indx] = xs;
103 }
104
105 /** @brief Set the cross section error corresponding to the weight
106 named \a wName.
107 */
108 void set_xsec_err(const std::string& wName, const double& xs_err) {
109 set_xsec_err(windx(wName), xs_err);
110 }
111
112 /** @brief Set the cross section error corresponding to the weight
113 with index \a indx.
114 */
115 void set_xsec_err(const int& indx, const double& xs_err) {
116 cross_section_errors[indx] = xs_err;
117 }
118
119 /** @brief Get the cross section corresponding to the weight named
120 \a wName.
121 */
122 double xsec(const std::string& wName) const {
123 return xsec(windx(wName));
124 }
125
126 /** @brief Get the cross section corresponding to the weight with index
127 \a indx.
128 */
129 double xsec(const int& indx = 0) const {
130 return cross_sections[indx];
131 }
132
133 /** @brief Get the cross section error corresponding to the weight
134 named \a wName.
135 */
136 double xsec_err(const std::string& wName) const {
137 return xsec_err(windx(wName));
138 }
139
140 /** @brief Get the cross section error corresponding to the weight
141 with index \a indx.
142 */
143 double xsec_err(const int& indx = 0) const {
144 return cross_section_errors[indx];
145 }
146
147 bool operator==( const GenCrossSection& ) const; ///< Operator ==
148 bool operator!=( const GenCrossSection& ) const; ///< Operator !=
149 bool is_valid() const; ///< Verify that the instance contains non-zero information
150
151private:
152
153 /** @brief get the weight index given a weight name. */
154 int windx(std::string wName) const;
155
156};
157
158
159} // namespace HepMC3
160
161#endif
Note: See TracBrowser for help on using the repository browser.