ATLAS Offline Software
Loading...
Searching...
No Matches
MuonCalib::VariableBinwidthHistogram Class Reference

A histogram where every bin has the same number of entries. More...

#include <VariableBinwidthHistogram.h>

Collaboration diagram for MuonCalib::VariableBinwidthHistogram:

Public Member Functions

 VariableBinwidthHistogram ()
 Default constructor.
 ~VariableBinwidthHistogram ()
 destructor
bool Initialize (TH1F *hist, double binc_r, double max_bin_width, double min_x=-9e9, double max_x=9e9)
 Initialize with new input histogram Returns on error false.
bool Smooth (double width)
 Removes steps that origin in a binning effekt.
const VariableBinwidthHistogramBinGetBin (unsigned int bin) const
 Get a bin.
const VariableBinwidthHistogramBinGetSortedBin (unsigned int bin)
 Get a bin sorted by content.
unsigned int GetNumberOfBins () const
 Get the number of bins.
TGraph * DenistyGraph () const
 create density graph - density vs bin center
TGraph * BinWidthGraph () const
 Plot binwidth graph - binwidth versus bin center.
TGraph * BinContentGraph () const
 Plot bin content graph - bin content vs bin center.
TGraph * DiffDensityGraph () const
 Plot graph with differential density.
TGraph * DiffBinwidthGraph () const
 Plot graph with differential binwidth.

Private Member Functions

double sign (double val) const

Private Attributes

std::vector< VariableBinwidthHistogramBin * > m_bins
 vector containing histogram bins
std::vector< VBHBinPtrSrtm_sort_bins
 bins sorted by content
unsigned int m_binc
 number of entries per bin
double m_max_bin_width
 maximum bin width
bool m_error
 error flag
bool m_sorted

Detailed Description

A histogram where every bin has the same number of entries.

The density is represented by the bin width.

Author
Felix.nosp@m..Rau.nosp@m.scher.nosp@m.@Phy.nosp@m.sik.U.nosp@m.ni-M.nosp@m.uench.nosp@m.en.D.nosp@m.e
Date
February 2006

Definition at line 40 of file VariableBinwidthHistogram.h.

Constructor & Destructor Documentation

◆ VariableBinwidthHistogram()

MuonCalib::VariableBinwidthHistogram::VariableBinwidthHistogram ( )
inline

◆ ~VariableBinwidthHistogram()

MuonCalib::VariableBinwidthHistogram::~VariableBinwidthHistogram ( )
inline

destructor

Definition at line 47 of file VariableBinwidthHistogram.h.

47 {
48 for (unsigned int i = 0; i < m_bins.size(); i++) { delete m_bins[i]; }
49 }
std::vector< VariableBinwidthHistogramBin * > m_bins
vector containing histogram bins

Member Function Documentation

◆ BinContentGraph()

TGraph * MuonCalib::VariableBinwidthHistogram::BinContentGraph ( ) const

Plot bin content graph - bin content vs bin center.

Definition at line 141 of file VariableBinwidthHistogram.cxx.

141 {
142 Double_t *x = new Double_t[m_bins.size()], *y = new Double_t[m_bins.size()];
143 for (unsigned int i = 0; i < m_bins.size(); i++) {
144 x[i] = m_bins[i]->Center();
145 y[i] = m_bins[i]->Entries();
146 }
147 TGraph *gr = new TGraph(m_bins.size(), x, y);
148 return gr;
149 }
#define gr
#define y
#define x

◆ BinWidthGraph()

TGraph * MuonCalib::VariableBinwidthHistogram::BinWidthGraph ( ) const

Plot binwidth graph - binwidth versus bin center.

Definition at line 127 of file VariableBinwidthHistogram.cxx.

127 {
128 Double_t *x = new Double_t[m_bins.size()], *y = new Double_t[m_bins.size()];
129 for (unsigned int i = 0; i < m_bins.size(); i++) {
130 x[i] = m_bins[i]->Center();
131 y[i] = m_bins[i]->Width();
132 }
133 TGraph *gr = new TGraph(m_bins.size(), x, y);
134 return gr;
135 }

◆ DenistyGraph()

TGraph * MuonCalib::VariableBinwidthHistogram::DenistyGraph ( ) const

create density graph - density vs bin center

Definition at line 113 of file VariableBinwidthHistogram.cxx.

113 {
114 Double_t *x = new Double_t[m_bins.size()], *y = new Double_t[m_bins.size()];
115 for (unsigned int i = 0; i < m_bins.size(); i++) {
116 x[i] = m_bins[i]->Center();
117 y[i] = m_bins[i]->Density();
118 }
119 TGraph *gr = new TGraph(m_bins.size(), x, y);
120 return gr;
121 }

◆ DiffBinwidthGraph()

TGraph * MuonCalib::VariableBinwidthHistogram::DiffBinwidthGraph ( ) const

Plot graph with differential binwidth.

Definition at line 166 of file VariableBinwidthHistogram.cxx.

166 {
167 if (m_bins.size() < 2) {
168 MsgStream log(Athena::getMessageSvc(), "VariableBinwidthHistogram");
169 log << MSG::WARNING << "DiffBinwidthGraph() - Need at alst 2 bins for differential density!" << endmsg;
170 return new TGraph();
171 }
172 Double_t *x = new Double_t[m_bins.size() - 1], *y = new Double_t[m_bins.size() - 1];
173 for (unsigned int i = 0; i < m_bins.size() - 1; i++) {
174 x[i] = m_bins[i]->Center();
175 y[i] = m_bins[i + 1]->Width() - m_bins[i]->Width();
176 }
177 TGraph *gr = new TGraph(m_bins.size() - 1, x, y);
178 return gr;
179 }
#define endmsg
IMessageSvc * getMessageSvc(bool quiet=false)

◆ DiffDensityGraph()

TGraph * MuonCalib::VariableBinwidthHistogram::DiffDensityGraph ( ) const

Plot graph with differential density.

Definition at line 151 of file VariableBinwidthHistogram.cxx.

151 {
152 if (m_bins.size() < 2) {
153 MsgStream log(Athena::getMessageSvc(), "VariableBinwidthHistogram");
154 log << MSG::WARNING << "DiffDensityGraph() - Need at alst 2 bins for differential density!" << endmsg;
155 return new TGraph();
156 }
157 Double_t *x = new Double_t[m_bins.size() - 1], *y = new Double_t[m_bins.size() - 1];
158 for (unsigned int i = 0; i < m_bins.size() - 1; i++) {
159 x[i] = m_bins[i]->Center();
160 y[i] = m_bins[i + 1]->Density() - m_bins[i]->Density();
161 }
162 TGraph *gr = new TGraph(m_bins.size() - 1, x, y);
163 return gr;
164 }

◆ GetBin()

const VariableBinwidthHistogramBin & MuonCalib::VariableBinwidthHistogram::GetBin ( unsigned int bin) const
inline

Get a bin.

Parameters
binThe bin index of the bin. bin < VariableBinwidthHistogram :: GetNumberOfBins()

Definition at line 71 of file VariableBinwidthHistogram.h.

71{ return *(m_bins[bin]); }

◆ GetNumberOfBins()

unsigned int MuonCalib::VariableBinwidthHistogram::GetNumberOfBins ( ) const
inline

Get the number of bins.

Definition at line 85 of file VariableBinwidthHistogram.h.

85{ return m_bins.size(); }

◆ GetSortedBin()

const VariableBinwidthHistogramBin & MuonCalib::VariableBinwidthHistogram::GetSortedBin ( unsigned int bin)
inline

Get a bin sorted by content.

Parameters
binThe bin index of the bin. bin < VariableBinwidthHistogram :: GetNumberOfBins()

Definition at line 76 of file VariableBinwidthHistogram.h.

76 {
77 if (!m_sorted) {
78 sort(m_sort_bins.begin(), m_sort_bins.end());
79 m_sorted = true;
80 }
81 return m_sort_bins[bin].Bin();
82 }
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
std::vector< VBHBinPtrSrt > m_sort_bins
bins sorted by content

◆ Initialize()

bool MuonCalib::VariableBinwidthHistogram::Initialize ( TH1F * hist,
double binc_r,
double max_bin_width,
double min_x = -9e9,
double max_x = 9e9 )

Initialize with new input histogram Returns on error false.

Parameters
histthe input histogram
binc_rnumber of entries per bin / max number of entries in hist
max_bin_widthThe binwidth will not exceed the maximum_bin_width
min_xonly the range between min_x and max_x will be used
max_xonly the range between min_x and max_x will be used

Definition at line 22 of file VariableBinwidthHistogram.cxx.

22 {
23 if (binc_rate < 1.0)
24 throw std::runtime_error(
25 Form("File: %s, Line: %d\nVariableBinwidthHistogram::Initialize - binc_rate must be greater than 1!", __FILE__, __LINE__));
26
27 // get maximum, firest bin and last bin
28 double max;
29 int first_bin(hist->FindBin(min_x)), last_bin(hist->FindBin(max_x));
30 bool max_valid(true);
31 if (first_bin == 0)
32 first_bin = 1;
33 else
34 max_valid = false;
35 if (last_bin > hist->GetNbinsX())
36 last_bin = hist->GetNbinsX();
37 else
38 max_valid = false;
39 if (max_valid)
40 max = hist->GetMaximum();
41 else {
42 max = -9e9;
43 for (int i = first_bin; i <= last_bin; i++) {
44 if (max < hist->GetBinContent(i)) max = hist->GetBinContent(i);
45 }
46 }
47 // get number of entries per bin
48 m_binc = static_cast<unsigned int>(ceil(max * binc_rate));
49 m_max_bin_width = max_bin_width;
50 // create first bin
51 m_bins.clear();
52 m_bins.push_back(new VariableBinwidthHistogramBin(hist->GetBinLowEdge(first_bin) + 0.5 * max_bin_width, max_bin_width, 0));
53 VariableBinwidthHistogramBin *current_bin(m_bins[0]);
54 // loop on histogram bins
55 for (int i = first_bin; i <= last_bin; i++) {
56 // maximum binwidth reached
57 if (hist->GetBinCenter(i) > current_bin->Right()) {
58 m_bins.push_back(new VariableBinwidthHistogramBin(current_bin->Right() + 0.5 * m_max_bin_width, max_bin_width, 0));
59 current_bin = m_bins[m_bins.size() - 1];
60 }
61 // will the bin be full
62 if (current_bin->Entries() + static_cast<unsigned int>(hist->GetBinContent(i)) > m_binc) {
63 // this will be filled in the next bin
64 unsigned int overflow = current_bin->Entries() + static_cast<unsigned int>(hist->GetBinContent(i)) - m_binc;
65 // the current bin ends here.
66 current_bin->MoveRight(hist->GetBinLowEdge(i));
67 current_bin->SetContent(m_binc);
68 // create new bin
69 m_bins.push_back(new VariableBinwidthHistogramBin(current_bin->Right() + 0.5 * m_max_bin_width, m_max_bin_width, overflow));
70 current_bin = m_bins[m_bins.size() - 1];
71 continue;
72 }
73 // add to content of current bin
74 (*current_bin) += static_cast<unsigned int>(hist->GetBinContent(i));
75 }
76 // sort bins by content
77 for (auto & bin : m_bins) { m_sort_bins.push_back(VBHBinPtrSrt(bin)); }
78 m_sorted = false;
79 m_error = false;
80 return true;
81 }
#define max(a, b)
Definition cfImp.cxx:41

◆ sign()

double MuonCalib::VariableBinwidthHistogram::sign ( double val) const
inlineprivate

Definition at line 122 of file VariableBinwidthHistogram.h.

122 {
123 if (val > 0.0) return 1.0;
124 if (val < 0.0) return -1.0;
125 if (val == 0.0) return 0;
126 return 0;
127 }

◆ Smooth()

bool MuonCalib::VariableBinwidthHistogram::Smooth ( double width)

Removes steps that origin in a binning effekt.

Parameters
widthsmoothin parameter (typically the binwidth of the input-hisatogram)

Definition at line 87 of file VariableBinwidthHistogram.cxx.

87 {
88 // needs at last 3 bins to smooth
89 if (m_bins.size() < 3) {
90 MsgStream log(Athena::getMessageSvc(), "VariableBinwidthHistogram");
91 log << MSG::WARNING << "Smooth() - VBH has less than 3 bins" << endmsg;
92 return false;
93 }
94 for (unsigned int i = 0; i < m_bins.size() - 3; i++) {
95 Double_t sl1 = m_bins[i + 1]->Width() - m_bins[i]->Width();
96 Double_t sl2 = m_bins[i + 2]->Width() - m_bins[i + 1]->Width();
97 // slopes must be oposit sign
98 if (sign(sl1) == sign(sl2)) continue;
99 // prevents numerical effects
100 if (std::abs(sl1) < width / 2 || std::abs(sl2) < width / 2) continue;
101 // move bin boarder
102 m_bins[i]->MoveRight(m_bins[i]->Right() - width / 2 * sign(sl2));
103 m_bins[i + 1]->MoveLeft(m_bins[i]->Right() - width / 2 * sign(sl2));
104 }
105 m_sorted = false;
106 return true;
107 }
const double width

Member Data Documentation

◆ m_binc

unsigned int MuonCalib::VariableBinwidthHistogram::m_binc
private

number of entries per bin

Definition at line 111 of file VariableBinwidthHistogram.h.

◆ m_bins

std::vector<VariableBinwidthHistogramBin *> MuonCalib::VariableBinwidthHistogram::m_bins
private

vector containing histogram bins

Definition at line 105 of file VariableBinwidthHistogram.h.

◆ m_error

bool MuonCalib::VariableBinwidthHistogram::m_error
private

error flag

Definition at line 117 of file VariableBinwidthHistogram.h.

◆ m_max_bin_width

double MuonCalib::VariableBinwidthHistogram::m_max_bin_width
private

maximum bin width

Definition at line 114 of file VariableBinwidthHistogram.h.

◆ m_sort_bins

std::vector<VBHBinPtrSrt> MuonCalib::VariableBinwidthHistogram::m_sort_bins
private

bins sorted by content

Definition at line 108 of file VariableBinwidthHistogram.h.

◆ m_sorted

bool MuonCalib::VariableBinwidthHistogram::m_sorted
private

Definition at line 119 of file VariableBinwidthHistogram.h.


The documentation for this class was generated from the following files: