ATLAS Offline Software
generate.h
Go to the documentation of this file.
1 /* emacs: this is -*- c++ -*- */
11 #ifndef GENERATE_H
12 #define GENERATE_H
13 
14 
15 #include <cmath>
16 #include <vector>
17 #include <iostream>
18 
19 #include "TH1D.h"
20 
21 #include "BasicRandom.h"
22 
23 
24 namespace generate {
25 
27 double get_mean( const std::vector<double>& m );
28 double get_rms( const std::vector<double>& m );
29 
30 
32 TH1D* integrate( TH1D* hin );
33 
35 void ScaleIntegral(TH1D* hin );
36 
38 TH1D* PDF(TH1D* hin);
39 
42 TH1D* smooth( TH1D* hin, bool sym=true );
43 
44 
47 public:
48  virtual ~generator_base() { }
49  virtual double generate() const = 0;
50 };
51 
52 
53 
56 
57  public:
58 
59  virtual double generate() const {
60  static BasicRandom r;
61  double x = std::tan((r.uniform()-0.5)*M_PI);
62  return x;
63  }
64 
65 };
66 
67 
68 
69 
72 class hist_generator : public generator_base {
73 
74 public:
75 
76  hist_generator(TH1D* h, bool _smooth=true );
77 
78  hist_generator(const hist_generator &) = delete;
80 
81  virtual ~hist_generator() { delete m_s; if ( m_random ) delete m_random; }
82 
84  double generate() const {
85  return invert(m_random->uniform());
86  }
87 
88  TH1D* histogram() { return m_s; }
89  TH1D* rawhistogram() { return m_raw; }
90  TH1D* smoothhistogram() { return m_smooth; }
91 
92 private:
93 
94  int getbin( double y ) const {
95  for ( unsigned i=0 ; i<m_y.size() ; i++ ) if ( y<=m_y[i] ) return i;
96  return m_y.size()-1;
97  }
98 
99  double invert( double y ) const {
100  int i = getbin(y);
101  if ( i==0 ) return 0;
102  else return (y-m_y[i-1])*m_dxdy[i-1]+m_x[i-1];
103  }
104 
105 
106 private:
107 
108  std::vector<double> m_y;
109  std::vector<double> m_x;
110 
111 
112  std::vector<double> m_dx;
113  std::vector<double> m_dy;
114 
115  std::vector<double> m_dxdy;
116 
117  TH1D* m_s; // accumated histogram
118  TH1D* m_raw; // raw histogram
119  TH1D* m_smooth; // smoothed raw histogram
120 
122 
123 };
124 
125 
126 
127 
131 
132 class experiment {
133 
134 public:
135 
136  experiment( TH1D* h, int Nexperiments=10, int fevents=0 );
137  experiment(const experiment& ) = delete;
138  experiment operator=(const experiment& ) = delete;
139 
140  double hmean() const { return m_hmean; }
141  double hrms() const { return m_hrms; }
142 
143  double mean() const { return m_global_mean; }
144  double mean_error() const { return m_global_mean_error; }
145 
146  double rms() const { return m_global_rms; }
147  double rms_error() const { return m_global_rms_error; }
148 
149  generator_base* gen() { return m_gen; }
150 
151  TH1D* rmshisto() { return m_THrms; }
152 
153 private:
154 
155  int m_N;
156 
157  std::vector<double> m_mean;
158  std::vector<double> m_rms;
159 
160  double m_hmean;
161  double m_hrms;
162 
166 
167  double m_global_rms;
169 
171 
173 };
174 
175 
176 
177 } // namespace generate
178 
179 #endif // GENERATE_H
generate::experiment::hmean
double hmean() const
Definition: generate.h:140
beamspotman.r
def r
Definition: beamspotman.py:676
generate::experiment::m_THrms
TH1D * m_THrms
Definition: generate.h:172
generate::experiment::operator=
experiment operator=(const experiment &)=delete
python.SystemOfUnits.m
int m
Definition: SystemOfUnits.py:91
generate::breit_generator
breit-wigner generator
Definition: generate.h:55
generate::hist_generator::generate
double generate() const
actually generate a random number from the distribution
Definition: generate.h:84
generate::experiment::experiment
experiment(const experiment &)=delete
TH1D
Definition: rootspy.cxx:342
generate::experiment::m_global_mean
double m_global_mean
actual histogram rms95
Definition: generate.h:164
M_PI
#define M_PI
Definition: ActiveFraction.h:11
generate::experiment
given a histogram, get a generator for the distribution in that histogram, then run some pseudo-exper...
Definition: generate.h:132
generate::hist_generator::hist_generator
hist_generator(const hist_generator &)=delete
generate::experiment::rmshisto
TH1D * rmshisto()
Definition: generate.h:151
generate::experiment::m_global_rms_error
double m_global_rms_error
Definition: generate.h:168
generate::hist_generator::histogram
TH1D * histogram()
Definition: generate.h:88
x
#define x
generate::hist_generator::m_dx
std::vector< double > m_dx
Definition: generate.h:112
generate::get_mean
double get_mean(const std::vector< double > &m)
mean and rms
Definition: generate.cxx:47
generate::generator_base::generate
virtual double generate() const =0
generate::experiment::rms_error
double rms_error() const
Definition: generate.h:147
generate::hist_generator::~hist_generator
virtual ~hist_generator()
Definition: generate.h:81
generate::generator_base::~generator_base
virtual ~generator_base()
Definition: generate.h:48
generate::experiment::gen
generator_base * gen()
Definition: generate.h:149
generate::hist_generator::m_s
TH1D * m_s
Definition: generate.h:117
generate::hist_generator::smoothhistogram
TH1D * smoothhistogram()
Definition: generate.h:90
lumiFormat.i
int i
Definition: lumiFormat.py:92
generate::experiment::m_hmean
double m_hmean
rms from each experiment
Definition: generate.h:160
generate::experiment::m_mean
std::vector< double > m_mean
number of experiments
Definition: generate.h:157
generate::hist_generator::getbin
int getbin(double y) const
Definition: generate.h:94
generate::PDF
TH1D * PDF(TH1D *hin)
generate a (normalised) pdf from a distribution
Definition: generate.cxx:112
drawFromPickle.tan
tan
Definition: drawFromPickle.py:36
generate::experiment::mean
double mean() const
Definition: generate.h:143
generate::experiment::m_N
int m_N
Definition: generate.h:155
generate::experiment::mean_error
double mean_error() const
Definition: generate.h:144
BasicRandom.h
generate::hist_generator::m_random
BasicRandom * m_random
Definition: generate.h:121
generate::hist_generator::operator=
hist_generator operator=(const hist_generator &)=delete
generate::hist_generator::invert
double invert(double y) const
Definition: generate.h:99
generate::smooth
TH1D * smooth(TH1D *hin, bool sym)
smooth a distribution about the maximum NB: maximum is not smoothed
Definition: generate.cxx:158
generate::experiment::rms
double rms() const
Definition: generate.h:146
generate::ScaleIntegral
void ScaleIntegral(TH1D *hin)
scale a distribution by the integral
Definition: generate.cxx:92
generate::experiment::m_global_rms
double m_global_rms
Definition: generate.h:167
generate::experiment::hrms
double hrms() const
Definition: generate.h:141
generate::hist_generator::m_raw
TH1D * m_raw
Definition: generate.h:118
generate
Definition: generate.cxx:28
generate::experiment::m_hrms
double m_hrms
actual histogram mean95
Definition: generate.h:161
generate::get_rms
double get_rms(const std::vector< double > &m)
Definition: generate.cxx:54
y
#define y
h
generate::experiment::experiment
experiment(TH1D *h, int Nexperiments=10, int fevents=0)
given a distribution, gnerate some number of pseudo experiments to estimate the uncertainties in the ...
Definition: generate.cxx:302
generate::experiment::m_gen
generator_base * m_gen
Definition: generate.h:170
BasicRandom
Definition: BasicRandom.h:32
generate::breit_generator::generate
virtual double generate() const
Definition: generate.h:59
generate::hist_generator::m_smooth
TH1D * m_smooth
Definition: generate.h:119
generate::hist_generator::m_dy
std::vector< double > m_dy
Definition: generate.h:113
generate::integrate
TH1D * integrate(TH1D *hin)
generate an integrated distribution
Definition: generate.cxx:68
generate::hist_generator::rawhistogram
TH1D * rawhistogram()
Definition: generate.h:89
generate::hist_generator
generate a distribution according to an input histogram (after smoothing)
Definition: generate.h:72
BasicRandom::uniform
double uniform()
Definition: BasicRandom.h:54
generate::hist_generator::m_x
std::vector< double > m_x
Definition: generate.h:109
generate::experiment::m_global_mean_error
double m_global_mean_error
Definition: generate.h:165
generate::experiment::m_rms
std::vector< double > m_rms
means from each experiment
Definition: generate.h:158
generate::hist_generator::m_y
std::vector< double > m_y
Definition: generate.h:108
generate::hist_generator::hist_generator
hist_generator(TH1D *h, bool _smooth=true)
Definition: generate.cxx:231
generate::generator_base
base class for a "generator" class
Definition: generate.h:46
generate::hist_generator::m_dxdy
std::vector< double > m_dxdy
Definition: generate.h:115