ATLAS Offline Software
Loading...
Searching...
No Matches
TFCS1DFunction.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include "TH1.h"
8#include "TCanvas.h"
9#include "TRandom.h"
10#include <string>
11#include <iostream>
12#include <cmath>
13
14//=============================================
15//======= TFCS1DFunction =========
16//=============================================
17
18void TFCS1DFunction::rnd_to_fct(float value[], const float rnd[]) const {
19 value[0] = rnd_to_fct(rnd[0]);
20}
21
22double TFCS1DFunction::get_maxdev(TH1 *h_input1, TH1 *h_approx1) {
23 TH1D *h_input = (TH1D *)h_input1->Clone("h_input");
24 TH1D *h_approx = (TH1D *)h_approx1->Clone("h_approx");
25
26 double maxdev = 0.0;
28 // normalize the histos to the same area:
29 double integral_input = h_input->Integral();
30 double integral_approx = 0.0;
31 for (int b = 1; b <= h_input->GetNbinsX(); b++){
32 integral_approx +=
33 h_approx->GetBinContent(h_approx->FindBin(h_input->GetBinCenter(b)));
34 }
35 if (integral_approx == 0.0) return 0.;
36 h_approx->Scale(integral_input / integral_approx);
37
38 double ymax = h_approx->GetBinContent(h_approx->GetNbinsX()) -
39 h_approx->GetBinContent(h_approx->GetMinimumBin());
40 for (int i = 1; i <= h_input->GetNbinsX(); i++) {
41 double val = std::abs(h_approx->GetBinContent(
42 h_approx->FindBin(h_input->GetBinCenter(i))) -
43 h_input->GetBinContent(i)) /
44 ymax;
45 if (val > maxdev)
46 maxdev = val;
47 }
48
49 delete h_input;
50 delete h_approx;
51
52 return maxdev * 100.0;
53}
54
56 const TH1 *hist, std::vector<double> &integral_vec, int &first, int &last) {
58 Int_t nbins = hist->GetNbinsX();
59
60 float integral = 0;
61 integral_vec.resize(nbins);
62 for (int ix = 1; ix <= nbins; ix++) {
63 float binval = hist->GetBinContent(ix);
64 if (binval < 0) {
65 // Can't work if a bin is negative, forcing bins to 0 in this case
66 double fraction = binval / hist->Integral();
67 if (std::abs(fraction) > 1e-5) {
68 ATH_MSG_NOCLASS(logger, "Warning: bin content is negative in histogram "
69 << hist->GetName() << " : "
70 << hist->GetTitle() << " binval=" << binval
71 << " " << fraction * 100
72 << "% of integral=" << hist->Integral()
73 << ". Forcing bin to 0.");
74 }
75 binval = 0;
76 }
77 integral += binval;
78 integral_vec[ix - 1] = integral;
79 }
80
81 for (first = 0; first < nbins; first++)
82 if (integral_vec[first] != 0)
83 break;
84 for (last = nbins - 1; last > 0; last--)
85 if (integral_vec[last] != integral)
86 break;
87 last++;
88
89 if (integral <= 0) {
90 ATH_MSG_NOCLASS(logger, "Error: histogram "
91 << hist->GetName() << " : " << hist->GetTitle()
92 << " integral=" << integral << " is <=0");
93 }
94 return integral;
95}
96
97TH1 *TFCS1DFunction::generate_histogram_random_slope
98ATLAS_NOT_THREAD_SAFE(int nbinsx, double xmin, double xmax,
99 double zerothreshold) {
100 TH1 *hist = new TH1D("test_slope1D", "test_slope1D", nbinsx, xmin, xmax);
101 hist->Sumw2();
102 for (int ix = 1; ix <= nbinsx; ++ix) {
103 double val = (0.5 + gRandom->Rndm()) * (nbinsx + ix);
104 if (gRandom->Rndm() < zerothreshold)
105 val = 0;
106 hist->SetBinContent(ix, val);
107 hist->SetBinError(ix, 0);
108 }
109 return hist;
110}
111
112TH1 *TFCS1DFunction::generate_histogram_random_gauss
113ATLAS_NOT_THREAD_SAFE(int nbinsx, int ntoy, double xmin, double xmax,
114 double xpeak, double sigma) {
115 TH1 *hist = new TH1D("test_gauss1D", "test_gauss1D", nbinsx, xmin, xmax);
116 hist->Sumw2();
117 for (int i = 1; i <= ntoy; ++i) {
118 double x = gRandom->Gaus(xpeak, sigma);
119 if (x >= xmin && x < xmax)
120 hist->Fill(x);
121 }
122 return hist;
123}
124
125void TFCS1DFunction::unit_test ATLAS_NOT_THREAD_SAFE(TH1 *hist,
126 TFCS1DFunction *rtof,
127 int nrnd, TH1 *histfine) {
129 ATH_MSG_NOCLASS(logger, "========= " << hist->GetName() << " funcsize="
130 << rtof->MemorySize() << " ========");
131 int nbinsx = hist->GetNbinsX();
132 double integral = hist->Integral();
133
134 float value[2];
135 float rnd[2];
136 //cppcheck-suppress uninitvar
137 for (rnd[0] = 0; rnd[0] < 0.9999; rnd[0] += 0.25) {
138 rtof->rnd_to_fct(value, rnd);
139 ATH_MSG_NOCLASS(logger, "rnd0=" << rnd[0] << " -> x=" << value[0]);
140 }
141
142 TH1 *hist_val = nullptr;
143 if (histfine)
144 hist_val = (TH1 *)histfine->Clone(TString(hist->GetName()) + "hist_val");
145 else
146 hist_val = (TH1 *)hist->Clone(TString(hist->GetName()) + "hist_val");
147 double weightfine = hist_val->Integral() / nrnd;
148 hist_val->SetTitle("toy simulation");
149 hist_val->Reset();
150 hist_val->SetLineColor(2);
151 hist_val->Sumw2();
152
153 TH1 *hist_diff = (TH1 *)hist->Clone(TString(hist->GetName()) + "_difference");
154 hist_diff->SetTitle("cut efficiency difference");
155 hist_diff->Reset();
156 hist_diff->Sumw2();
157
158 double weight = integral / nrnd;
159 for (int i = 0; i < nrnd; ++i) {
160 rnd[0] = gRandom->Rndm();
161 rtof->rnd_to_fct(value, rnd);
162 hist_val->Fill(value[0], weightfine);
163 hist_diff->Fill(value[0], weight);
164 }
165 hist_diff->Add(hist, -1);
166 hist_diff->Scale(1.0 / integral);
167
168 TH1F *hist_pull =
169 new TH1F(TString(hist->GetName()) + "_pull", "pull", 200, -10, 10);
170 for (int ix = 1; ix <= nbinsx; ++ix) {
171 float val = hist_diff->GetBinContent(ix);
172 float err = hist_diff->GetBinError(ix);
173 if (err > 0)
174 hist_pull->Fill(val / err);
175 }
176
177// Screen output in athena won't make sense and would require linking of
178// additional libraries
179#if defined(__FastCaloSimStandAlone__)
180 TCanvas *c = new TCanvas(hist->GetName(), hist->GetName());
181 c->Divide(2, 2);
182
183 c->cd(1);
184 if (histfine) {
185 histfine->SetLineColor(kGray);
186 histfine->DrawClone("hist");
187 hist->DrawClone("same");
188 } else {
189 hist->DrawClone();
190 }
191 hist_val->DrawClone("sameshist");
192
193 c->cd(2);
194 if (histfine) {
195 histfine->SetLineColor(kGray);
196 histfine->DrawClone("hist");
197 hist->DrawClone("same");
198 } else {
199 hist->DrawClone();
200 }
201 hist_val->DrawClone("sameshist");
202 gPad->SetLogy();
203
204 c->cd(3);
205 hist_diff->Draw();
206
207 c->cd(4);
208 hist_pull->Draw();
209
210 c->SaveAs(".png");
211#endif
212}
static TRandom * rnd
#define ATH_MSG_NOCLASS(logger_name, x)
Definition MLogging.h:52
#define x
#define ATLAS_NOT_THREAD_SAFE
getNoisyStrip() Find noisy strips from hitmaps and write out into xml/db formats
Cut down AthMessaging.
Definition MLogging.h:176
static double get_maxdev(TH1 *, TH1 *)
static double CheckAndIntegrate1DHistogram(const TH1 *hist, std::vector< double > &integral_vec, int &first, int &last)
virtual void rnd_to_fct(float value[], const float rnd[]) const
Function gets array of random numbers rnd[] in the range [0,1) as arguments and returns function valu...
virtual std::size_t MemorySize() const
Gives the total memory size, including the size of additional memory allocated inside the class.
double integral(TH1 *h)
Definition computils.cxx:59
static Root::TMsgLogger logger("iLumiCalc")
double xmax
Definition listroot.cxx:61
double xmin
Definition listroot.cxx:60
double ymax
Definition listroot.cxx:64