[95a917c] | 1 | // -*- C++ -*-
|
---|
| 2 | //
|
---|
| 3 | // This file is part of HepMC
|
---|
| 4 | // Copyright (C) 2014-2019 The HepMC collaboration (see AUTHORS for details)
|
---|
| 5 | //
|
---|
| 6 | /**
|
---|
| 7 | * @file GenCrossSection.cc
|
---|
| 8 | * @brief Implementation of \b class GenCrossSection
|
---|
| 9 | *
|
---|
| 10 | */
|
---|
| 11 | #include "HepMC3/GenCrossSection.h"
|
---|
| 12 | #include "HepMC3/GenEvent.h"
|
---|
| 13 | #include <cstring> // memcmp
|
---|
| 14 | #include <cstdlib> // atoi
|
---|
| 15 | #include <sstream>
|
---|
| 16 | #include <iomanip>
|
---|
| 17 |
|
---|
| 18 | namespace HepMC3 {
|
---|
| 19 |
|
---|
| 20 |
|
---|
| 21 | int GenCrossSection::windx(std::string wName) const {
|
---|
| 22 | if ( !event() || !event()->run_info() ) return 0;
|
---|
| 23 | return event()->run_info()->weight_index(wName);
|
---|
| 24 | }
|
---|
| 25 |
|
---|
| 26 | void GenCrossSection::set_cross_section(const double& xs, const double& xs_err,const long& n_acc, const long& n_att ) {
|
---|
| 27 | double cross_section = xs;
|
---|
| 28 | double cross_section_error = xs_err;
|
---|
| 29 | accepted_events = n_acc;
|
---|
| 30 | attempted_events = n_att;
|
---|
| 31 | size_t N=1;
|
---|
| 32 | if ( event() ) N=std::max(event()->weights().size(),N);
|
---|
| 33 | cross_sections = std::vector<double>(N, cross_section);
|
---|
| 34 | cross_section_errors = std::vector<double>(N, cross_section_error);
|
---|
| 35 | }
|
---|
| 36 |
|
---|
| 37 |
|
---|
| 38 | bool GenCrossSection::from_string(const std::string &att) {
|
---|
| 39 | const char *cursor = att.data();
|
---|
| 40 | cross_sections.clear();
|
---|
| 41 | cross_section_errors.clear();
|
---|
| 42 |
|
---|
| 43 |
|
---|
| 44 | double cross_section = atof(cursor);
|
---|
| 45 | cross_sections.push_back(cross_section);
|
---|
| 46 |
|
---|
| 47 | if( !(cursor = strchr(cursor+1,' ')) ) return false;
|
---|
| 48 | double cross_section_error = atof(cursor);
|
---|
| 49 | cross_section_errors.push_back(cross_section_error);
|
---|
| 50 |
|
---|
| 51 | if( !(cursor = strchr(cursor+1,' ')) ) {accepted_events = -1; attempted_events = -1;}
|
---|
| 52 | else
|
---|
| 53 | {
|
---|
| 54 | accepted_events = atoi(cursor);
|
---|
| 55 | if( !(cursor = strchr(cursor+1,' ')) ) attempted_events = -1;
|
---|
| 56 | else attempted_events = atoi(cursor);
|
---|
| 57 | }
|
---|
| 58 | size_t N=1;
|
---|
| 59 | if ( event() ) N=std::max(event()->weights().size(),N);
|
---|
| 60 | const size_t max_n_cross_sections=1000;
|
---|
| 61 | while (cross_sections.size()<max_n_cross_sections) {
|
---|
| 62 | if( !(cursor = strchr(cursor+1,' ')) ) break;
|
---|
| 63 | cross_sections.push_back(atof(cursor));
|
---|
| 64 | if( !(cursor = strchr(cursor+1,' ')) ) break;
|
---|
| 65 | cross_section_errors.push_back(atof(cursor));
|
---|
| 66 | }
|
---|
| 67 | if (cross_sections.size()>=max_n_cross_sections)
|
---|
| 68 | HEPMC3_WARNING( "GenCrossSection::from_string: too many optional cross-sections N="<<cross_sections.size()<<" or ill-formed input:"<<att )
|
---|
| 69 | // if (cross_sections.size()!=N)
|
---|
| 70 | // So far it is not clear if there should be a warning or not
|
---|
| 71 | // Frank suggests no. HEPMC3_WARNING( "GenCrossSection::from_string: optional cross-sections are available not for all weights")
|
---|
| 72 | size_t oldsize=cross_sections.size();
|
---|
| 73 | for (size_t i=oldsize; i<N; i++) {cross_sections.push_back(cross_section); cross_section_errors.push_back(cross_section_error);}
|
---|
| 74 |
|
---|
| 75 | return true;
|
---|
| 76 | }
|
---|
| 77 |
|
---|
| 78 | bool GenCrossSection::to_string(std::string &att) const {
|
---|
| 79 | std::ostringstream os;
|
---|
| 80 |
|
---|
| 81 | os << std::setprecision(8) << std::scientific
|
---|
| 82 | << cross_sections.at(0) << " "
|
---|
| 83 | << cross_section_errors.at(0) << " "
|
---|
| 84 | << accepted_events << " "
|
---|
| 85 | << attempted_events;
|
---|
| 86 |
|
---|
| 87 | for (size_t i = 1; i < cross_sections.size(); ++i )
|
---|
| 88 | os << " " << cross_sections.at(i)
|
---|
| 89 | << " " << cross_section_errors.at(i);
|
---|
| 90 |
|
---|
| 91 | att = os.str();
|
---|
| 92 |
|
---|
| 93 | return true;
|
---|
| 94 | }
|
---|
| 95 |
|
---|
| 96 | bool GenCrossSection::operator==( const GenCrossSection& a ) const {
|
---|
| 97 | return ( memcmp( (void*)this, (void*) &a, sizeof(class GenCrossSection) ) == 0 );
|
---|
| 98 | }
|
---|
| 99 |
|
---|
| 100 | bool GenCrossSection::operator!=( const GenCrossSection& a ) const {
|
---|
| 101 | return !( a == *this );
|
---|
| 102 | }
|
---|
| 103 |
|
---|
| 104 | bool GenCrossSection::is_valid() const {
|
---|
| 105 | if( cross_sections.size() == 0 ) return false;
|
---|
| 106 | if( cross_section_errors.size() == 0 ) return false;
|
---|
| 107 | if( cross_section_errors.size()!=cross_sections.size() ) return false;
|
---|
| 108 | if( cross_sections.at(0) != 0 ) return true;
|
---|
| 109 | if( cross_section_errors.at(0) != 0 ) return true;
|
---|
| 110 | return false;
|
---|
| 111 | }
|
---|
| 112 |
|
---|
| 113 | } // namespace HepMC3
|
---|