CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
th1fmorph.cc File Reference
#include "../interface/th1fmorph.h"
#include "TROOT.h"
#include "TAxis.h"
#include "TArrayD.h"
#include <iostream>
#include <cmath>
#include <set>

Go to the source code of this file.

Functions

TH1F * th1fmorph (const char *chname, const char *chtitle, TH1F *hist1, TH1F *hist2, Double_t par1, Double_t par2, Double_t parinterp, Double_t morphedhistnorm, Int_t idebug)
 
TH1D * th1fmorph (const char *chname, const char *chtitle, TH1D *hist1, TH1D *hist2, Double_t par1, Double_t par2, Double_t parinterp, Double_t morphedhistnorm, Int_t idebug)
 
template<typename TH1_t , typename Value_t >
TH1_t * th1fmorph_ (const char *chname, const char *chtitle, TH1_t *hist1, TH1_t *hist2, Double_t par1, Double_t par2, Double_t parinterp, Double_t morphedhistnorm, Int_t idebug)
 

Function Documentation

TH1F* th1fmorph ( const char *  chname,
const char *  chtitle,
TH1F *  hist1,
TH1F *  hist2,
Double_t  par1,
Double_t  par2,
Double_t  parinterp,
Double_t  morphedhistnorm,
Int_t  idebug 
)

Definition at line 489 of file th1fmorph.cc.

495 { return th1fmorph_<TH1F, Float_t>(chname, chtitle, hist1, hist2, par1, par2, parinterp, morphedhistnorm, idebug); }
TH1D* th1fmorph ( const char *  chname,
const char *  chtitle,
TH1D *  hist1,
TH1D *  hist2,
Double_t  par1,
Double_t  par2,
Double_t  parinterp,
Double_t  morphedhistnorm,
Int_t  idebug 
)

Definition at line 497 of file th1fmorph.cc.

503 { return th1fmorph_<TH1D, Double_t>(chname, chtitle, hist1, hist2, par1, par2, parinterp, morphedhistnorm, idebug); }
template<typename TH1_t , typename Value_t >
TH1_t* th1fmorph_ ( const char *  chname,
const char *  chtitle,
TH1_t *  hist1,
TH1_t *  hist2,
Double_t  par1,
Double_t  par2,
Double_t  parinterp,
Double_t  morphedhistnorm,
Int_t  idebug 
)

Definition at line 13 of file th1fmorph.cc.

References gather_cfg::cout, alignCSCRings::e, i, pileupDistInMC::total, x, and detailsBasic3DVector::y.

19 {
20  //--------------------------------------------------------------------------
21  // Author : Alex Read
22  // Version 0.2 of ROOT implementation, 08.05.2011
23  // *
24  // * Perform a linear interpolation between two histograms as a function
25  // * of the characteristic parameter of the distribution.
26  // *
27  // * The algorithm is described in Read, A. L., "Linear Interpolation
28  // * of Histograms", NIM A 425 (1999) 357-360.
29  // *
30  // * This ROOT-based CINT implementation is based on the FORTRAN77
31  // * implementation used by the DELPHI experiment at LEP (d_pvmorph.f).
32  // * The use of double precision allows pdf's to be accurately
33  // * interpolated down to something like 10**-15.
34  // *
35  // * The input histograms don't have to be identical, the binning is also
36  // * interpolated.
37  // *
38  // * Extrapolation is allowed (a warning is given) but the extrapolation
39  // * is not as well-defined as the interpolation and the results should
40  // * be used with great care.
41  // *
42  // * Data in the underflow and overflow bins are completely ignored.
43  // * They are neither interpolated nor do they contribute to the
44  // * normalization of the histograms.
45  // *
46  // * Input arguments:
47  // * ================
48  // * chname, chtitle : The ROOT name and title of the interpolated histogram.
49  // * Defaults for the name and title are "THF1-interpolated"
50  // * and "Interpolated histogram", respectively.
51  // *
52  // * hist1, hist2 : The two input histograms.
53  // *
54  // * par1,par2 : The values of the linear parameter that characterises
55  // * the histograms (e.g. a particle mass).
56  // *
57  // * parinterp : The value of the linear parameter we wish to
58  // * interpolate to.
59  // *
60  // * morphedhistnorm : The normalization of the interpolated histogram
61  // * (default is 1.0).
62  // *
63  // * idebug : Default is zero, no internal information displayed.
64  // * Values between 1 and increase the verbosity of
65  // * informational output which may prove helpful to
66  // * understand errors and pathalogical results.
67  // *
68  // * The routine returns a pointer (TH1_t *) to a new histogram which is
69  // * the interpolated result.
70  // *
71  // *------------------------------------------------------------------------
72  // Changes from 0.1 to 0.2:
73  // o Treatment of empty and non-existing histograms now well-defined.
74  // o The tricky regions of the first and last bins are improved (and
75  // well-tested).
76  // *------------------------------------------------------------------------
77 
78  // Return right away if one of the input histograms doesn't exist.
79  if(!hist1) {
80  cout << "ERROR! th1morph says first input histogram doesn't exist." << endl;
81  return(0);
82  }
83  if(!hist2) {
84  cout << "ERROR! th1morph says second input histogram doesn't exist." << endl;
85  return(0);
86  }
87 
88  // Extract bin parameters of input histograms 1 and 2.
89  // Supports the cases of non-equidistant as well as equidistant binning
90  // and also the case that binning of histograms 1 and 2 is different.
91  TAxis* axis1 = hist1->GetXaxis();
92  Int_t nb1 = axis1->GetNbins();
93  TAxis* axis2 = hist2->GetXaxis();
94  Int_t nb2 = axis2->GetNbins();
95 
96  std::set<Double_t> bedgesn_tmp;
97  for(Int_t i = 1; i <= nb1; ++i){
98  bedgesn_tmp.insert(axis1->GetBinLowEdge(i));
99  bedgesn_tmp.insert(axis1->GetBinUpEdge(i));
100  }
101  for(Int_t i = 1; i <= nb2; ++i){
102  bedgesn_tmp.insert(axis2->GetBinLowEdge(i));
103  bedgesn_tmp.insert(axis2->GetBinUpEdge(i));
104  }
105  Int_t nbn = bedgesn_tmp.size() - 1;
106  TArrayD bedgesn(nbn+1);
107  Int_t idx = 0;
108  for (std::set<Double_t>::const_iterator bedge = bedgesn_tmp.begin();
109  bedge != bedgesn_tmp.end(); ++bedge){
110  bedgesn[idx]=(*bedge);
111  ++idx;
112  }
113  Double_t xminn = bedgesn[0];
114  Double_t xmaxn = bedgesn[nbn];
115 
116  // ......The weights (wt1,wt2) are the complements of the "distances" between
117  // the values of the parameters at the histograms and the desired
118  // interpolation point. For example, wt1=0, wt2=1 means that the
119  // interpolated histogram should be identical to input histogram 2.
120  // Check that they make sense. If par1=par2 then we can choose any
121  // valid set of wt1,wt2 so why not take the average?
122 
123  Double_t wt1,wt2;
124  if (par2 != par1) {
125  wt1 = 1. - (parinterp-par1)/(par2-par1);
126  wt2 = 1. + (parinterp-par2)/(par2-par1);
127  }
128  else {
129  wt1 = 0.5;
130  wt2 = 0.5;
131  }
132 
133  //......Give a warning if this is an extrapolation.
134 
135  if (wt1 < 0 || wt1 > 1. || wt2 < 0. || wt2 > 1. || fabs(1-(wt1+wt2))
136  > 1.0e-4) {
137  cout << "Warning! th1fmorph: This is an extrapolation!! Weights are "
138  << wt1 << " and " << wt2 << " (sum=" << wt1+wt2 << ")" << endl;
139  }
140  if (idebug >= 1) cout << "th1morph - Weights: " << wt1 << " " << wt2 << endl;
141 
142  if (idebug >= 1) cout << "New hist: " << nbn << " " << xminn << " "
143  << xmaxn << endl;
144 
145  // Treatment for empty histograms: Return an empty histogram
146  // with interpolated bins.
147 
148  if (hist1->GetSum() <= 0 || hist2->GetSum() <=0 ) {
149  cout << "Warning! th1morph detects an empty input histogram. Empty interpolated histogram returned: "
150  <<endl << " " << chname << " - " << chtitle << endl;
151  TH1_t *morphedhist = (TH1_t *)gROOT->FindObject(chname);
152  if (morphedhist) delete morphedhist;
153  morphedhist = new TH1_t(chname,chtitle,nbn,xminn,xmaxn);
154  return(morphedhist);
155  }
156  if (idebug >= 1) cout << "Input histogram content sums: "
157  << hist1->GetSum() << " " << hist2->GetSum() << endl;
158  // *
159  // *......Extract the single precision histograms into double precision arrays
160  // * for the interpolation computation. The offset is because sigdis(i)
161  // * describes edge i (there are nbins+1 of them) while dist1/2
162  // * describe bin i. Be careful, ROOT does not use C++ convention to
163  // * number bins: dist1[ibin] is content of bin ibin where ibin runs from
164  // * 1 to nbins. We allocate some extra space for the derived distributions
165  // * because there may be as many as nb1+nb2+2 edges in the intermediate
166  // * interpolated cdf described by xdisn[i] (position of edge i) and
167  // * sigdisn[i] (cummulative probability up this edge) before we project
168  // * into the final binning.
169 
170  Value_t *dist1=hist1->GetArray();
171  Value_t *dist2=hist2->GetArray();
172  Double_t *sigdis1 = new Double_t[1+nb1];
173  Double_t *sigdis2 = new Double_t[1+nb2];
174  Double_t *sigdisn = new Double_t[2+nb1+nb2];
175  Double_t *xdisn = new Double_t[2+nb1+nb2];
176  Double_t *sigdisf = new Double_t[nbn+1];
177 
178  for(Int_t i=0;i<2+nb1+nb2;i++) xdisn[i] = 0; // Start with empty edges
179  sigdis1[0] = 0; sigdis2[0] = 0; // Start with cdf=0 at left edge
180 
181  for(Int_t i=1;i<nb1+1;i++) { // Remember, bin i has edges at i-1 and
182  sigdis1[i] = dist1[i]; // i and i runs from 1 to nb.
183  }
184  for(Int_t i=1;i<nb2+1;i++) {
185  sigdis2[i] = dist2[i];
186  }
187 
188  if (idebug >= 3) {
189  for(Int_t i=0;i<nb1+1;i++) {
190  cout << i << " dist1" << dist1[i] << endl;
191  }
192  for(Int_t i=0;i<nb2+1;i++) {
193  cout << i << " dist2" << dist2[i] << endl;
194  }
195  }
196 
197  //......Normalize the distributions to 1 to obtain pdf's and integrate
198  // (sum) to obtain cdf's.
199 
200  Double_t total = 0;
201  for(Int_t i=0;i<nb1+1;i++) {
202  total += sigdis1[i];
203  }
204  if (idebug >=1) cout << "Total histogram 1: " << total << endl;
205  for(Int_t i=1;i<nb1+1;i++) {
206  sigdis1[i] = sigdis1[i]/total + sigdis1[i-1];
207  }
208 
209  total = 0.;
210  for(Int_t i=0;i<nb2+1;i++) {
211  total += sigdis2[i];
212  }
213  if (idebug >=1) cout << "Total histogram 22: " << total << endl;
214  for(Int_t i=1;i<nb2+1;i++) {
215  sigdis2[i] = sigdis2[i]/total + sigdis2[i-1];
216  }
217 
218  // *
219  // *......We are going to step through all the edges of both input
220  // * cdf's ordered by increasing value of y. We start at the
221  // * lower edge, but first we should identify the upper ends of the
222  // * curves. These (ixl1, ixl2) are the first point in each cdf from
223  // * above that has the same integral as the last edge.
224  // *
225 
226  Int_t ix1l = nb1;
227  Int_t ix2l = nb2;
228  while(sigdis1[ix1l-1] >= sigdis1[ix1l]) {
229  ix1l = ix1l - 1;
230  }
231  while(sigdis2[ix2l-1] >= sigdis2[ix2l]) {
232  ix2l = ix2l - 1;
233  }
234 
235  // *
236  // *......Step up to the beginnings of the curves. These (ix1, ix2) are the
237  // * first non-zero points from below.
238 
239  Int_t ix1 = -1;
240  do {
241  ix1 = ix1 + 1;
242  } while(sigdis1[ix1+1] <= sigdis1[0]);
243 
244  Int_t ix2 = -1;
245  do {
246  ix2 = ix2 + 1;
247  } while(sigdis2[ix2+1] <= sigdis2[0]);
248 
249  if (idebug >= 1) {
250  cout << "First and last edge of hist1: " << ix1 << " " << ix1l << endl;
251  cout << " " << sigdis1[ix1] << " " << sigdis1[ix1+1] << endl;
252  cout << "First and last edge of hist2: " << ix2 << " " << ix2l << endl;
253  cout << " " << sigdis2[ix2] << " " << sigdis2[ix2+1] << endl;
254  }
255 
256  //......The first interpolated point should be computed now.
257 
258  Int_t nx3 = 0;
259  Double_t x1,x2,x;
260  x1 = axis1->GetBinLowEdge(ix1+1);
261  x2 = axis2->GetBinLowEdge(ix2+1);
262  x = wt1*x1 + wt2*x2;
263  xdisn[nx3] = x;
264  sigdisn[nx3] = 0;
265  if(idebug >= 1) {
266  cout << "First interpolated point: " << xdisn[nx3] << " "
267  << sigdisn[nx3] << endl;
268  cout << " " << x1 << " <= " << x << " <= "
269  << x2 << endl;
270  }
271 
272  //......Loop over the remaining point in both curves. Getting the last
273  // points may be a bit tricky due to limited floating point
274  // precision.
275 
276  if (idebug >= 1) {
277  cout << "----BEFORE while with ix1=" << ix1 << ", ix1l=" << ix1l
278  << ", ix2=" << ix2 << ", ix2l=" << ix2l << endl;
279  }
280 
281  Double_t yprev = -1; // The probability y of the previous point, it will
282  //get updated and used in the loop.
283  Double_t y = 0;
284  while((ix1 < ix1l) | (ix2 < ix2l)) {
285  if (idebug >= 1 ) cout << "----Top of while with ix1=" << ix1
286  << ", ix1l=" << ix1l << ", ix2=" << ix2
287  << ", ix2l=" << ix2l << endl;
288 
289  //......Increment to the next lowest point. Step up to the next
290  // kink in case there are several empty (flat in the integral)
291  // bins.
292 
293  Int_t i12type = -1; // Tells which input distribution we need to
294  // see next point of.
295  if ((sigdis1[ix1+1] <= sigdis2[ix2+1] || ix2 == ix2l) && ix1 < ix1l) {
296  ix1 = ix1 + 1;
297  while(sigdis1[ix1+1] <= sigdis1[ix1] && ix1 < ix1l) {
298  ix1 = ix1 + 1;
299  }
300  i12type = 1;
301  } else if (ix2 < ix2l) {
302  ix2 = ix2 + 1;
303  while(sigdis2[ix2+1] <= sigdis2[ix2] && ix2 < ix2l) {
304  ix2 = ix2 + 1;
305  }
306  i12type = 2;
307  }
308  if (i12type == 1) {
309  if (idebug >= 3) {
310  cout << "Pair for i12type=1: " << sigdis2[ix2] << " "
311  << sigdis1[ix1] << " " << sigdis2[ix2+1] << endl;
312  }
313  x1 = axis1->GetBinLowEdge(ix1+1);
314  y = sigdis1[ix1];
315  Double_t x20 = axis2->GetBinLowEdge(ix2+1);
316  Double_t x21 = axis2->GetBinUpEdge(ix2+1);
317  Double_t y20 = sigdis2[ix2];
318  Double_t y21 = sigdis2[ix2+1];
319 
320  //......Calculate where the cummulative probability y in distribution 1
321  // intersects between the 2 points from distribution 2 which
322  // bracket it.
323 
324  if (y21 > y20) {
325  x2 = x20 + (x21-x20)*(y-y20)/(y21-y20);
326  }
327  else {
328  x2 = x20;
329  }
330  } else {
331  if (idebug >= 3) {
332  cout << "Pair for i12type=2: " << sigdis1[ix1] << " " << sigdis2[ix2]
333  << " " << sigdis1[ix1+1] << endl;
334  }
335  x2 = axis2->GetBinLowEdge(ix2+1);
336  y = sigdis2[ix2];
337  Double_t x10 = axis1->GetBinLowEdge(ix1+1);
338  Double_t x11 = axis1->GetBinUpEdge(ix1+1);
339  Double_t y10 = sigdis1[ix1];
340  Double_t y11 = sigdis1[ix1+1];
341 
342  //......Calculate where the cummulative probability y in distribution 2
343  // intersects between the 2 points from distribution 1 which
344  // brackets it.
345 
346  if (y11 > y10) {
347  x1 = x10 + (x11-x10)*(y-y10)/(y11-y10);
348  } else {
349  x1 = x10;
350  }
351  }
352 
353  //......Interpolate between the x's in the 2 distributions at the
354  // cummulative probability y. Store the (x,y) for provisional
355  // edge nx3 in (xdisn[nx3],sigdisn[nx3]). nx3 grows for each point
356  // we add the the arrays. Note: Should probably turn the pair into
357  // a structure to make the code more object-oriented and readable.
358 
359  x = wt1*x1 + wt2*x2;
360  if (y > yprev) {
361  nx3 = nx3+1;
362  if (idebug >= 1) {
363  cout << " ---> y > yprev: i12type=" << i12type << ", nx3="
364  << nx3 << ", x= " << x << ", y=" << y << ", yprev=" << yprev
365  << endl;
366  }
367  yprev = y;
368  xdisn[nx3] = x;
369  sigdisn[nx3] = y;
370  if(idebug >= 1) {
371  cout << " ix1=" << ix1 << ", ix2= " << ix2 << ", i12type= "
372  << i12type << ", sigdis1[ix1]=" << sigdis1[ix1] << endl;
373  cout << " " << ", nx3=" << nx3 << ", x=" << x << ", y= "
374  << sigdisn[nx3] << endl;
375  }
376  }
377  }
378  if (idebug >=3) for (Int_t i=0;i<nx3;i++) {
379  cout << " nx " << i << " " << xdisn[i] << " " << sigdisn[i] << endl;
380  }
381 
382  // *......Now we loop over the edges of the bins of the interpolated
383  // * histogram and find out where the interpolated cdf 3
384  // * crosses them. This projection defines the result and will
385  // * be stored (after differention and renormalization) in the
386  // * output histogram.
387  // *
388  // *......We set all the bins following the final edge to the value
389  // * of the final edge.
390 
391  x = xmaxn;
392  Int_t ix = nbn;
393 
394  if (idebug >= 1) cout << "------> Any final bins to set? " << x << " "
395  << xdisn[nx3] << endl;
396  while(x >= xdisn[nx3]) {
397  sigdisf[ix] = sigdisn[nx3];
398  if (idebug >= 2) cout << " Setting final bins" << ix << " " << x
399  << " " << sigdisf[ix] << endl;
400  ix = ix-1;
401  x = bedgesn[ix];
402  }
403  Int_t ixl = ix + 1;
404  if (idebug >= 1) cout << " Now ixl=" << ixl << " ix=" << ix << endl;
405 
406  // *
407  // *......The beginning may be empty, so we have to step up to the first
408  // * edge where the result is nonzero. We zero the bins which have
409  // * and upper (!) edge which is below the first point of the
410  // * cummulative distribution we are going to project to this
411  // * output histogram binning.
412  // *
413 
414  ix = 0;
415  x = bedgesn[ix+1];
416  if (idebug >= 1) cout << "Start setting initial bins at x=" << x << endl;
417  while(x <= xdisn[0]) {
418  sigdisf[ix] = sigdisn[0];
419  if (idebug >= 1) cout << " Setting initial bins " << ix << " " << x
420  << " " << xdisn[1] << " " << sigdisf[ix] << endl;
421  ix = ix+1;
422  x = bedgesn[ix+1];
423  }
424  Int_t ixf = ix;
425 
426  if (idebug >= 1)
427  cout << "Bins left to loop over:" << ixf << "-" << ixl << endl;
428 
429  // *......Also the end (from y to 1.0) often comes before the last edge
430  // * so we have to set the following to 1.0 as well.
431 
432  Int_t ix3 = 0; // Problems with initial edge!!!
433  for(ix=ixf;ix<ixl;ix++) {
434  x = bedgesn[ix];
435  if (x < xdisn[0]) {
436  y = 0;
437  } else if (x > xdisn[nx3]) {
438  y = 1.;
439  } else {
440  while(xdisn[ix3+1] <= x && ix3 < 2*nbn) {
441  ix3 = ix3 + 1;
442  }
443  Double_t dx2=axis2->GetBinWidth(axis2->FindBin(x));
444  if (xdisn[ix3+1]-x > 1.1*dx2) { // Empty bin treatment
445  y = sigdisn[ix3+1];
446  }
447  else if (xdisn[ix3+1] > xdisn[ix3]) { // Normal bins
448  y = sigdisn[ix3] + (sigdisn[ix3+1]-sigdisn[ix3])
449  *(x-xdisn[ix3])/(xdisn[ix3+1]-xdisn[ix3]);
450  } else { // Is this ever used?
451  y = 0;
452  cout << "Warning - th1fmorph: This probably shoudn't happen! "
453  << endl;
454  cout << "Warning - th1fmorph: Zero slope solving x(y)" << endl;
455  }
456  }
457  sigdisf[ix] = y;
458  if (idebug >= 3) {
459  cout << ix << ", ix3=" << ix3 << ", xdisn=" << xdisn[ix3] << ", x="
460  << x << ", next xdisn=" << xdisn[ix3+1] << endl;
461  cout << " cdf n=" << sigdisn[ix3] << ", y=" << y << ", next point="
462  << sigdisn[ix3+1] << endl;
463  }
464  }
465 
466  //......Differentiate interpolated cdf and return renormalized result in
467  // new histogram.
468 
469  TH1_t *morphedhist = (TH1_t *)gROOT->FindObject(chname);
470  if (morphedhist) delete morphedhist;
471  morphedhist = new TH1_t(chname,chtitle,nbn,bedgesn.GetArray());
472 
473  for(ix=nbn-1;ix>-1;ix--) {
474  y = sigdisf[ix+1]-sigdisf[ix];
475  morphedhist->SetBinContent(ix+1,y*morphedhistnorm);
476  }
477 
478  //......Clean up the temporary arrays we allocated.
479 
480  delete sigdis1; delete sigdis2;
481  delete sigdisn; delete xdisn; delete sigdisf;
482 
483  //......All done, return the result.
484 
485  return(morphedhist);
486 }
int i
Definition: DBlmapReader.cc:9
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10