ATLAS Offline Software
Loading...
Searching...
No Matches
Resplot.cxx File Reference
#include <sstream>
#include <cmath>
#include "Directory.h"
#include "tagname.h"
#include "Resplot.h"
#include "TCanvas.h"
#include "TStyle.h"
#include "TMath.h"
#include "generate.h"
#include "rmsFrac.h"

Go to the source code of this file.

Functions

void binwidth (TH1D *h)
void ZeroErrors (TH1D *h)
void unZeroErrors (TH1D *h)
double getWeights (TH1D *h)
double getWeights (TH2D *h)
void GetStats (TH1D *h, std::vector< double > &stats)
std::string number (const double &x)
double NCounter (TH1D *h)
double FindMean (TH1D *s, double frac=0.95)
int findmax (TH1D *s)
Double_t langaufun (Double_t *x_par, Double_t *par)

Variables

int nexperiments_ = 0
int nexperiments_max = 40
int nexperiments_min = 20

Detailed Description

Author
M.Sutton
Date
Mon Jun 21 18:35:22 BST 2004

Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration

Definition in file Resplot.cxx.

Function Documentation

◆ binwidth()

void binwidth ( TH1D * h)

Definition at line 49 of file Resplot.cxx.

49 {
50
51 for (int i=1 ; i<=h->GetNbinsX() ; i++ ) {
52 double d = h->GetBinLowEdge(i+1) - h->GetBinLowEdge(i);
53
54 h->SetBinContent(i, h->GetBinContent(i)/d );
55 h->SetBinError(i, h->GetBinError(i)/d );
56 }
57
58}
Header file for AthHistogramAlgorithm.

◆ findmax()

int findmax ( TH1D * s)

find maximum bin

Definition at line 1395 of file Resplot.cxx.

1395 {
1396 int imax = 1;
1397 for ( int i=2 ; i<=s->GetNbinsX() ; i++ ) {
1399 if ( s->GetBinContent(i)>s->GetBinContent(imax) ) imax = i;
1400 }
1401 return imax;
1402}
int imax(int i, int j)

◆ FindMean()

double FindMean ( TH1D * s,
double frac = 0.95 )

Definition at line 1308 of file Resplot.cxx.

1308 {
1309
1310 // double entries = s->GetEntries() + s->GetBinContent(0) + s->GetBinContent(s->GetNbinsX()+1);
1311 double entries = generate::GetEntries(s,0,s->GetNbinsX()+1);
1312
1313 // std::cout << s->GetName() << "\tentries: " << entries << std::endl;
1314
1315 if ( entries==0 ) return 0;
1316
1317 // s->GetXaxis()->SetRangeUser(-100, 100);
1318 // s->GetXaxis()->UnZoom();
1319 s->GetXaxis()->SetRange(1,s->GetNbinsX());
1320
1321 double mean = s->GetMean();
1322 double meane = s->GetMeanError();
1323
1324 int upperbin = 0;
1325 int lowerbin = 0;
1326
1327
1328 for ( int it=0 ; it<20 ; it++ ) {
1329
1330 // std::cout << it << "\tmean " << mean << " +- " << meane << std::endl;
1331
1332 int imax = s->GetXaxis()->FindBin(mean);
1333
1334 double sumn = s->GetBinContent(imax);
1335
1336 // double uppersum = 0;
1337 // double lowersum = 0;
1338
1339 upperbin = imax;
1340 lowerbin = imax;
1341
1342
1343 int i=1;
1344 while ( true ) {
1345
1346 const int upperbin_i = imax+i;
1347 const int lowerbin_i = imax-i;
1348
1349 if ( upperbin_i>s->GetNbinsX() || lowerbin_i<1 ) break;
1350
1351 double tsum = sumn + s->GetBinContent(upperbin_i) + s->GetBinContent(lowerbin_i);
1352
1353 if ( tsum>entries*frac ) {
1354 // uppersum = tsum/entries;
1355 break;
1356 }
1357
1358 // lowersum = tsum/entries;
1359
1360 sumn = tsum;
1361
1362 upperbin = upperbin_i;
1363 lowerbin = lowerbin_i;
1364
1365 i++;
1366 }
1367
1368 s->GetXaxis()->SetRange(lowerbin, upperbin);
1369
1370 double m = s->GetMean();
1371 double me = s->GetMeanError();
1372
1373 if ( it>0 ) {
1374 if ( mean==m ||
1375 std::fabs(mean-m)<me*1e-5 ||
1376 std::fabs(mean-m)<meane*1e-5 ) {
1377 mean = m;
1378 meane = me;
1379 break;
1380 }
1381 }
1382
1383 mean = m;
1384 meane = me;
1385
1386 }
1387
1388 // std::cout << "mean " << mean << " +- " << meane << "\tentries: " << generate::GetEntries(s) << "\tlower " << lowerbin << " - " << upperbin << std::endl;
1389
1390
1391 return mean;
1392}
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
double entries
Definition listroot.cxx:49
double GetEntries(TH1D *h, int ilow, int ihi)
Definition rmsFrac.cxx:20

◆ GetStats()

void GetStats ( TH1D * h,
std::vector< double > & stats )

Definition at line 112 of file Resplot.cxx.

112 {
113
114 stats.clear();
115 stats.resize(4);
116
117 TAxis& fXaxis = *h->GetXaxis();
118
119 double s = 0;
120 double sx = 0;
121 double sx2 = 0;
122 double sx4 = 0;
123
124
125 // if ( fXaxis.TestBit(TAxis::kAxisRange) ) {
126 // for (bin=0;bin<4;bin++) stats[bin] = 0;
127
128 Int_t firstBinX = fXaxis.GetFirst();
129 Int_t lastBinX = fXaxis.GetLast();
130 // include underflow/overflow if TH1::StatOverflows(kTRUE) in case no range is set on the axis
131 // if (fgStatOverflows && !fXaxis.TestBit(TAxis::kAxisRange)) {
132 // if (firstBinX == 1) firstBinX = 0;
133 // if (lastBinX == fXaxis.GetNbins() ) lastBinX += 1;
134 // }
135 for (int binx = firstBinX; binx <= lastBinX; binx++) {
136 double x = fXaxis.GetBinCenter(binx);
137 double w = TMath::Abs(h->GetBinContent(binx));
138 // err = TMath::Abs(GetBinError(binx));
139 s += w;
140 // stats[1] += err*err;
141 sx += w*x;
142 sx2 += w*x*x;
143 sx4 += w*x*x*x*x;
144 }
145
146
147 double mu = sx/s;
148
149 double v = sx2/s-mu*mu;
150
151 double rms = std::sqrt(v);
152 double rmserror = std::sqrt(0.25*(sx4/v-v*s))/s;
153 // double rmserror1 = std::sqrt(0.4*(sx4/s-v*v)/s); << is this the correct one? it seems about 2* too large
154
155 stats[0] = mu;
156 stats[1] = std::sqrt(v/s);
157
158 stats[2] = rms;
159 stats[3] = rmserror;
160
161 // double duff = std::sqrt(0.5*v/s);
162
163#if 0
164 std::cout << "GetStats() "
165 // << "\tmean " << stats[0] << " +- " << stats[1]
166 << "\trms " << stats[2] << " +- " << stats[3] << "\t(" << rmserror1 << ")" << "\t root " << duff << std::endl;
167#endif
168
169}
#define x

◆ getWeights() [1/2]

double getWeights ( TH1D * h)

Definition at line 78 of file Resplot.cxx.

78 {
79 double sum = 0;
80 double sume2 = 0;
81 for ( int i=0 ; i<h->GetNbinsX()+1 ; i++ ) {
82 sum += h->GetBinContent(i);
83 sume2 += h->GetBinError(i)*h->GetBinError(i);
84 }
85 if ( sume2!=0 ) return sum*sum/sume2;
86 return 0;
87}

◆ getWeights() [2/2]

double getWeights ( TH2D * h)

Definition at line 90 of file Resplot.cxx.

90 {
91
92 double sum = 0;
93 double sume2 = 0;
94
95 for ( int i=1 ; i<=h->GetNbinsX() ; i++ ) {
96 for ( int j=1 ; j<=h->GetNbinsY() ; j++ ) {
97 sum += h->GetBinContent(i,j);
98 sume2 += h->GetBinError(i,j)*h->GetBinError(i,j);
99 }
100 }
101
102 if ( sume2!=0 ) return sum*sum/sume2;
103
104 return 0;
105}

◆ langaufun()

Double_t langaufun ( Double_t * x_par,
Double_t * par )

Definition at line 2017 of file Resplot.cxx.

2017 {
2018
2019 //Fit parameters:
2020 //par[0]=Total area (integral -inf to inf, normalization constant)
2021 //par[1]=Most Probable (MP, location) parameter of Landau density
2022 //par[2]=Width (scale) parameter of Landau density
2023 //par[3]=Width (sigma) of convoluted Gaussian function
2024 //
2025 //In the Landau distribution (represented by the CERNLIB approximation),
2026 //the maximum is located at x=-0.22278298 with the location parameter=0.
2027 //This shift is corrected within this function, so that the actual
2028 //maximum is identical to the MP parameter.
2029
2030 // Numeric constants
2031 Double_t invsq2pi = 0.3989422804014; // (2 pi)^(-1/2)
2032 Double_t mpshift = -0.22278298; // Landau maximum location
2033
2034 // Control constants
2035 Double_t np = 500; // number of convolution steps
2036 Double_t sc = 5; // convolution extends to +-sc Gaussian sigmas
2037
2038 double& x = x_par[0];
2039
2040 // Variables
2041 Double_t xx;
2042 Double_t mpc;
2043 Double_t fland;
2044 Double_t sum = 0.0;
2045 Double_t xlow,xupp;
2046 Double_t step;
2047 Double_t i;
2048
2049 // if ( par[2]<0.1*par[3] ) return 0;
2050
2051 // correct the most probable bin position as described above
2052 mpc = par[1] - mpshift*par[2];
2053
2054 // Range of convolution integral
2055
2056 // if ( par[3]<par[2] ) np = np*par[2]/par[3];
2057
2058 xlow = x - sc * par[3];
2059 xupp = x + sc * par[3];
2060
2061 step = (xupp-xlow) / np;
2062
2063 // Convolution integral of Landau and Gaussian by sum
2064 for(i=1.0; i<=np/2; i++) {
2065 xx = xlow + (i-0.5)*step;
2066 fland = TMath::Landau(xx,mpc,par[2]) / par[2];
2067 sum += fland * TMath::Gaus(x,xx,par[3]);
2068
2069 xx = xupp - (i-0.5)*step;
2070 fland = TMath::Landau(xx,mpc,par[2]) / par[2];
2071 sum += fland * TMath::Gaus(x,xx,par[3]);
2072 }
2073
2074 return (par[0] * step * sum * invsq2pi / par[3]);
2075}
static Double_t sc

◆ NCounter()

double NCounter ( TH1D * h)

Definition at line 379 of file Resplot.cxx.

379 {
380 double N=0;
381 for ( int i=1 ; i<=h->GetNbinsX() ; i++ ) N += h->GetBinContent(i);
382 return N;
383}

◆ number()

std::string number ( const double & x)

Definition at line 372 of file Resplot.cxx.

372 {
373 std::ostringstream s;
374 s << x;
375 return s.str();
376}

◆ unZeroErrors()

void unZeroErrors ( TH1D * h)

Definition at line 72 of file Resplot.cxx.

72 {
73 for (int i=1 ; i<=h->GetNbinsX() ; i++ ) if (h->GetBinContent(i)==0) h->SetBinError(i,0);
74}

◆ ZeroErrors()

void ZeroErrors ( TH1D * h)

Definition at line 66 of file Resplot.cxx.

66 {
67 for (int i=1 ; i<=h->GetNbinsX() ; i++ ) if (h->GetBinContent(i)==0) h->SetBinError(i,1);
68}

Variable Documentation

◆ nexperiments_

int nexperiments_ = 0

Definition at line 1680 of file Resplot.cxx.

◆ nexperiments_max

int nexperiments_max = 40

Definition at line 1681 of file Resplot.cxx.

◆ nexperiments_min

int nexperiments_min = 20

Definition at line 1682 of file Resplot.cxx.