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)
27 __author__ =
"Vanya BELYAEV ibelyaev@physics.syr.edu"
56 Helper private auxiliary function to get Application Manager
58 gaudi = kwargs.get (
'gaudi' ,
None )
59 if not gaudi : gaudi =
AppMgr()
60 if not gaudi :
raise RuntimeError,
'Unable to get valid ApplicationMgr'
62 state = gaudi._isvc.FSMState()
63 if state < cpp.Gaudi.StateMachine.CONFIGURED : gaudi.config ()
64 state = gaudi._isvc.FSMState()
65 if state < cpp.Gaudi.StateMachine.INITIALIZED : gaudi.initialize ()
73 Helper private auxiliary function to get iHistogramSvs
75 svc = kwargs.get (
'service' ,
None )
76 if not svc : svc = kwargs.get (
'svc' ,
None )
78 gaudi = _getAppMgr ( **kwargs )
79 return gaudi.histsvc ()
85 Helper private auxiliary function to get iDataSvs
87 svc = kwargs.get (
'service' ,
None )
88 if not svc : svc = kwargs.get (
'svc' ,
None )
90 gaudi = _getAppMgr ( **kwargs )
95 def book ( *args , **kwargs ) :
97 The trivial function to book the various 1D,2D&3D-histograms
99 (1) book the trivial 1D histogram with full path
101 >>> h1D = book ( 'path/to/my/histo' , ## path in Histogram Transient Store
102 'cosine of decay angle ' , ## histogram title
103 100 , ## number of bins
107 (2) book the trivial 1D histogram with directory path and string ID :
109 >>> h1D = book ( 'path/to/directory' , ## the path to directory in HTS
110 'H1' , ## string histogram identifier
111 'cosine of decay angle ' , ## histogram title
112 100 , ## number of bins
116 (3) book the trivial 1D histogram with directory path and integer ID :
118 >>> h1D = book ( 'path/to/directory' , ## the path to directory in HTS
119 124 , ## integer histogram identifier
120 'cosine of decay angle ' , ## histogram title
121 100 , ## number of bins
125 (4) book the trivial 2D histogram with full path
127 >>> h1D = book ( 'path/to/my/histo' , ## path in Histogram Transient Store
128 'm12**2 versus m13**2' , ## histogram title
129 50 , ## number of X-bins
132 50 , ## number of Y-bins
136 (5) book the trivial 2D histogram with directory path and literal ID
138 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store
139 'Dalitz1' , ## literal histogram identifier
140 'm12**2 versus m13**2' , ## histogram title
141 50 , ## number of X-bins
144 50 , ## number of Y-bins
148 (6) book the trivial 2D histogram with directory path and integer ID
150 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store
151 854 , ## integer histogram identifier
152 'm12**2 versus m13**2' , ## histogram title
153 50 , ## number of X-bins
156 50 , ## number of Y-bins
160 (7) book the trivial 3D histogram with full path
162 >>> h1D = book ( 'path/to/my/histo' , ## path in Histogram Transient Store
163 'x vs y vs z ' , ## histogram title
164 10 , ## number of X-bins
167 10 , ## number of Y-bins
170 10 , ## number of Z-bins
174 (8) book the trivial 3D histogram with directory path and literal ID
176 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store
177 'xyz' , ## literal histogram identifier
178 'x vs y vs z' , ## histogram title
179 10 , ## number of X-bins
182 10 , ## number of Y-bins
185 10 , ## number of Z-bins
189 (9) book the trivial 3D histogram with directory path and integer ID
191 >>> h1D = book ( 'path/to/directory' , ## path in Histogram Transient Store
192 888 , ## integer histogram identifier
193 'x vs y vs z' , ## histogram title
194 10 , ## number of X-bins
197 10 , ## number of Y-bins
200 10 , ## number of Z-bins
204 Many other booking methods are available,
205 e.g. for the histograms with non-equidistant bins, see IHistogamSvc::book
208 if useROOT
or kwargs.get(
'useROOT' ,
False )
or not kwargs.get(
'useAIDA' ,
True ) :
209 from ROOT
import TH1D
213 if not str
is type(a1) :
216 return TH1D ( a0+a1 , a2 , *args[3:] )
218 return TH1D ( a0 , a1 , *args[2:] )
220 svc = _getHistoSvc ( **kwargs )
221 if not svc :
raise RuntimeError,
'Unable to get valid HistogramService '
223 return svc.book(*args)
225 book.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.book) : \n\n' \
226 + iHistogramSvc.book . __doc__
227 book.__doc__ +=
'\n\n' +
'\thelp(IHistogramSvc::book) : \n\n' \
228 + cpp.IHistogramSvc.book . __doc__
235 The trivial function to book 1D&2D profile histograms:
237 (1) book 1D-profile histogram with full path in Histogram Transient Store:
238 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store
239 'Energy Correction' , ## the histogram title
240 100 , ## number of X-bins
244 (2) book 1D-profile histogram with directory path and literal ID
245 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store
246 'Calibration' , ## the histogram literal identifier
247 'Energy Correction' , ## the histogram title
248 100 , ## number of X-bins
252 (3) book 1D-profile histogram with directory path and integer ID
253 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store
254 418 , ## the histogram integer identifier
255 'Energy Correction' , ## the histogram title
256 100 , ## number of X-bins
260 (4) book 2D-profile histogram with full path in Histogram Transient Store:
261 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store
262 'Energy Correction' , ## the histogram title
263 50 , ## number of X-bins
266 50 , ## number of Y-bins
270 (5) book 2D-profile histogram with directory path and literal ID
271 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store
272 'Calibration' , ## the histogram literal identifier
273 'Energy Correction' , ## the histogram title
274 50 , ## number of X-bins
277 50 , ## number of Y-bins
281 (6) book 2D-profile histogram with directory path and integer ID
282 >>> histo = bookProf ( 'path/to/my/profile' , ## path in Histogram Transient Store
283 418 , ## the histogram integer identifier
284 'Energy Correction' , ## the histogram title
285 50 , ## number of X-bins
288 50 , ## number of Y-bins
292 Many other booking methods are available,
293 e.g. for the histograms with non-equidistant bins, see IHistogamSvc::book
296 svc = _getHistoSvc ( **kwargs )
297 if not svc :
raise RuntimeError,
'Unable to get valid HistogramService '
299 return svc.bookProf(*args)
302 bookProf.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.bookProf) : \n\n' \
303 + iHistogramSvc.bookProf . __doc__
304 bookProf.__doc__ +=
'\n\n' +
'\thelp(IHistogramSvc::bookProf) : \n\n' \
305 + cpp.IHistogramSvc.bookProf . __doc__
312 The most trivial function to retrieve the histogram from Histogram Transient Store
313 The histogram is returned by reference to its AIDA-representation (if possible)
315 >>> h = getAsAIDA ( 'some/path/to/my/histogram' )
318 svc = _getHistoSvc ( **kwargs )
319 if not svc :
raise RuntimeError,
'Unable to get valid HistogramService '
321 return svc.getAsAIDA( path )
323 getAsAIDA.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.getAsAIDA) : \n\n' \
324 + iHistogramSvc.getAsAIDA . __doc__
325 getAsAIDA.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.retrieve) : \n\n' \
326 + iHistogramSvc.retrieve . __doc__
333 The most trivial function to retrieve the histogram from Histogram Transient Store
334 The histogram is returned by reference to its underlying native ROOT-representation (if possible)
336 >>> h = getAsROOT ( 'some/path/to/my/histogram' )
339 svc = _getHistoSvc ( **kwargs )
340 if not svc :
raise RuntimeError,
'Unable to get valid HistogramService '
342 return svc.getAsROOT( path )
344 getAsROOT.__doc__ +=
'\n\n' +
'\thelp(iHistogramSvc.getAsROOT) : \n\n' \
345 + iHistogramSvc.getAsROOT . __doc__
353 cut =
lambda x :
True ,
357 The function which allows 'the smart fill' of 1D-histogram
361 The most trivial case, fill with the value
362 >>> fill ( histo , 1.0 )
364 Fill from any iterable object (sequence)
365 >>> fill ( histo , [1,,2,3,4,5,10] )
367 Fill from iterable object and apply the function:
368 >>> fill ( histo , [1,2,3,4,5] , math.sin )
371 >>> fill ( histo , [1,2,3,4,5] , lambda x : x*x )
374 >>> fill ( histo , [1,2,3,4,5] , fun = lambda x : x*x )
376 Use internal attributes:
377 >>> tracks = evtSvc['Rec/Track/Best'] ## iterable container of tracks
378 >>> fill ( histo , tracks , lambda t : t.pt() )
380 Apply the predicate: fill only even numbers:
381 >>> fill ( histo , [1,2,3,4,5,6,7] , lambda x : x , lambda y : y%2 )
383 The same (omit the trivial function) :
384 >>> fill ( histo , [1,2,3,4,5,6,7] , cut = lambda y : y%2 )
386 Apply the predicate: fill only pt for positively charged tracks:
387 >>> tracks = evtSvc['Rec/Track/Best']
388 >>> fill ( histo , tracks , lambda t : t.pt() , lambda t : 0<t.charge() )
391 >>> fill ( histo , tracks ,
392 fun = lambda t : t.pt() ,
393 cut = lambda t : 0<t.charge() )
395 Ordinary functions are also fine:
396 >>> def myfun ( track ) : return sin( track.pt() + track.p() )
397 >>> def mycut ( track ) : return track.p() > 100 * GeV
398 >>> fill ( histo , tracks , myfun , mycut )
400 The 'data' could be the address in TES, in this case the object
401 is retrieved from TES and the method is applied to the objects,
403 >>> fill ( histo , ## the reference to the histogram
404 'Rec/Track/Best' , ## the location of objects in TES
405 lambda t : t.pt() ) ## function to be used for histogram fill
406 >>> fill ( histo , ## the reference to the histogram
407 'Rec/Track/Best' , ## the address of objects in TES
408 lambda t : t.pt() , ## the function to be used for histogram fill
409 lambda t : t.charge()>0 ) ## the criteria to select tracks
411 The arguments 'fun' and 'cut' could be strings, in this case
412 they are evaluated by python before usage.
413 This option could be often very useful.
418 if type ( data ) == str :
419 svc = _getEvtSvc ( **kwargs )
421 return fill ( histo , data , fun , cut , **kwargs )
424 if type ( fun ) == str : fun = eval ( fun , globals() )
427 if type ( cut ) == str : cut = eval ( cut , globals() )
429 if not hasattr ( data ,
'__iter__' ) : data = [ data ]
431 if not hasattr ( histo ,
'fill' )
and hasattr ( histo ,
'Fill' ) :
432 setattr ( histo ,
'fill' , getattr ( histo ,
'Fill' ) )
435 if not cut ( item ) :
continue
436 histo.fill ( fun ( item ) )
442 aida2root = cpp.Gaudi.Utils.Aida2ROOT.aida2root
449 >>> aida = ... ## get AIDA histogram
450 >>> root = aida.toROOT() ## convert it to ROOT
453 return aida2root ( self )
455 _to_root_ . __doc__ += aida2root . __doc__
457 for t
in ( cpp.AIDA.IHistogram3D ,
458 cpp.AIDA.IHistogram2D ,
459 cpp.AIDA.IHistogram1D ,
460 cpp.AIDA.IProfile2D ,
461 cpp.AIDA.IProfile1D ) :
462 if not hasattr ( t ,
'Fill' )
and hasattr ( t ,
'fill' ) :
463 setattr ( t ,
'Fill' , getattr ( t ,
'fill' ) )
464 for attr
in (
'toROOT' ,
'toRoot' ,
465 'asROOT' ,
'asRoot' ,
466 'AsROOT' ,
'AsRoot' ) :
467 if not hasattr ( t , attr ) : setattr ( t , attr , _to_root_ )
469 cpp.AIDA.IHistogram3D. __repr__ =
lambda s : cpp.GaudiAlg.Print3D.toString( s ,
HID( s.title() ) )
470 cpp.AIDA.IHistogram3D. __str__ = cpp.AIDA.IHistogram3D. __repr__
473 HistoStats = cpp.Gaudi.Utils.HistoStats
479 Evaluate 'bin-by-bin' momentum of order 'order' around the value 'value'
483 >>> print h1.moment ( 5 )
486 return HistoStats.moment ( self , order , value )
492 Evaluate error for 'bin-by-bin' momentum of order 'order' around the value 'value'
496 >>> print h1.momentErr ( 5 )
499 return HistoStats.momentErr ( self , order )
504 Evaluate 'bin-by-bin' central momentum (around mean value)
508 >>> print h1.centralMoment ( 5 )
511 return HistoStats.centralMoment ( self , order )
517 Evaluate error for 'bin-by-bin' central momentum (around mean value)
521 >>> print h1.centralMomentErr ( 5 )
524 return HistoStats.centralMomentErr ( self , order )
530 Evaluate 'bin-by-bin' skewness for 1D AIDA histogram
533 >>> print h1.skewness()
536 return HistoStats.skewness ( self )
542 Evaluate error for 'bin-by-bin' skewness
545 >>> print h1.skewnessErr()
548 return HistoStats.skewnessErr ( self )
554 Evaluate 'bin-by-bin' kurtosis
557 >>> print h1.kurtosis ()
560 return HistoStats.kurtosis ( self )
566 Evaluate error for 'bin-by-bin' kurtotis for 1D AIDA histogram
569 >>> print h1.kurtotisErr()
572 return HistoStats.kurtosisErr ( self )
577 Number of equivalent entries
579 return HistoStats.nEff ( self )
583 Evaluate the MEAN value
585 return HistoStats.mean ( self )
589 Evaluate the error for MEAN estimate
591 return HistoStats.meanErr ( self )
596 Evaluate the RMS for AIDA histogram
598 return HistoStats.rms ( self )
602 Evaluate the error for RMS estimate
604 return HistoStats.rmsErr ( self )
609 Get an error in the sum bin height ('in-range integral')
611 return HistoStats.sumBinHeightErr ( self )
615 """ Get an error in the sum bin height ('in-range integral') """
616 return HistoStats.sumAllBinHeightErr ( self )
621 The fraction of overflow entries (useful for shape comparison)
623 return HistoStats.overflowEntriesFrac ( self )
627 The error for fraction of overflow entries (useful for shape comparison)
629 return HistoStats.overflowEntriesFracErr ( self )
633 The fraction of underflow entries (useful for shape comparison)
635 return HistoStats.underflowEntriesFrac ( self )
639 The error for fraction of underflow entries (useful for shape comparison)
641 return HistoStats.underflowEntriesFracErr ( self )
646 The fraction of overflow integral (useful for shape comparison)
648 return HistoStats.overflowIntegralFrac ( self )
652 The error for fraction of overflow integral (useful for shape comparison)
654 return HistoStats.overflowIntegralFracErr ( self )
658 The fraction of underflow integral (useful for shape comparison)
660 return HistoStats.underflowIntegralFrac ( self )
664 The error for fraction of underflow integral (useful for shape comparison)
666 return HistoStats.underflowIntegralFracErr ( self )
674 Get number of entries in histogram up to the certain bin (not-included)
676 attention: underflow bin is included!
679 >>> print h1.nEntries ( 10 )
681 Get number of entries in histogram form the certain
682 minimal bin up to the certain maximal bin (not-included)
685 >>> print h1.nEntries ( 10 , 15 )
688 if i2 < i1
or i2 < 0 :
return HistoStats.nEntries ( self , i1 )
689 return HistoStats.nEntries ( self , i1 , i2 )
693 Get the fraction of entries in histogram up to the certain bin (not-included)
695 attention: underflow bin is included!
698 >>> print h1.nEntriesFrac ( 10 )
700 Get the fraction of entries in histogram form the certain
701 minimal bin up to the certain maximal bin (not-included)
704 >>> print h1.nEntriesFrac ( 10 , 15 )
707 if i2 < i1
or i2 < 0 :
return HistoStats.nEntriesFrac ( self , i1 )
708 return HistoStats.nEntriesFrac ( self , i1 , i2 )
712 Get error for fraction of entries in histogram up to the certain bin (not-included)
714 attention: underflow bin is included!
717 >>> print h1.nEntriesFracErr( 10 )
719 Get error fraction of entries in histogram form the certain
720 minimal bin up to the certain maximal bin (not-included)
723 >>> print h1.nEntriesFracErr ( 10 , 15 )
726 if i2 < i1
or i2 < 0 :
return HistoStats.nEntriesFracErr ( self , i1 )
727 return HistoStats.nEntriesFracErr ( self , i1 , i2 )
730 i1DH = cpp.AIDA.IHistogram1D
732 if not hasattr ( i1DH ,
'moment' ) : i1DH.moment = _moment_
733 if not hasattr ( i1DH ,
'momentErr' ) : i1DH.momentErr = _momentErr_
734 if not hasattr ( i1DH ,
'centralMoment' ) : i1DH.centralMoment = _centralMoment_
735 if not hasattr ( i1DH ,
'momentMomentErr' ) : i1DH.centralMomentErr = _centralMomentErr_
736 if not hasattr ( i1DH ,
'nEff' ) : i1DH.nEff = _nEff_
737 if not hasattr ( i1DH ,
'mean' ) : i1DH.mean = _mean_
738 if not hasattr ( i1DH ,
'meanErr' ) : i1DH.meanErr = _meanErr_
739 if not hasattr ( i1DH ,
'rms' ) : i1DH.rms = _rms_
740 if not hasattr ( i1DH ,
'rmsErr' ) : i1DH.rmsErr = _rmsErr_
741 if not hasattr ( i1DH ,
'skewness' ) : i1DH.skewness = _skewness_
742 if not hasattr ( i1DH ,
'skewnessErr' ) : i1DH.skewnessErr = _skewnessErr_
743 if not hasattr ( i1DH ,
'kurtosis' ) : i1DH.kurtosis = _kurtosis_
744 if not hasattr ( i1DH ,
'kurtosisErr' ) : i1DH.kurtosisErr = _kurtosisErr_
746 if not hasattr ( i1DH ,
'overflowEntriesFrac' ) : i1DH.overflowEntriesFrac = _overflowEntriesFrac_
747 if not hasattr ( i1DH ,
'overflowEntriesFracErr' ) : i1DH.overflowEntriesFracErr = _overflowEntriesFracErr_
748 if not hasattr ( i1DH ,
'underflowEntriesFrac' ) : i1DH.underflowEntriesFrac = _underflowEntriesFrac_
749 if not hasattr ( i1DH ,
'underflowEntriesFracErr' ) : i1DH.underflowEntriesFracErr = _underflowEntriesFracErr_
751 if not hasattr ( i1DH ,
'overflowIntegralFrac' ) : i1DH.overflowIntegralFrac = _overflowIntegralFrac_
752 if not hasattr ( i1DH ,
'overflowIntegralFracErr' ) : i1DH.overflowIntegralFracErr = _overflowIntegralFracErr_
753 if not hasattr ( i1DH ,
'underflowIntegralFrac' ) : i1DH.underflowIntegralFrac = _underflowIntegralFrac_
754 if not hasattr ( i1DH ,
'underflowIntegralFracErr' ) : i1DH.underflowIntegralFracErr = _underflowIntegralFracErr_
756 if not hasattr ( i1DH ,
'nEntries' ) : i1DH.nEntries = _nEntries_
757 if not hasattr ( i1DH ,
'nEntriesFrac' ) : i1DH.nEntriesFrac = _nEntriesFrac_
758 if not hasattr ( i1DH ,
'nEntriesFracErr' ) : i1DH.nEntriesFracErr = _nEntriesFracErr_
763 Get the path in THS for the given AIDA object:
766 >>> print aida.path()
769 return cpp.Gaudi.Utils.Histos.path ( self )
771 iBH = cpp.AIDA.IBaseHistogram
772 if not hasattr ( iBH ,
'path' ) : iBH.path = _path_
773 if not hasattr ( iBH ,
'TESpath' ) : iBH.TESpath = _path_
774 if not hasattr ( iBH ,
'location' ) : iBH.location = _path_
781 Dump the histogram/profile in text format (a'la HBOOK)
784 >>> print dumpHisto ( histo )
786 >>> print histo.dump()
787 >>> print histo.dump( 20 , 20 )
788 >>> print histo.dump( 20 , 20 , True )
793 return cpp.Gaudi.Utils.Histos.histoDump ( histo , *args )
795 __dumpHisto__ .__doc__ =
'\n' + cpp.Gaudi.Utils.Histos.histoDump . __doc__
799 histoDump = __dumpHisto__
800 dumpHisto = __dumpHisto__
802 for t
in ( cpp.AIDA.IHistogram1D ,
803 cpp.AIDA.IProfile1D ,
808 for method
in (
'dump' ,
811 if not hasattr ( t , method ) : setattr ( t , method , __dumpHisto__ )
816 Class to write histograms to a ROOT file.
817 hFile = HistoFile("myFile.root")
823 hFile.save(myHisto0, myHisto1, myHisto2)
824 histoList = [h0, h1, h2, h3]
825 hFile.save(histoList)
829 __author__ =
"Juan Palacios juan.palacios@nikhef.nl"
832 self.
file = ROOT.TFile(fileName,
"RECREATE")
833 from GaudiPython
import gbl
836 gbl.AIDA.IHistogram2D,
837 gbl.AIDA.IHistogram3D,
840 gbl.AIDA.IHistogram ]
843 histoType =
type(histo)
845 if histoType == t :
return True
850 This function stores histograms on the file for future saving.
851 It takes an arbitrary number of AIDA or ROOT histograms or
855 if type(args[0])==list :
868 if '__main__' == __name__ :
873 print sys.modules[__name__].__dict__[o].__doc__
def _skewness_
Evaluate 'bin-by-bin' skewness for 1D histogram.
def _underflowIntegralFracErr_
def getAsAIDA
The most trivial function to retrieve the histogram from Histogram Transient Store.
def __dumpHisto__
=============================================================================
def book
The trivial function to book the various 1D,2D&3D-histograms.
def _getHistoSvc
Helper private auxiliary function to get iHistogramSvs.
def _to_root_
Convert AIDA to ROOT.
def _moment_
Evaluate 'bin-by-bin' momentum of certain order around the value.
def getAsROOT
The most trivial function to retrieve the histogram from Histogram Transient Store.
def _overflowIntegralFrac_
def _centralMomentErr_
Evaluate error in 'bin-by-bin' momentum of certain order around the value.
def _getEvtSvc
Helper private auxiliary function to get iDataSvs.
def _skewnessErr_
Evaluate error for 'bin-by-bin' skewness for 1D histogram.
def _underflowEntriesFracErr_
def _overflowEntriesFrac_
def fill
The function which allows 'the smart fill' of 1D-histogram.
def _kurtosisErr_
Evaluate error for 'bin-by-bin' kurtosis for 1D histogram.
def _overflowEntriesFracErr_
def _path_
============================================================================
def _overflowIntegralFracErr_
def _getAppMgr
Helper private auxiliary function to get Application Manager.
def _underflowIntegralFrac_
def _nEntries_
get number of entries in histogram up to the certain bin (not-included) get number of entries in hist...
def _underflowEntriesFrac_
def bookProf
The trivial function to book 1D&2D profile histograms:
def _centralMoment_
Evaluate 'bin-by-bin' central momentum (around mean value)
def _kurtosis_
Evaluate 'bin-by-bin' kurtosis for 1D histogram.
def _momentErr_
Evaluate error in 'bin-by-bin' momentum of certain order around the value.