9 This module contains set of simple and useful utilities for booking and 10 manipulations with Gaudi-AIDA histograms, inspired by Thomas' request 12 The module contains following public symbols: 14 - book for booking of various 1D,2D&3D-histograms 15 - bookProf for booking of various 1D&2D-profiles 16 - getAsAIDA for retrieval of histograms/profiles from HTS in AIDA format 17 - getAsROOT for retrieval of histograms/profiles from HTS in ROOT format 18 - fill for smart filling of 1D-histogram (AIDA or ROOT) 19 - aida2root for conversion of AIDA histogram to ROOT 20 - HistoStats for many statistical information 21 - HistoFile class for storing histograms to a file 22 - histoDump for dumping of the histogram in text format (a'la HBOOK) 23 - dumpHisto for dumping of the histogram in text format (a'la HBOOK) 26 from __future__
import print_function
28 __author__ =
"Vanya BELYAEV ibelyaev@physics.syr.edu" 59 Helper private auxiliary function to get Application Manager 61 gaudi = kwargs.get(
'gaudi',
None)
65 raise RuntimeError(
'Unable to get valid ApplicationMgr')
67 state = gaudi._isvc.FSMState()
68 if state < cpp.Gaudi.StateMachine.CONFIGURED:
70 state = gaudi._isvc.FSMState()
71 if state < cpp.Gaudi.StateMachine.INITIALIZED:
83 Helper private auxiliary function to get iHistogramSvs 85 svc = kwargs.get(
'service',
None)
87 svc = kwargs.get(
'svc',
None)
91 return gaudi.histsvc()
100 Helper private auxiliary function to get iDataSvs 102 svc = kwargs.get(
'service',
None)
104 svc = kwargs.get(
'svc',
None)
108 return gaudi.evtsvc()
117 The trivial function to book the various 1D,2D&3D-histograms 119 (1) book the trivial 1D histogram with full path 121 >>> h1D = book ( 'path/to/my/histo' , ## path in Histogram Transient Store 122 'cosine of decay angle ' , ## histogram title 123 100 , ## number of bins 127 (2) book the trivial 1D histogram with directory path and string ID : 129 >>> h1D = book ( 'path/to/directory' , ## the path to directory in HTS 130 'H1' , ## string histogram identifier 131 'cosine of decay angle ' , ## histogram title 132 100 , ## number of bins 136 (3) book the trivial 1D histogram with directory path and integer ID : 138 >>> h1D = book ( 'path/to/directory' , ## the path to directory in HTS 139 124 , ## integer histogram identifier 140 'cosine of decay angle ' , ## histogram title 141 100 , ## number of bins 145 (4) book the trivial 2D histogram with full path 147 >>> h1D = book ( 'path/to/my/histo' , ## path in Histogram Transient Store 148 'm12**2 versus m13**2' , ## histogram title 149 50 , ## number of X-bins 152 50 , ## number of Y-bins 156 (5) book the trivial 2D histogram with directory path and literal ID 158 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store 159 'Dalitz1' , ## literal histogram identifier 160 'm12**2 versus m13**2' , ## histogram title 161 50 , ## number of X-bins 164 50 , ## number of Y-bins 168 (6) book the trivial 2D histogram with directory path and integer ID 170 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store 171 854 , ## integer histogram identifier 172 'm12**2 versus m13**2' , ## histogram title 173 50 , ## number of X-bins 176 50 , ## number of Y-bins 180 (7) book the trivial 3D histogram with full path 182 >>> h1D = book ( 'path/to/my/histo' , ## path in Histogram Transient Store 183 'x vs y vs z ' , ## histogram title 184 10 , ## number of X-bins 187 10 , ## number of Y-bins 190 10 , ## number of Z-bins 194 (8) book the trivial 3D histogram with directory path and literal ID 196 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store 197 'xyz' , ## literal histogram identifier 198 'x vs y vs z' , ## histogram title 199 10 , ## number of X-bins 202 10 , ## number of Y-bins 205 10 , ## number of Z-bins 209 (9) book the trivial 3D histogram with directory path and integer ID 211 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store 212 888 , ## integer histogram identifier 213 'x vs y vs z' , ## histogram title 214 10 , ## number of X-bins 217 10 , ## number of Y-bins 220 10 , ## number of Z-bins 224 Many other booking methods are available, 225 e.g. for the histograms with non-equidistant bins, see IHistogamSvc::book 228 if useROOT
or kwargs.get(
'useROOT',
229 False)
or not kwargs.get(
'useAIDA',
True):
230 from ROOT
import TH1D
234 if not str
is type(a1):
237 return TH1D(a0 + a1, a2, *args[3:])
239 return TH1D(a0, a1, *args[2:])
243 raise RuntimeError(
'Unable to get valid HistogramService ')
245 return svc.book(*args)
248 book.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.book) : \n\n' \
249 + iHistogramSvc.book . __doc__
250 book.__doc__ +=
'\n\n' +
'\thelp(IHistogramSvc::book) : \n\n' \
251 + cpp.IHistogramSvc.book . __doc__
260 The trivial function to book 1D&2D profile histograms: 262 (1) book 1D-profile histogram with full path in Histogram Transient Store: 263 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store 264 'Energy Correction' , ## the histogram title 265 100 , ## number of X-bins 269 (2) book 1D-profile histogram with directory path and literal ID 270 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store 271 'Calibration' , ## the histogram literal identifier 272 'Energy Correction' , ## the histogram title 273 100 , ## number of X-bins 277 (3) book 1D-profile histogram with directory path and integer ID 278 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store 279 418 , ## the histogram integer identifier 280 'Energy Correction' , ## the histogram title 281 100 , ## number of X-bins 285 (4) book 2D-profile histogram with full path in Histogram Transient Store: 286 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store 287 'Energy Correction' , ## the histogram title 288 50 , ## number of X-bins 291 50 , ## number of Y-bins 295 (5) book 2D-profile histogram with directory path and literal ID 296 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store 297 'Calibration' , ## the histogram literal identifier 298 'Energy Correction' , ## the histogram title 299 50 , ## number of X-bins 302 50 , ## number of Y-bins 306 (6) book 2D-profile histogram with directory path and integer ID 307 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store 308 418 , ## the histogram integer identifier 309 'Energy Correction' , ## the histogram title 310 50 , ## number of X-bins 313 50 , ## number of Y-bins 317 Many other booking methods are available, 318 e.g. for the histograms with non-equidistant bins, see IHistogamSvc::book 323 raise RuntimeError(
'Unable to get valid HistogramService ')
325 return svc.bookProf(*args)
328 bookProf.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.bookProf) : \n\n' \
329 + iHistogramSvc.bookProf . __doc__
330 bookProf.__doc__ +=
'\n\n' +
'\thelp(IHistogramSvc::bookProf) : \n\n' \
331 + cpp.IHistogramSvc.bookProf . __doc__
340 The most trivial function to retrieve the histogram from Histogram Transient Store 341 The histogram is returned by reference to its AIDA-representation (if possible) 343 >>> h = getAsAIDA ( 'some/path/to/my/histogram' ) 348 raise RuntimeError(
'Unable to get valid HistogramService ')
350 return svc.getAsAIDA(path)
353 getAsAIDA.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.getAsAIDA) : \n\n' \
354 + iHistogramSvc.getAsAIDA . __doc__
355 getAsAIDA.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.retrieve) : \n\n' \
356 + iHistogramSvc.retrieve . __doc__
365 The most trivial function to retrieve the histogram from Histogram Transient Store 366 The histogram is returned by reference to its underlying native ROOT-representation (if possible) 368 >>> h = getAsROOT ( 'some/path/to/my/histogram' ) 373 raise RuntimeError(
'Unable to get valid HistogramService ')
375 return svc.getAsROOT(path)
378 getAsROOT.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.getAsROOT) : \n\n' \
379 + iHistogramSvc.getAsROOT . __doc__
392 The function which allows 'the smart fill' of 1D-histogram 396 The most trivial case, fill with the value 397 >>> fill ( histo , 1.0 ) 399 Fill from any iterable object (sequence) 400 >>> fill ( histo , [1,,2,3,4,5,10] ) 402 Fill from iterable object and apply the function: 403 >>> fill ( histo , [1,2,3,4,5] , math.sin ) 406 >>> fill ( histo , [1,2,3,4,5] , lambda x : x*x ) 409 >>> fill ( histo , [1,2,3,4,5] , fun = lambda x : x*x ) 411 Use internal attributes: 412 >>> tracks = evtSvc['Rec/Track/Best'] ## iterable container of tracks 413 >>> fill ( histo , tracks , lambda t : t.pt() ) 415 Apply the predicate: fill only even numbers: 416 >>> fill ( histo , [1,2,3,4,5,6,7] , lambda x : x , lambda y : y%2 ) 418 The same (omit the trivial function) : 419 >>> fill ( histo , [1,2,3,4,5,6,7] , cut = lambda y : y%2 ) 421 Apply the predicate: fill only pt for positively charged tracks: 422 >>> tracks = evtSvc['Rec/Track/Best'] 423 >>> fill ( histo , tracks , lambda t : t.pt() , lambda t : 0<t.charge() ) 426 >>> fill ( histo , tracks , 427 fun = lambda t : t.pt() , 428 cut = lambda t : 0<t.charge() ) 430 Ordinary functions are also fine: 431 >>> def myfun ( track ) : return sin( track.pt() + track.p() ) 432 >>> def mycut ( track ) : return track.p() > 100 * GeV 433 >>> fill ( histo , tracks , myfun , mycut ) 435 The 'data' could be the address in TES, in this case the object 436 is retrieved from TES and the method is applied to the objects, 438 >>> fill ( histo , ## the reference to the histogram 439 'Rec/Track/Best' , ## the location of objects in TES 440 lambda t : t.pt() ) ## function to be used for histogram fill 441 >>> fill ( histo , ## the reference to the histogram 442 'Rec/Track/Best' , ## the address of objects in TES 443 lambda t : t.pt() , ## the function to be used for histogram fill 444 lambda t : t.charge()>0 ) ## the criteria to select tracks 446 The arguments 'fun' and 'cut' could be strings, in this case 447 they are evaluated by python before usage. 448 This option could be often very useful. 453 if type(data) == str:
456 return fill(histo, data, fun, cut, **kwargs)
460 fun = eval(fun, globals())
464 cut = eval(cut, globals())
466 if not hasattr(data,
'__iter__'):
469 if not hasattr(histo,
'fill')
and hasattr(histo,
'Fill'):
470 setattr(histo,
'fill', getattr(histo,
'Fill'))
475 histo.fill(fun(item))
482 aida2root = cpp.Gaudi.Utils.Aida2ROOT.aida2root
492 >>> aida = ... ## get AIDA histogram 493 >>> root = aida.toROOT() ## convert it to ROOT 499 _to_root_.__doc__ += aida2root.__doc__
501 for t
in (cpp.AIDA.IHistogram3D, cpp.AIDA.IHistogram2D, cpp.AIDA.IHistogram1D,
502 cpp.AIDA.IProfile2D, cpp.AIDA.IProfile1D):
503 if not hasattr(t,
'Fill')
and hasattr(t,
'fill'):
504 setattr(t,
'Fill', getattr(t,
'fill'))
505 for attr
in (
'toROOT',
'toRoot',
'asROOT',
'asRoot',
'AsROOT',
'AsRoot'):
506 if not hasattr(t, attr):
507 setattr(t, attr, _to_root_)
509 cpp.AIDA.IHistogram3D. __repr__ =
lambda s: cpp.GaudiAlg.Print3D.toString(
511 cpp.AIDA.IHistogram3D.__str__ = cpp.AIDA.IHistogram3D.__repr__
513 HistoStats = cpp.Gaudi.Utils.HistoStats
521 Evaluate 'bin-by-bin' momentum of order 'order' around the value 'value' 525 >>> print h1.moment ( 5 ) 528 return HistoStats.moment(self, order, value)
537 Evaluate error for 'bin-by-bin' momentum of order 'order' around the value 'value' 541 >>> print h1.momentErr ( 5 ) 544 return HistoStats.momentErr(self, order)
553 Evaluate 'bin-by-bin' central momentum (around mean value) 557 >>> print h1.centralMoment ( 5 ) 560 return HistoStats.centralMoment(self, order)
569 Evaluate error for 'bin-by-bin' central momentum (around mean value) 573 >>> print h1.centralMomentErr ( 5 ) 576 return HistoStats.centralMomentErr(self, order)
585 Evaluate 'bin-by-bin' skewness for 1D AIDA histogram 588 >>> print h1.skewness() 591 return HistoStats.skewness(self)
600 Evaluate error for 'bin-by-bin' skewness 603 >>> print h1.skewnessErr() 606 return HistoStats.skewnessErr(self)
615 Evaluate 'bin-by-bin' kurtosis 618 >>> print h1.kurtosis () 621 return HistoStats.kurtosis(self)
630 Evaluate error for 'bin-by-bin' kurtotis for 1D AIDA histogram 633 >>> print h1.kurtotisErr() 636 return HistoStats.kurtosisErr(self)
644 Number of equivalent entries 646 return HistoStats.nEff(self)
654 Evaluate the MEAN value 656 return HistoStats.mean(self)
664 Evaluate the error for MEAN estimate 666 return HistoStats.meanErr(self)
674 Evaluate the RMS for AIDA histogram 676 return HistoStats.rms(self)
684 Evaluate the error for RMS estimate 686 return HistoStats.rmsErr(self)
694 Get an error in the sum bin height ('in-range integral') 696 return HistoStats.sumBinHeightErr(self)
703 """ Get an error in the sum bin height ('in-range integral') """ 704 return HistoStats.sumAllBinHeightErr(self)
712 The fraction of overflow entries (useful for shape comparison) 714 return HistoStats.overflowEntriesFrac(self)
722 The error for fraction of overflow entries (useful for shape comparison) 724 return HistoStats.overflowEntriesFracErr(self)
732 The fraction of underflow entries (useful for shape comparison) 734 return HistoStats.underflowEntriesFrac(self)
742 The error for fraction of underflow entries (useful for shape comparison) 744 return HistoStats.underflowEntriesFracErr(self)
752 The fraction of overflow integral (useful for shape comparison) 754 return HistoStats.overflowIntegralFrac(self)
762 The error for fraction of overflow integral (useful for shape comparison) 764 return HistoStats.overflowIntegralFracErr(self)
772 The fraction of underflow integral (useful for shape comparison) 774 return HistoStats.underflowIntegralFrac(self)
782 The error for fraction of underflow integral (useful for shape comparison) 784 return HistoStats.underflowIntegralFracErr(self)
795 Get number of entries in histogram up to the certain bin (not-included) 797 attention: underflow bin is included! 800 >>> print h1.nEntries ( 10 ) 802 Get number of entries in histogram form the certain 803 minimal bin up to the certain maximal bin (not-included) 806 >>> print h1.nEntries ( 10 , 15 ) 809 if i2 < i1
or i2 < 0:
810 return HistoStats.nEntries(self, i1)
811 return HistoStats.nEntries(self, i1, i2)
819 Get the fraction of entries in histogram up to the certain bin (not-included) 821 attention: underflow bin is included! 824 >>> print h1.nEntriesFrac ( 10 ) 826 Get the fraction of entries in histogram form the certain 827 minimal bin up to the certain maximal bin (not-included) 830 >>> print h1.nEntriesFrac ( 10 , 15 ) 833 if i2 < i1
or i2 < 0:
834 return HistoStats.nEntriesFrac(self, i1)
835 return HistoStats.nEntriesFrac(self, i1, i2)
843 Get error for fraction of entries in histogram up to the certain bin (not-included) 845 attention: underflow bin is included! 848 >>> print h1.nEntriesFracErr( 10 ) 850 Get error fraction of entries in histogram form the certain 851 minimal bin up to the certain maximal bin (not-included) 854 >>> print h1.nEntriesFracErr ( 10 , 15 ) 857 if i2 < i1
or i2 < 0:
858 return HistoStats.nEntriesFracErr(self, i1)
859 return HistoStats.nEntriesFracErr(self, i1, i2)
863 i1DH = cpp.AIDA.IHistogram1D
865 if not hasattr(i1DH,
'moment'):
866 i1DH.moment = _moment_
867 if not hasattr(i1DH,
'momentErr'):
868 i1DH.momentErr = _momentErr_
869 if not hasattr(i1DH,
'centralMoment'):
870 i1DH.centralMoment = _centralMoment_
871 if not hasattr(i1DH,
'momentMomentErr'):
872 i1DH.centralMomentErr = _centralMomentErr_
873 if not hasattr(i1DH,
'nEff'):
875 if not hasattr(i1DH,
'mean'):
877 if not hasattr(i1DH,
'meanErr'):
878 i1DH.meanErr = _meanErr_
879 if not hasattr(i1DH,
'rms'):
881 if not hasattr(i1DH,
'rmsErr'):
882 i1DH.rmsErr = _rmsErr_
883 if not hasattr(i1DH,
'skewness'):
884 i1DH.skewness = _skewness_
885 if not hasattr(i1DH,
'skewnessErr'):
886 i1DH.skewnessErr = _skewnessErr_
887 if not hasattr(i1DH,
'kurtosis'):
888 i1DH.kurtosis = _kurtosis_
889 if not hasattr(i1DH,
'kurtosisErr'):
890 i1DH.kurtosisErr = _kurtosisErr_
892 if not hasattr(i1DH,
'overflowEntriesFrac'):
893 i1DH.overflowEntriesFrac = _overflowEntriesFrac_
894 if not hasattr(i1DH,
'overflowEntriesFracErr'):
895 i1DH.overflowEntriesFracErr = _overflowEntriesFracErr_
896 if not hasattr(i1DH,
'underflowEntriesFrac'):
897 i1DH.underflowEntriesFrac = _underflowEntriesFrac_
898 if not hasattr(i1DH,
'underflowEntriesFracErr'):
899 i1DH.underflowEntriesFracErr = _underflowEntriesFracErr_
901 if not hasattr(i1DH,
'overflowIntegralFrac'):
902 i1DH.overflowIntegralFrac = _overflowIntegralFrac_
903 if not hasattr(i1DH,
'overflowIntegralFracErr'):
904 i1DH.overflowIntegralFracErr = _overflowIntegralFracErr_
905 if not hasattr(i1DH,
'underflowIntegralFrac'):
906 i1DH.underflowIntegralFrac = _underflowIntegralFrac_
907 if not hasattr(i1DH,
'underflowIntegralFracErr'):
908 i1DH.underflowIntegralFracErr = _underflowIntegralFracErr_
910 if not hasattr(i1DH,
'nEntries'):
911 i1DH.nEntries = _nEntries_
912 if not hasattr(i1DH,
'nEntriesFrac'):
913 i1DH.nEntriesFrac = _nEntriesFrac_
914 if not hasattr(i1DH,
'nEntriesFracErr'):
915 i1DH.nEntriesFracErr = _nEntriesFracErr_
922 Get the path in THS for the given AIDA object: 925 >>> print aida.path() 928 return cpp.Gaudi.Utils.Histos.path(self)
931 iBH = cpp.AIDA.IBaseHistogram
932 if not hasattr(iBH,
'path'):
934 if not hasattr(iBH,
'TESpath'):
936 if not hasattr(iBH,
'location'):
937 iBH.location = _path_
944 Dump the histogram/profile in text format (a'la HBOOK) 947 >>> print dumpHisto ( histo ) 949 >>> print histo.dump() 950 >>> print histo.dump( 20 , 20 ) 951 >>> print histo.dump( 20 , 20 , True ) 956 return cpp.Gaudi.Utils.Histos.histoDump(histo, *args)
959 __dumpHisto__.__doc__ =
'\n' + cpp.Gaudi.Utils.Histos.histoDump.__doc__
963 histoDump = __dumpHisto__
964 dumpHisto = __dumpHisto__
966 for t
in (cpp.AIDA.IHistogram1D, cpp.AIDA.IProfile1D, ROOT.TH1D, ROOT.TH1F,
967 ROOT.TH1, ROOT.TProfile):
968 for method
in (
'dump',
'dumpHisto',
'dumpAsText'):
969 if not hasattr(t, method):
970 setattr(t, method, __dumpHisto__)
977 Class to write histograms to a ROOT file. 978 hFile = HistoFile("myFile.root") 984 hFile.save(myHisto0, myHisto1, myHisto2) 985 histoList = [h0, h1, h2, h3] 986 hFile.save(histoList) 990 __author__ =
"Juan Palacios juan.palacios@nikhef.nl" 993 self.
file = ROOT.TFile(fileName,
"RECREATE")
994 from GaudiPython
import gbl
997 gbl.AIDA.IHistogram1D, gbl.AIDA.IHistogram2D,
998 gbl.AIDA.IHistogram3D, gbl.AIDA.IProfile1D, gbl.AIDA.IProfile2D,
1003 histoType =
type(histo)
1011 This function stores histograms on the file for future saving. 1012 It takes an arbitrary number of AIDA or ROOT histograms or 1016 if type(args[0]) == list:
1031 if '__main__' == __name__:
1036 print(sys.modules[__name__].__dict__[o].__doc__)
def _underflowIntegralFrac_(self)
def _nEntriesFracErr_(self, i1, i2=-10000000)
def _sumBinHeightErr_(self)
def __init__(self, fileName)
def _sumAllBinHeightErr_(self)
def getAsROOT(path, **kwargs)
def _overflowIntegralFracErr_(self)
def _getHistoSvc(**kwargs)
def _overflowEntriesFrac_(self)
def book(*args, **kwargs)
def _moment_(self, order, value=0)
def _nEntriesFrac_(self, i1, i2=-10000000)
def getAsAIDA(path, **kwargs)
def _underflowEntriesFrac_(self)
def _overflowIntegralFrac_(self)
def _overflowEntriesFracErr_(self)
def __dumpHisto__(histo, *args)
def _underflowEntriesFracErr_(self)
def _centralMomentErr_(self, order)
def bookProf(*args, **kwargs)
def _nEntries_(self, i1, i2=-10000000)
def fill(histo, data, fun=lambda x:x, cut=lambda x:True, **kwargs)
def _momentErr_(self, order)
def _centralMoment_(self, order)
def _underflowIntegralFracErr_(self)
def __convertibleType(self, histo)