ATLAS Offline Software
TFCS1DFunctionInt32Histogram.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 #include <algorithm>
7 #include <iostream>
8 #include "TMath.h"
9 #include "TCanvas.h"
10 #include "TH2F.h"
11 #include "TRandom.h"
12 #include "TFile.h"
13 #include "TClass.h"
14 
15 //=============================================
16 //======= TFCS1DFunctionInt32Histogram =========
17 //=============================================
18 
21 
23  Int_t nbinsx = hist->GetNbinsX();
24  Int_t nbins = nbinsx;
25 
26  float integral = 0;
27  m_HistoBorders.resize(nbinsx + 1);
28  m_HistoContents.resize(nbins);
29  std::vector<double> temp_HistoContents(nbins);
30  int ibin = 0;
31  for (int ix = 1; ix <= nbinsx; ix++) {
32  float binval = hist->GetBinContent(ix);
33  if (binval < 0) {
34  // Can't work if a bin is negative, forcing bins to 0 in this case
35  double fraction = binval / hist->Integral();
36  if (TMath::Abs(fraction) > 1e-5) {
37  ATH_MSG_WARNING("bin content is negative in histogram "
38  << hist->GetName() << " : " << hist->GetTitle()
39  << " binval=" << binval << " " << fraction * 100
40  << "% of integral=" << hist->Integral()
41  << ". Forcing bin to 0.");
42  }
43  binval = 0;
44  }
45  integral += binval;
46  temp_HistoContents[ibin] = integral;
47  ++ibin;
48  }
49  if (integral <= 0) {
50  ATH_MSG_ERROR("histogram " << hist->GetName() << " : " << hist->GetTitle()
51  << " integral=" << integral << " is <=0");
52  m_HistoBorders.resize(0);
53  m_HistoContents.resize(0);
54  return;
55  }
56 
57  for (int ix = 1; ix <= nbinsx; ix++)
58  m_HistoBorders[ix - 1] = hist->GetXaxis()->GetBinLowEdge(ix);
59  m_HistoBorders[nbinsx] = hist->GetXaxis()->GetXmax();
60 
61  for (ibin = 0; ibin < nbins; ++ibin) {
62  m_HistoContents[ibin] = s_MaxValue * (temp_HistoContents[ibin] / integral);
63  }
64 }
65 
66 double TFCS1DFunctionInt32Histogram::rnd_to_fct(double rnd) const {
67  if (m_HistoContents.empty()) {
68  return 0;
69  }
70  HistoContent_t int_rnd = s_MaxValue * rnd;
71  auto it =
72  std::upper_bound(m_HistoContents.begin(), m_HistoContents.end(), int_rnd);
73  int ibin = std::distance(m_HistoContents.begin(), it);
74  if (ibin >= (int)m_HistoContents.size())
75  ibin = m_HistoContents.size() - 1;
76  Int_t binx = ibin;
77 
78  HistoContent_t basecont = 0;
79  if (ibin > 0)
80  basecont = m_HistoContents[ibin - 1];
81 
82  HistoContent_t dcont = m_HistoContents[ibin] - basecont;
83  if (dcont > 0) {
84  return m_HistoBorders[binx] +
85  ((m_HistoBorders[binx + 1] - m_HistoBorders[binx]) *
86  (int_rnd - basecont)) /
87  dcont;
88  } else {
89  return m_HistoBorders[binx] +
90  (m_HistoBorders[binx + 1] - m_HistoBorders[binx]) / 2;
91  }
92 }
93 
95  if (IsA() != ref.IsA())
96  return false;
97  const TFCS1DFunctionInt32Histogram &ref_typed =
98  static_cast<const TFCS1DFunctionInt32Histogram &>(ref);
99 
100  if (m_HistoBorders != ref_typed.m_HistoBorders)
101  return false;
102  if (m_HistoContents != ref_typed.m_HistoContents)
103  return false;
104  return true;
105 }
106 
107 void TFCS1DFunctionInt32Histogram::unit_test ATLAS_NOT_THREAD_SAFE(TH1 *hist) {
109  int nbinsx;
110  if (hist == nullptr) {
111  nbinsx = 400;
112  hist = new TH1D("test1D", "test1D", nbinsx, 0, 1);
113  hist->Sumw2();
114  for (int ix = 1; ix <= nbinsx; ++ix) {
115  double val = (0.5 + gRandom->Rndm()) * (nbinsx + ix);
116  if (gRandom->Rndm() < 0.1)
117  val = 0;
118  hist->SetBinContent(ix, val);
119  hist->SetBinError(ix, 0);
120  }
121  }
123  nbinsx = hist->GetNbinsX();
124 
125  float value[2];
126  float rnd[2];
127  //cppcheck-suppress uninitvar
128  for (rnd[0] = 0; rnd[0] < 0.9999; rnd[0] += 0.25) {
129  rtof.rnd_to_fct(value, rnd);
130  ATH_MSG_NOCLASS(logger, "rnd0=" << rnd[0] << " -> x=" << value[0]);
131  }
132 
133  TH1 *hist_val = (TH1 *)hist->Clone("hist_val");
134  hist_val->SetTitle("difference");
135  hist_val->Reset();
136  int nrnd = 10000000;
137  double weight = hist->Integral() / nrnd;
138  hist_val->Sumw2();
139  for (int i = 0; i < nrnd; ++i) {
140  rnd[0] = gRandom->Rndm();
141  rtof.rnd_to_fct(value, rnd);
142  hist_val->Fill(value[0], weight);
143  }
144  hist_val->Add(hist, -1);
145 
146  TH1F *hist_pull = new TH1F("pull", "pull", 200, -10, 10);
147  for (int ix = 1; ix <= nbinsx; ++ix) {
148  float val = hist_val->GetBinContent(ix);
149  float err = hist_val->GetBinError(ix);
150  if (err > 0)
151  hist_pull->Fill(val / err);
152  ATH_MSG_NOCLASS(logger, "val=" << val << " err=" << err);
153  }
154 
155 // Screen output in athena won't make sense and would require linking of
156 // additional libraries
157 #if defined(__FastCaloSimStandAlone__)
158  new TCanvas("input", "Input");
159  hist->Draw();
160 
161  new TCanvas("validation", "Validation");
162  hist_val->Draw();
163 
164  new TCanvas("pull", "Pull");
165  hist_pull->Draw();
166 #endif
167 }
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
ISF_FCS::MLogging
Cut down AthMessaging.
Definition: MLogging.h:176
IsA
#define IsA
Declare the TObject style functions.
Definition: xAODTEventBranch.h:59
integral
double integral(TH1 *h)
Definition: computils.cxx:57
plotmaker.hist
hist
Definition: plotmaker.py:148
skel.it
it
Definition: skel.GENtoEVGEN.py:396
TFCS1DFunctionInt32Histogram::HistoContent_t
uint32_t HistoContent_t
Definition: TFCS1DFunctionInt32Histogram.h:27
athena.value
value
Definition: athena.py:124
SCT_CalibAlgs::nbins
@ nbins
Definition: SCT_CalibNumbers.h:10
TFCS1DFunctionInt32Histogram::rnd_to_fct
virtual double rnd_to_fct(double rnd) const
Function gets random number rnd in the range [0,1) as argument and returns function value according t...
Definition: TFCS1DFunctionInt32Histogram.cxx:66
dqt_zlumi_pandas.weight
int weight
Definition: dqt_zlumi_pandas.py:189
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
TFCS1DFunctionInt32Histogram::m_HistoBorders
std::vector< float > m_HistoBorders
Definition: TFCS1DFunctionInt32Histogram.h:48
dqt_zlumi_pandas.err
err
Definition: dqt_zlumi_pandas.py:182
lumiFormat.i
int i
Definition: lumiFormat.py:85
TFCS1DFunctionInt32Histogram
Definition: TFCS1DFunctionInt32Histogram.h:15
TFCS1DFunctionInt32Histogram::operator==
virtual bool operator==(const TFCS1DFunction &ref) const
The == operator compares the content of instances.
Definition: TFCS1DFunctionInt32Histogram.cxx:94
TFCS1DFunctionInt32Histogram::s_MaxValue
static const HistoContent_t s_MaxValue
Definition: TFCS1DFunctionInt32Histogram.h:28
ATH_MSG_NOCLASS
#define ATH_MSG_NOCLASS(logger_name, x)
Definition: MLogging.h:52
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
ref
const boost::regex ref(r_ef)
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
python.TrigEgammaMonitorHelper.TH1F
def TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
Definition: TrigEgammaMonitorHelper.py:24
TFCS1DFunctionInt32Histogram::m_HistoContents
std::vector< HistoContent_t > m_HistoContents
Definition: TFCS1DFunctionInt32Histogram.h:49
TFCS1DFunction
Definition: TFCS1DFunction.h:17
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
python.iconfTool.gui.pad.logger
logger
Definition: pad.py:14
TFCS1DFunctionInt32Histogram::Initialize
void Initialize(const TH1 *hist)
Definition: TFCS1DFunctionInt32Histogram.cxx:22
ATLAS_NOT_THREAD_SAFE
void TFCS1DFunctionInt32Histogram::unit_test ATLAS_NOT_THREAD_SAFE(TH1 *hist)
Definition: TFCS1DFunctionInt32Histogram.cxx:107
TFCS1DFunctionInt32Histogram.h