Loading [MathJax]/extensions/tex2jax.js
The Gaudi Framework  v31r0 (aeb156f0)
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Gaudi::Utils::Histos Namespace Reference

Collection of useful utilities for manipulations with AIDA hisgograms. More...

Namespaces

 Formats
 The 15 format fields are predefined now:
 

Classes

class  HistoStrings
 Helper class to produce "good" dictionaries. More...
 
class  Table
 Simple class for the customizeble printout of the histogram tables. More...
 

Typedefs

typedef std::vector< std::stringLabels
 Typedef for a list of labels. More...
 
typedef std::pair< unsigned, std::stringBinLabel
 Typedef for a bin number and its associated label. More...
 
typedef std::vector< BinLabelBinLabels
 Typedef for a list of bin numbers and their associated label. More...
 

Functions

GAUDI_API void fill (AIDA::IHistogram1D *histo, const double value, const double weight=1.0)
 simple function to fill AIDA::IHistogram1D objects More...
 
GAUDI_API void fill (AIDA::IHistogram2D *histo, const double valueX, const double valueY, const double weight=1.0)
 simple function to fill AIDA::IHistogram2D objects More...
 
GAUDI_API void fill (AIDA::IHistogram3D *histo, const double valueX, const double valueY, const double valueZ, const double weight=1.0)
 simple function to fill AIDA::IHistogram3D objects More...
 
GAUDI_API void fill (AIDA::IProfile1D *histo, const double valueX, const double valueY, const double weight=1.0)
 simple function to fill AIDA::IProfile1D objects More...
 
GAUDI_API void fill (AIDA::IProfile2D *histo, const double valueX, const double valueY, const double valueZ, const double weight=1.0)
 simple function to fill AIDA::IProfile2D objects More...
 
GAUDI_API std::string htitle (const AIDA::IBaseHistogram *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IHistogram *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IHistogram1D *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IHistogram2D *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IHistogram3D *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IProfile *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IProfile1D *histo, const std::string &title="")
 get the title More...
 
GAUDI_API std::string htitle (const AIDA::IProfile2D *histo, const std::string &title="")
 get the title More...
 
GAUDI_API AIDA::IBaseHistogram * toBase (AIDA::IHistogram1D *histo)
 
GAUDI_API AIDA::IBaseHistogram * toBase (AIDA::IHistogram2D *histo)
 
GAUDI_API AIDA::IBaseHistogram * toBase (AIDA::IHistogram3D *histo)
 
GAUDI_API AIDA::IBaseHistogram * toBase (AIDA::IProfile1D *histo)
 
GAUDI_API AIDA::IBaseHistogram * toBase (AIDA::IProfile2D *histo)
 
GAUDI_API std::ostreamhistoDump_ (const AIDA::IHistogram1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
 dump the text representation of the histogram More...
 
GAUDI_API std::string histoDump (const AIDA::IHistogram1D *histo, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
 dump the text representation of the histogram More...
 
GAUDI_API std::ostreamhistoDump_ (const AIDA::IProfile1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool spread=true)
 dump the text representation of 1D-profile More...
 
GAUDI_API std::string histoDump (const AIDA::IProfile1D *histo, const std::size_t width=80, const std::size_t height=50, const bool spread=true)
 dump the text representation of the 1D-profile More...
 
GAUDI_API std::ostreamhistoDump_ (const TProfile *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50)
 dump the text representation of the Profile More...
 
GAUDI_API std::string histoDump (const TProfile *histo, const std::size_t width=80, const std::size_t height=50)
 dump the text representation of the histogram More...
 
GAUDI_API std::ostreamhistoDump_ (const TH1 *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
 dump the text representation of the histogram More...
 
GAUDI_API std::string histoDump (const TH1 *histo, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
 dump the text representation of the histogram More...
 
GAUDI_API bool setBinLabels (AIDA::IHistogram1D *hist, const Labels &labels)
 Set the Bin labels for a given 1D histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IHistogram1D *hist, const BinLabels &labels)
 Set the Bin labels for a given 1D histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IProfile1D *hist, const Labels &labels)
 Set the Bin labels for a given 1D profile histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IProfile1D *hist, const BinLabels &labels)
 Set the Bin labels for a given 1D profile histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IHistogram2D *hist, const Labels &xlabels, const Labels &ylabels)
 Set the Bin labels for a given 2D histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IHistogram2D *hist, const BinLabels &xlabels, const BinLabels &ylabels)
 Set the Bin labels for a given 2D histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IProfile2D *hist, const Labels &xlabels, const Labels &ylabels)
 Set the Bin labels for a given 2D profile histogram. More...
 
GAUDI_API bool setBinLabels (AIDA::IProfile2D *hist, const BinLabels &xlabels, const BinLabels &ylabels)
 Set the Bin labels for a given 2D profile histogram. More...
 
GAUDI_API bool setAxisLabels (AIDA::IHistogram1D *hist, const std::string &xAxis, const std::string &yAxis)
 Set the axis labels for the given 1D histogram. More...
 
GAUDI_API bool setAxisLabels (AIDA::IProfile1D *hist, const std::string &xAxis, const std::string &yAxis)
 Set the axis labels for the given 1D profile histogram. More...
 
GAUDI_API bool setAxisLabels (AIDA::IHistogram2D *hist, const std::string &xAxis, const std::string &yAxis)
 Set the axis labels for the given 2D histogram. More...
 
GAUDI_API bool setAxisLabels (AIDA::IProfile2D *hist, const std::string &xAxis, const std::string &yAxis)
 Set the axis labels for the given 2D profile histogram. More...
 
GAUDI_API std::string path (const AIDA::IBaseHistogram *aida)
 get the path in THS for AIDA histogram More...
 
GAUDI_API std::string format (const AIDA::IHistogram1D *histo, const std::string &fmt)
 Make the string representation of the histogram according to the specified format. More...
 
GAUDI_API std::string format (const AIDA::IProfile1D *histo, const std::string &fmt)
 Make the string representation of the profile histogram according to the specified format. More...
 
GAUDI_API std::string format (const AIDA::IHistogram1D *histo, const std::string &ID, const std::string &fmt1, const std::string &fmt2)
 format a full row in table, including ID, label, path or any other "extra" identifier in string form More...
 
GAUDI_API std::string format (const AIDA::IProfile1D *histo, const std::string &ID, const std::string &fmt1, const std::string &fmt2)
 format a full row in table, including ID, label, path or any other "extra" identifier in string form More...
 
template<class HISTO , class STREAM , class TERMINATOR >
STREAM & printList (HISTO first, HISTO last, const std::string &fmt, STREAM &stream, TERMINATOR term)
 print the simple sequence (list-like) of histograms as table More...
 
template<class LIST , class STREAM , class TERMINATOR >
STREAM & printList (const LIST &histos, const std::string &fmt, STREAM &stream, TERMINATOR term)
 print the simple container of histograms as table More...
 
template<class HISTO , class STREAM , class TERMINATOR >
STREAM & printMap (HISTO begin, HISTO end, const std::string &fmt1, const std::string &fmt2, STREAM &stream, TERMINATOR term)
 Print the "associative sequence" (e.g. More...
 
template<class MAP , class STREAM , class TERMINATOR >
STREAM & printMap (const MAP &histos, const std::string &fmt1, const std::string &fmt2, STREAM &stream, TERMINATOR term)
 Print the "associative sequence" (e.g. More...
 
GAUDI_API std::string format (const std::string &val1, const std::string &val2, const std::string &fmt)
 helper method to merge the headers for short format table More...
 
GAUDI_API std::ostreamtoXml (const TH1D &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TH2D &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TH3D &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TProfile &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TProfile2D &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const AIDA::IHistogram1D &histo, std::ostream &stream)
 stream the AIDA histogram into the output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const AIDA::IHistogram2D &histo, std::ostream &stream)
 stream the AIDA histogram into the output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const AIDA::IHistogram3D &histo, std::ostream &stream)
 stream the AIDA histogram into the output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const AIDA::IProfile1D &histo, std::ostream &stream)
 stream the AIDA histogram into the output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TH1F &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TH2F &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const TH3F &histo, std::ostream &stream)
 stream the ROOT histogram into output stream as XML More...
 
GAUDI_API std::ostreamtoXml (const AIDA::IProfile2D &histo, std::ostream &stream)
 stream the AIDA histogram into the output stream as XML More...
 
GAUDI_API StatusCode fromXml (TH1D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH2D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH3D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TProfile &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TProfile2D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH1F &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH2F &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH3F &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH1D *&result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH2D *&result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TH3D *&result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TProfile *&result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (TProfile2D *&result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (AIDA::IHistogram1D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (AIDA::IHistogram2D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (AIDA::IHistogram3D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (AIDA::IProfile1D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 
GAUDI_API StatusCode fromXml (AIDA::IProfile2D &result, const std::string &input)
 parse the histogram from standard ROOT XML More...
 

Detailed Description

Collection of useful utilities for manipulations with AIDA hisgograms.

Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-08

Typedef Documentation

typedef std::pair<unsigned, std::string> Gaudi::Utils::Histos::BinLabel

Typedef for a bin number and its associated label.

Definition at line 25 of file HistoLabels.h.

typedef std::vector<BinLabel> Gaudi::Utils::Histos::BinLabels

Typedef for a list of bin numbers and their associated label.

Definition at line 27 of file HistoLabels.h.

typedef std::vector<std::string> Gaudi::Utils::Histos::Labels

Typedef for a list of labels.

Definition at line 23 of file HistoLabels.h.

Function Documentation

void Gaudi::Utils::Histos::fill ( AIDA::IHistogram1D *  histo,
const double  value,
const double  weight = 1.0 
)

simple function to fill AIDA::IHistogram1D objects

See also
AIDA::IHistogram1D
Parameters
histopointer to the histogram
valuevalue to be added to the histogram
weightthe "weight" assciated with this entry
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-10-02

Definition at line 35 of file Fill.cpp.

35  {
36  if ( histo ) { histo->fill( value, weight ); }
37 }
void Gaudi::Utils::Histos::fill ( AIDA::IHistogram2D *  histo,
const double  valueX,
const double  valueY,
const double  weight = 1.0 
)

simple function to fill AIDA::IHistogram2D objects

See also
AIDA::IHistogram2D
Parameters
histopointer to the histogram
valueXvalue to be added to the histogram
valueYvalue to be added to the histogram
weightthe "weight" assciated with this entry
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-10-02

Definition at line 49 of file Fill.cpp.

50  {
51  if ( histo ) { histo->fill( valueX, valueY, weight ); }
52 }
void Gaudi::Utils::Histos::fill ( AIDA::IHistogram3D *  histo,
const double  valueX,
const double  valueY,
const double  valueZ,
const double  weight = 1.0 
)

simple function to fill AIDA::IHistogram3D objects

See also
AIDA::IHistogram3D
Parameters
histopointer to the histogram
valueXvalue to be added to the histogram
valueYvalue to be added to the histogram
valueZvalue to be added to the histogram
weightthe "weight" assciated with this entry
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-10-02

Definition at line 65 of file Fill.cpp.

66  {
67  if ( histo ) { histo->fill( valueX, valueY, valueZ, weight ); }
68 }
void Gaudi::Utils::Histos::fill ( AIDA::IProfile1D *  histo,
const double  valueX,
const double  valueY,
const double  weight = 1.0 
)

simple function to fill AIDA::IProfile1D objects

See also
AIDA::IProfile1D
Parameters
histopointer to the histogram
valueXvalue to be added to the histogram
valueYvalue to be added to the histogram
weightthe "weight" assciated with this entry
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-10-02

Definition at line 80 of file Fill.cpp.

81  {
82  if ( histo ) { histo->fill( valueX, valueY, weight ); }
83 }
void Gaudi::Utils::Histos::fill ( AIDA::IProfile2D *  histo,
const double  valueX,
const double  valueY,
const double  valueZ,
const double  weight = 1.0 
)

simple function to fill AIDA::IProfile2D objects

See also
AIDA::IProfile2D
Parameters
histopointer to the histogram
valueXvalue to be added to the histogram
valueYvalue to be added to the histogram
valueZvalue to be added to the histogram
weightthe "weight" assciated with this entry
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-10-02

Definition at line 96 of file Fill.cpp.

97  {
98  if ( histo ) { histo->fill( valueX, valueY, valueZ, weight ); }
99 }
std::string Gaudi::Utils::Histos::format ( const AIDA::IHistogram1D *  histo,
const std::string fmt 
)

Make the string representation of the histogram according to the specified format.

The method could be used to access/print various quantities

using namespace Gaudi::Utils::Histos ;
const AIDA::IHistogram1D* histo = ... ;
// print title in a free format:
std::cout << format ( histo , "%2%" ) << std::endl ;
// print the path in HTS in a free format:
std::cout << format ( histo , " path in HTS: %1% " ) << std::endl ;
// print the formatted nEntries/Overflow/Underflow:
format ( histo , " #nEntries/Overflow/Underflow=%3%/%4%/%5%" )
<< std::endl ;
// print the formatted Mean+-ErrorMean:
format ( histo , " Mean+-Error=(%8%+-%9%)" )
<< std::endl ;
// print the skewness and kurtosis:
format ( histo , " Skewness/Kurtosis=%12%/%14% " )
<< std::endl ;
Parameters
historeference to the histogram
fmtthe printout format
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-08

Definition at line 234 of file HistoTableFormat.cpp.

234  {
235  if ( !histo ) { return "<NULL>"; }
236  using namespace Gaudi::Utils;
237  using namespace boost::io;
238  boost::format _fmt( fmt );
239  // allow various number of arguments
240  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
241 
242  _fmt % ( "\"" + path( histo ) + "\"" ) // 1) histogram path
243  % ( "\"" + _title( histo ) + "\"" ) // 2) title
244  % histo->allEntries() // 3) # entries
245  % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) // 4) # underflow
246  % histo->binEntries( AIDA::IAxis::OVERFLOW_BIN ) // 5) # overflow
247  % HistoStats::nEff( histo ) // 6) equivalent entries
248  % histo->sumBinHeights() // 7) integral
249  % HistoStats::mean( histo ) // 8) mean value
250  % HistoStats::meanErr( histo ) // 9) error in mean
251  % HistoStats::rms( histo ) // 10) rms
252  % HistoStats::rmsErr( histo ) // 11) error in rms
253  % HistoStats::skewness( histo ) // 12) skewness
254  % HistoStats::skewnessErr( histo ) // 13) error in skewness
255  % HistoStats::kurtosis( histo ) // 14) kurtosis
256  % HistoStats::kurtosisErr( histo ) // 15) error in kurtosis
257  //
258  % histo->sumAllBinHeights() // 16) full integral (in and out range)
259  % HistoStats::sumAllBinHeightErr( histo ) // 17) error on 16
260  % HistoStats::sumBinHeightErr( histo ) // 18) error on 7
261  //
262  % ( 100 * HistoStats::underflowEntriesFrac( histo ) ) // 19) fraction of underflow entries
263  % ( 100 * HistoStats::underflowEntriesFracErr( histo ) ) // 20) error on 19
264  % ( 100 * HistoStats::overflowEntriesFrac( histo ) ) // 21) fraction of overflow entries
265  % ( 100 * HistoStats::overflowEntriesFracErr( histo ) ) // 22) error on 21
266  //
267  % HistoStats::underflowIntegralFrac( histo ) // 23) fraction of underflow integral
268  % HistoStats::underflowIntegralFracErr( histo ) // 24) error on 23
269  % HistoStats::overflowIntegralFrac( histo ) // 25) fraction of overflow intergal
270  % HistoStats::overflowIntegralFracErr( histo ) // 26) error on 25
271  ;
272 
273  return _fmt.str();
274 }
static double underflowIntegralFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of underflow integral
Definition: HistoStats.cpp:765
static double overflowIntegralFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of overflow intergal
Definition: HistoStats.cpp:753
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
static double underflowEntriesFrac(const AIDA::IHistogram1D *histo)
the fraction of underflow entries (useful for shape comparison)
Definition: HistoStats.cpp:693
static double rmsErr(const AIDA::IHistogram1D *histo)
get an error in the rms value
Definition: HistoStats.cpp:651
static double underflowIntegralFrac(const AIDA::IHistogram1D *histo)
the fraction of underflow integral (useful for shape comparison)
Definition: HistoStats.cpp:717
static double meanErr(const AIDA::IHistogram1D *histo)
get an error in the mean value
Definition: HistoStats.cpp:635
static double overflowIntegralFrac(const AIDA::IHistogram1D *histo)
the fraction of overflow intergal (useful for shape comparison)
Definition: HistoStats.cpp:705
static double rms(const AIDA::IHistogram1D *histo)
get the rms value for the histogram (just for completeness)
Definition: HistoStats.cpp:643
static double overflowEntriesFrac(const AIDA::IHistogram1D *histo)
the fraction of overflow entries (useful for shape comparison)
Definition: HistoStats.cpp:681
static double sumBinHeightErr(const AIDA::IHistogram1D *histo)
get an error in the sum bin height ("in-range integral")
Definition: HistoStats.cpp:659
static double overflowEntriesFracErr(const AIDA::IHistogram1D *histo)
error on fraction of overflow entries (useful for shape comparison)
Definition: HistoStats.cpp:729
static double sumAllBinHeightErr(const AIDA::IHistogram1D *histo)
get an error in the sum of all bin height ("integral")
Definition: HistoStats.cpp:669
static double kurtosisErr(const AIDA::IHistogram1D *histo)
get the error in kurtosis for the histogram
Definition: HistoStats.cpp:611
static double mean(const AIDA::IHistogram1D *histo)
get the mean value for the histogram (just for completeness)
Definition: HistoStats.cpp:627
static double skewnessErr(const AIDA::IHistogram1D *histo)
get the error in skewness for the histogram
Definition: HistoStats.cpp:595
static double underflowEntriesFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of underflow entries (useful for shape comparison)
Definition: HistoStats.cpp:741
static double kurtosis(const AIDA::IHistogram1D *histo)
get the kurtosis for the histogram
Definition: HistoStats.cpp:603
static double nEff(const AIDA::IHistogram1D *histo)
get the effective entries (just for completeness)
Definition: HistoStats.cpp:619
static double skewness(const AIDA::IHistogram1D *histo)
get the skewness for the histogram
Definition: HistoStats.cpp:587
std::string Gaudi::Utils::Histos::format ( const AIDA::IProfile1D *  histo,
const std::string fmt 
)

Make the string representation of the profile histogram according to the specified format.

The method could be used to access/print various quantities

using namespace Gaudi::Utils::Histos ;
const AIDA::IProfile1D* histo = ... ;
// print title in a free format:
std::cout << format ( histo , "%2%" ) << std::endl ;
// print the path in HTS in a free format:
std::cout << format ( histo , " path in HTS: %1% " ) << std::endl ;
Parameters
historeference to the histogram
fmtthe printout format
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-08

Definition at line 280 of file HistoTableFormat.cpp.

280  {
281  if ( !histo ) { return "<NULL>"; }
282  using namespace Gaudi::Utils;
283  using namespace boost::io;
284  boost::format _fmt( fmt );
285  // allow various number of arguments
286  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
287 
288  _fmt % ( "\"" + path( histo ) + "\"" ) // 1) histogram path
289  % ( "\"" + _title( histo ) + "\"" ) // 2) title
290  % histo->allEntries() // 3) # entries
291  % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) // 4) # underflow
292  % histo->binEntries( AIDA::IAxis::OVERFLOW_BIN ) // 5) # overflow
293  % HistoStats::nEff( histo ) // 6) equivalent entries
294  % histo->sumBinHeights() // 7) integral
295  % HistoStats::mean( histo ) // 8) mean value
296  % HistoStats::meanErr( histo ) // 9) error in mean
297  % HistoStats::rms( histo ) // 10) rms
298  % HistoStats::rmsErr( histo ) // 11) error in rms
299  % HistoStats::skewness( histo ) // 12) skewness
300  % HistoStats::skewnessErr( histo ) // 13) error in skewness
301  % HistoStats::kurtosis( histo ) // 14) kurtosis
302  % HistoStats::kurtosisErr( histo ) // 15) error in kurtosis
303  //
304  % histo->sumAllBinHeights() // 16) full integral (in and out range)
305  % HistoStats::sumAllBinHeightErr( histo ) // 17) error on 16
306  % HistoStats::sumBinHeightErr( histo ) // 18) error on 7
307  //
308  % ( 100 * HistoStats::underflowEntriesFrac( histo ) ) // 19) fraction of underflow entries
309  % ( 100 * HistoStats::underflowEntriesFracErr( histo ) ) // 20) error on 19
310  % ( 100 * HistoStats::overflowEntriesFrac( histo ) ) // 21) fraction of overflow entries
311  % ( 100 * HistoStats::overflowEntriesFracErr( histo ) ) // 22) error on 21
312  //
313  % HistoStats::underflowIntegralFrac( histo ) // 23) fraction of underflow integral
314  % HistoStats::underflowIntegralFracErr( histo ) // 24) error on 23
315  % HistoStats::overflowIntegralFrac( histo ) // 25) fraction of overflow intergal
316  % HistoStats::overflowIntegralFracErr( histo ) // 26) error on 25
317  ;
318 
319  return _fmt.str();
320 }
static double underflowIntegralFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of underflow integral
Definition: HistoStats.cpp:765
static double overflowIntegralFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of overflow intergal
Definition: HistoStats.cpp:753
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
static double underflowEntriesFrac(const AIDA::IHistogram1D *histo)
the fraction of underflow entries (useful for shape comparison)
Definition: HistoStats.cpp:693
static double rmsErr(const AIDA::IHistogram1D *histo)
get an error in the rms value
Definition: HistoStats.cpp:651
static double underflowIntegralFrac(const AIDA::IHistogram1D *histo)
the fraction of underflow integral (useful for shape comparison)
Definition: HistoStats.cpp:717
static double meanErr(const AIDA::IHistogram1D *histo)
get an error in the mean value
Definition: HistoStats.cpp:635
static double overflowIntegralFrac(const AIDA::IHistogram1D *histo)
the fraction of overflow intergal (useful for shape comparison)
Definition: HistoStats.cpp:705
static double rms(const AIDA::IHistogram1D *histo)
get the rms value for the histogram (just for completeness)
Definition: HistoStats.cpp:643
static double overflowEntriesFrac(const AIDA::IHistogram1D *histo)
the fraction of overflow entries (useful for shape comparison)
Definition: HistoStats.cpp:681
static double sumBinHeightErr(const AIDA::IHistogram1D *histo)
get an error in the sum bin height ("in-range integral")
Definition: HistoStats.cpp:659
static double overflowEntriesFracErr(const AIDA::IHistogram1D *histo)
error on fraction of overflow entries (useful for shape comparison)
Definition: HistoStats.cpp:729
static double sumAllBinHeightErr(const AIDA::IHistogram1D *histo)
get an error in the sum of all bin height ("integral")
Definition: HistoStats.cpp:669
static double kurtosisErr(const AIDA::IHistogram1D *histo)
get the error in kurtosis for the histogram
Definition: HistoStats.cpp:611
static double mean(const AIDA::IHistogram1D *histo)
get the mean value for the histogram (just for completeness)
Definition: HistoStats.cpp:627
static double skewnessErr(const AIDA::IHistogram1D *histo)
get the error in skewness for the histogram
Definition: HistoStats.cpp:595
static double underflowEntriesFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of underflow entries (useful for shape comparison)
Definition: HistoStats.cpp:741
static double kurtosis(const AIDA::IHistogram1D *histo)
get the kurtosis for the histogram
Definition: HistoStats.cpp:603
static double nEff(const AIDA::IHistogram1D *histo)
get the effective entries (just for completeness)
Definition: HistoStats.cpp:619
static double skewness(const AIDA::IHistogram1D *histo)
get the skewness for the histogram
Definition: HistoStats.cpp:587
std::string Gaudi::Utils::Histos::format ( const AIDA::IHistogram1D *  histo,
const std::string ID,
const std::string fmt1,
const std::string fmt2 
)

format a full row in table, including ID, label, path or any other "extra" identifier in string form

using namespace Gaudi::Utils::Histos
const AIDA::IHistogram1D* histo = ... ;
// define short format
const std::string fmt1 = " |%1$-30.30s %|33t| %2" ;
// define format for the histogram
const std::string fmt2 = ... ;
info () <<
format ( "My Histo" , histo , fmt1 , fmt2 )
<< endmsg ;
Parameters
histopointer to the histogram
IDhistorgam ID, title, label or other extra information
fmt1"short" format used for the table
fmt2format used for the histogram printout
Returns
formatted row
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-08

Definition at line 326 of file HistoTableFormat.cpp.

327  {
328  using namespace boost::io;
329  boost::format _fmt( fmt1 );
330  // allow various number of arguments
331  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
332  //
333  _fmt % ID // format ID
334  % format( histo, fmt2 ); // format the histogram
335  //
336  return _fmt.str();
337 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
std::string Gaudi::Utils::Histos::format ( const AIDA::IProfile1D *  histo,
const std::string ID,
const std::string fmt1,
const std::string fmt2 
)

format a full row in table, including ID, label, path or any other "extra" identifier in string form

using namespace Gaudi::Utils::Histos
const AIDA::IProfile1D* histo = ... ;
// define short format
const std::string fmt1 = " |%1$-30.30s %|33t| %2" ;
// define format for the histogram
const std::string fmt2 = ... ;
info () <<
format ( "My Histo" , histo , fmt1 , fmt2 )
<< endmsg ;
Parameters
histopointer to the histogram
IDhistorgam ID, title, label or other extra information
fmt1"short" format used for the table
fmt2format used for the histogram printout
Returns
formatted row
Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-08

Definition at line 343 of file HistoTableFormat.cpp.

344  {
345  using namespace boost::io;
346  boost::format _fmt( fmt1 );
347  // allow various number of arguments
348  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
349  //
350  _fmt % ID // format ID
351  % format( histo, fmt2 ); // format the histogram
352  //
353  return _fmt.str();
354 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
std::string Gaudi::Utils::Histos::format ( const std::string val1,
const std::string val2,
const std::string fmt 
)

helper method to merge the headers for short format table

Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-07

Definition at line 358 of file HistoTableFormat.cpp.

358  {
359  using namespace boost::io;
360  boost::format _fmt( fmt );
361  // allow various number of arguments
362  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
363  //
364  _fmt % val1 // format ID
365  % val2; // format the histogram
366  //
367  return _fmt.str();
368 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
StatusCode Gaudi::Utils::Histos::fromXml ( TH1D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 187 of file HistoXML.cpp.

187  {
188  //
189  result.Reset(); // RESET old histogram
190  //
191 
192  auto histo = _Xml<TH1D>( input );
193  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
194  //
195  histo->Copy( result );
196  //
197  return StatusCode::SUCCESS;
198 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH2D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 206 of file HistoXML.cpp.

206  {
207  //
208  result.Reset(); // RESET old histogram
209  //
210  auto histo = _Xml<TH2D>( input );
211  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
212  //
213  histo->Copy( result );
214  //
215  return StatusCode::SUCCESS;
216 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH3D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 224 of file HistoXML.cpp.

224  {
225  //
226  result.Reset(); // RESET old histogram
227  //
228  auto histo = _Xml<TH3D>( input );
229  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
230  //
231  histo->Copy( result );
232  //
233  return StatusCode::SUCCESS;
234 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TProfile &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 296 of file HistoXML.cpp.

296  {
297  //
298  result.Reset(); // RESET old histogram
299  //
300  auto histo = _Xml<TProfile>( input );
301  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
302  //
303  histo->Copy( result );
304  //
305  return StatusCode::SUCCESS;
306 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TProfile2D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 314 of file HistoXML.cpp.

314  {
315  //
316  result.Reset(); // RESET old histogram
317  //
318  auto histo = _Xml<TProfile2D>( input );
319  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
320  //
321  histo->Copy( result );
322  //
323  return StatusCode::SUCCESS;
324 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH1F &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 242 of file HistoXML.cpp.

242  {
243  //
244  result.Reset(); // RESET old histogram
245  //
246  auto histo = _Xml<TH1F>( input );
247  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
248  //
249  histo->Copy( result );
250  //
251  return StatusCode::SUCCESS;
252 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH2F &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 260 of file HistoXML.cpp.

260  {
261  //
262  result.Reset(); // RESET old histogram
263  //
264  auto histo = _Xml<TH2F>( input );
265  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
266  //
267  histo->Copy( result );
268  //
269  return StatusCode::SUCCESS;
270 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH3F &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 278 of file HistoXML.cpp.

278  {
279  //
280  result.Reset(); // RESET old histogram
281  //
282  auto histo = _Xml<TH3F>( input );
283  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
284  //
285  histo->Copy( result );
286  //
287  return StatusCode::SUCCESS;
288 }
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH1D *&  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 333 of file HistoXML.cpp.

333  {
334  if ( result ) { return fromXml( *result, input ); }
335  //
336  auto histo = _Xml<TH1D>( input );
337  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
338  //
339  result = histo.release(); // ASSIGN
340  //
341  return StatusCode::SUCCESS;
342 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH2D *&  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 350 of file HistoXML.cpp.

350  {
351  if ( result ) { return fromXml( *result, input ); }
352  //
353  auto histo = _Xml<TH2D>( input );
354  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
355  //
356  result = histo.release(); // ASSIGN
357  //
358  return StatusCode::SUCCESS;
359 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TH3D *&  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 367 of file HistoXML.cpp.

367  {
368  if ( result ) { return fromXml( *result, input ); }
369  //
370  auto histo = _Xml<TH3D>( input );
371  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
372  //
373  result = histo.release(); // ASSIGN
374  //
375  return StatusCode::SUCCESS;
376 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TProfile *&  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 385 of file HistoXML.cpp.

385  {
386  if ( result ) { return fromXml( *result, input ); }
387  //
388  auto histo = _Xml<TProfile>( input );
389  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
390  //
391  result = histo.release(); // ASSIGN
392  //
393  return StatusCode::SUCCESS;
394 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( TProfile2D *&  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 402 of file HistoXML.cpp.

402  {
403  if ( result ) { return fromXml( *result, input ); }
404  //
405  auto histo = _Xml<TProfile2D>( input );
406  if ( !histo ) { return StatusCode::FAILURE; } // RETURN
407  //
408  result = histo.release(); // ASSIGN
409  //
410  return StatusCode::SUCCESS;
411 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
constexpr static const auto SUCCESS
Definition: StatusCode.h:85
constexpr static const auto FAILURE
Definition: StatusCode.h:86
StatusCode Gaudi::Utils::Histos::fromXml ( AIDA::IHistogram1D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 419 of file HistoXML.cpp.

419  {
420  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
421  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
422 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
StatusCode Gaudi::Utils::Histos::fromXml ( AIDA::IHistogram2D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 430 of file HistoXML.cpp.

430  {
431  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
432  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
433 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
StatusCode Gaudi::Utils::Histos::fromXml ( AIDA::IHistogram3D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 441 of file HistoXML.cpp.

441  {
442  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
443  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
444 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
StatusCode Gaudi::Utils::Histos::fromXml ( AIDA::IProfile1D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 452 of file HistoXML.cpp.

452  {
453  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
454  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
455 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
StatusCode Gaudi::Utils::Histos::fromXml ( AIDA::IProfile2D &  result,
const std::string input 
)

parse the histogram from standard ROOT XML

See also
TBufferXML
TBufferXML::ConvertFromXML
Parameters
result(OUTPUT) the parsed histogram
input(INPUT) the input XML string
Returns
status code

Definition at line 463 of file HistoXML.cpp.

463  {
464  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
465  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
466 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:187
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
std::string Gaudi::Utils::Histos::histoDump ( const AIDA::IHistogram1D *  histo,
const std::size_t  width = 80,
const std::size_t  height = 50,
const bool  errors = false 
)

dump the text representation of the histogram

Parameters
histo(INPUT) the histogram
width(INPUT) the maximal column width
height(INPUT) the proposed column height
erorrs(INPUT) print/plot errors
Returns
string representation of the histogram
Author
Vanya BELYAEV Ivan..nosp@m.BEly.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 643 of file HistoDump.cpp.

644  {
645  std::ostringstream stream;
646  histoDump_( histo, stream, width, height, errors );
647  return stream.str();
648 }
GAUDI_API std::ostream & histoDump_(const AIDA::IHistogram1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
dump the text representation of the histogram
Definition: HistoDump.cpp:584
std::string Gaudi::Utils::Histos::histoDump ( const AIDA::IProfile1D *  histo,
const std::size_t  width = 80,
const std::size_t  height = 50,
const bool  spread = true 
)

dump the text representation of the 1D-profile

Parameters
histo(INPUT) the histogram
width(INPUT) the maximal column width
height(INPUT) the proposed column height
spread(INPUT) print/plto spread vs rms
Returns
string representation of the histogram
Author
Vanya BELYAEV Ivan..nosp@m.BEly.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 712 of file HistoDump.cpp.

713  {
714  std::ostringstream stream;
715  histoDump_( histo, stream, width, height, spread );
716  return stream.str();
717 }
GAUDI_API std::ostream & histoDump_(const AIDA::IHistogram1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
dump the text representation of the histogram
Definition: HistoDump.cpp:584
std::string Gaudi::Utils::Histos::histoDump ( const TProfile *  histo,
const std::size_t  width = 80,
const std::size_t  height = 50 
)

dump the text representation of the histogram

Parameters
histo(INPUT) the histogram
width(INPUT) the maximal column width
height(INPUT) the propsoed coulmn height
erorrs(INPUT) print/plot errors
Returns
string representation of the histogram
Author
Vanya BELYAEV Ivan..nosp@m.Bely.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 824 of file HistoDump.cpp.

825  {
826  std::ostringstream stream;
827  histoDump_( histo, stream, width, height );
828  return stream.str();
829 }
GAUDI_API std::ostream & histoDump_(const AIDA::IHistogram1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
dump the text representation of the histogram
Definition: HistoDump.cpp:584
std::string Gaudi::Utils::Histos::histoDump ( const TH1 *  histo,
const std::size_t  width = 80,
const std::size_t  height = 50,
const bool  errors = false 
)

dump the text representation of the histogram

Parameters
histo(INPUT) the histogram
width(INPUT) the maximal column width
height(INPUT) the propsoed coulmn height
erorrs(INPUT) print/plot errors
Returns
string representation of the histogram
Author
Vanya BELYAEV Ivan..nosp@m.Bely.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 808 of file HistoDump.cpp.

809  {
810  std::ostringstream stream;
811  histoDump_( histo, stream, width, height, errors );
812  return stream.str();
813 }
GAUDI_API std::ostream & histoDump_(const AIDA::IHistogram1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
dump the text representation of the histogram
Definition: HistoDump.cpp:584
std::ostream & Gaudi::Utils::Histos::histoDump_ ( const AIDA::IHistogram1D *  histo,
std::ostream stream,
const std::size_t  width = 80,
const std::size_t  height = 50,
const bool  errors = false 
)

dump the text representation of the histogram

Parameters
histo(INPUT) the histogram
stream(OUTUT) the stream
width(INPUT) the maximal column width
height(INPUT) the proposed column height
errors(INPUT) print/plot errors
Returns
the stream
Author
Vanya BELYAEV Ivan..nosp@m.BEly.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 584 of file HistoDump.cpp.

585  {
586  stream << std::endl;
587  if ( !histo ) { return stream; } // RETURN
588  Histo hist;
589  StatusCode sc = _getHisto( histo, hist );
590  if ( sc.isFailure() ) { return stream; } // RETURN
591  //
592  stream << boost::format( " Histo TES : \"%s\"" ) % path( histo ) << std::endl
593  << boost::format( " Histo Title : \"%s\"" ) % histo->title() << std::endl
594  << std::endl;
595  //
596  stream << boost::format( " Mean : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::mean( histo ) %
598  << std::endl
599  << boost::format( " Rms : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::rms( histo ) %
601  << std::endl
602  << boost::format( " Skewness : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::skewness( histo ) %
604  << std::endl
605  << boost::format( " Kurtosis : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::kurtosis( histo ) %
607  << std::endl
608  << std::endl;
609  //
610  stream << boost::format( " Entries :\n | %=9s | %=9s | %=9s | %9s | %=11s | %=11s | %=11s |" ) % "All" %
611  "In Range" % "Underflow" % "Overflow" % "#Equivalent" % "Integral" % "Total"
612  << std::endl
613  << boost::format( " | %=9d | %=9d | %=9d | %=9d | %=11.5g | %=11.5g | %=11.5g |" ) % histo->allEntries() %
614  histo->entries() % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) %
615  histo->binEntries( AIDA::IAxis::OVERFLOW_BIN ) % histo->equivalentBinEntries() %
616  histo->sumBinHeights() % histo->sumAllBinHeights()
617  << std::endl
618  << std::endl;
619  //
620  const AIDA::IAnnotation& a = histo->annotation();
621  if ( 0 != a.size() ) {
622  stream << " Annotation" << std::endl;
623  for ( int i = 0; i < a.size(); ++i ) {
624  stream << boost::format( " | %-25.25s : %-45.45s | " ) % a.key( i ) % a.value( i ) << std::endl;
625  }
626  stream << std::endl;
627  }
628  //
629  return dumpText( hist, width, height, errors, stream );
630 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
static double rmsErr(const AIDA::IHistogram1D *histo)
get an error in the rms value
Definition: HistoStats.cpp:651
T endl(T...args)
static double meanErr(const AIDA::IHistogram1D *histo)
get an error in the mean value
Definition: HistoStats.cpp:635
bool isFailure() const
Definition: StatusCode.h:130
static double rms(const AIDA::IHistogram1D *histo)
get the rms value for the histogram (just for completeness)
Definition: HistoStats.cpp:643
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
static double kurtosisErr(const AIDA::IHistogram1D *histo)
get the error in kurtosis for the histogram
Definition: HistoStats.cpp:611
static double mean(const AIDA::IHistogram1D *histo)
get the mean value for the histogram (just for completeness)
Definition: HistoStats.cpp:627
static double skewnessErr(const AIDA::IHistogram1D *histo)
get the error in skewness for the histogram
Definition: HistoStats.cpp:595
static double kurtosis(const AIDA::IHistogram1D *histo)
get the kurtosis for the histogram
Definition: HistoStats.cpp:603
static double skewness(const AIDA::IHistogram1D *histo)
get the skewness for the histogram
Definition: HistoStats.cpp:587
std::ostream & Gaudi::Utils::Histos::histoDump_ ( const AIDA::IProfile1D *  histo,
std::ostream stream,
const std::size_t  width = 80,
const std::size_t  height = 50,
const bool  spread = true 
)

dump the text representation of 1D-profile

Parameters
histo(INPUT) the 1D-profile
stream(OUTUT) the stream
width(INPUT) the maximal column width
height(INPUT) the proposed column height
spread(INPUT) print/plot spread/rms ?
Returns
the stream
Author
Vanya BELYAEV Ivan..nosp@m.BEly.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 661 of file HistoDump.cpp.

662  {
663  stream << std::endl;
664  if ( !histo ) { return stream; } // RETURN
665  Histo hist;
666  StatusCode sc = _getHisto( histo, hist, spread );
667  if ( sc.isFailure() ) { return stream; } // RETURN
668  //
669  stream << boost::format( " Histo TES : \"%s\"" ) % path( histo ) << std::endl
670  << boost::format( " Histo Title : \"%s\"" ) % histo->title() << std::endl
671  << std::endl;
672  //
673  stream << boost::format( " Mean : %11.5g " ) % histo->mean() << std::endl
674  << boost::format( " Rms : %11.5g " ) % histo->rms() << std::endl
675  << std::endl;
676  //
677  stream << boost::format( " Entries :\n | %=9s | %=9s | %=9s | %9s | %=11s | %=11s |" ) % "All" % "In Range" %
678  "Underflow" %
679  "Overflow"
680  // % "#Equivalent"
681  % "Integral" % "Total"
682  << std::endl
683  << boost::format( " | %=9d | %=9d | %=9d | %=9d | %=11.5g | %=11.5g |" ) % histo->allEntries() %
684  histo->entries() % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) %
685  histo->binEntries( AIDA::IAxis::OVERFLOW_BIN )
686  // % histo -> equivalentBinEntries ()
687  % histo->sumBinHeights() % histo->sumAllBinHeights()
688  << std::endl
689  << std::endl;
690  //
691  const AIDA::IAnnotation& a = histo->annotation();
692  if ( 0 != a.size() ) {
693  stream << " Annotation" << std::endl;
694  for ( int i = 0; i < a.size(); ++i ) {
695  stream << boost::format( " | %-25.25s : %-45.45s | " ) % a.key( i ) % a.value( i ) << std::endl;
696  }
697  stream << std::endl;
698  }
699  //
700  return dumpText( hist, width, height, true, stream );
701 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
T endl(T...args)
bool isFailure() const
Definition: StatusCode.h:130
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
std::ostream & Gaudi::Utils::Histos::histoDump_ ( const TProfile *  histo,
std::ostream stream,
const std::size_t  width = 80,
const std::size_t  height = 50 
)

dump the text representation of the Profile

Parameters
histo(INPUT) the histogram
stream(OUTUT) the stream
width(INPUT) the maximal column width
height(INPUT) the proposed coulmn height
spread(INPUT) print/plot rms versus erorr
Returns
the stream
Author
Vanya BELYAEV Ivan..nosp@m.BEly.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 772 of file HistoDump.cpp.

773  {
774  stream << std::endl;
775  if ( !histo ) { return stream; } // RETURN
776  Histo hist;
777  StatusCode sc = _getHisto( histo, hist, true );
778  if ( sc.isFailure() ) { return stream; } // RETURN
779  //
780  stream << boost::format( " Profile Name : \"%s\"" ) % histo->GetName() << std::endl
781  << boost::format( " Profile Title : \"%s\"" ) % histo->GetTitle() << std::endl
782  << std::endl;
783  //
784  stream << boost::format( " Mean : %11.5g " ) % histo->GetMean() << std::endl
785  << boost::format( " Rms : %11.5g " ) % histo->GetRMS() << std::endl
786  << std::endl;
787  //
788  stream << boost::format( " Entries :\n | %=11s | %=11s | %=11s | %=11s |" ) % "All" % "Underflow" % "Overflow" %
789  "Integral"
790  << std::endl
791  << boost::format( " | %=11.5g | %=11.5g | %=11.5g | %=11.5g |" ) % histo->GetEntries() %
792  histo->GetBinContent( 0 ) % histo->GetBinContent( histo->GetNbinsX() + 1 ) % histo->Integral()
793  << std::endl
794  << std::endl;
795  //
796  return dumpText( hist, width, height, true, stream );
797 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
T endl(T...args)
bool isFailure() const
Definition: StatusCode.h:130
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
std::ostream & Gaudi::Utils::Histos::histoDump_ ( const TH1 *  histo,
std::ostream stream,
const std::size_t  width = 80,
const std::size_t  height = 50,
const bool  errors = false 
)

dump the text representation of the histogram

Parameters
histo(INPUT) the histogram
stream(OUTUT) the stream
width(INPUT) the maximal column width
height(INPUT) the proposed coulmn height
errors(INPUT) print/plot errors
Returns
the stream
Author
Vanya BELYAEV Ivan..nosp@m.BEly.nosp@m.aev@n.nosp@m.ikhe.nosp@m.f.nl
Date
2009-09-19

Definition at line 729 of file HistoDump.cpp.

730  {
731  const TProfile* profile = dynamic_cast<const TProfile*>( histo );
732  if ( profile ) { return histoDump_( profile, stream, width, height ); }
733  //
734  stream << std::endl;
735  if ( !histo ) { return stream; } // RETURN
736  Histo hist;
737  StatusCode sc = _getHisto( histo, hist );
738  if ( sc.isFailure() ) { return stream; } // RETURN
739  //
740  stream << boost::format( " Histo Name : \"%s\"" ) % histo->GetName() << std::endl
741  << boost::format( " Histo Title : \"%s\"" ) % histo->GetTitle() << std::endl
742  << std::endl;
743  //
744  stream << boost::format( " Mean : %11.5g +- %-10.4g " ) % histo->GetMean() % histo->GetMeanError() << std::endl
745  << boost::format( " Rms : %11.5g +- %-10.4g " ) % histo->GetRMS() % histo->GetRMSError() << std::endl
746  << boost::format( " Skewness : %11.5g " ) % histo->GetSkewness() << std::endl
747  << boost::format( " Kurtosis : %11.5g " ) % histo->GetKurtosis() << std::endl
748  << std::endl;
749  //
750  stream << boost::format( " Entries :\n | %=11s | %=11s | %=11s | %=11s | %=11s |" ) % "All" % "Underflow" %
751  "Overflow" % "#Equivalent" % "Integral"
752  << std::endl
753  << boost::format( " | %=11.5g | %=11.5g | %=11.5g | %=11.5g | %=11.5g |" ) % histo->GetEntries() %
754  histo->GetBinContent( 0 ) % histo->GetBinContent( histo->GetNbinsX() + 1 ) %
755  histo->GetEffectiveEntries() % histo->Integral()
756  << std::endl
757  << std::endl;
758  //
759  return dumpText( hist, width, height, errors, stream );
760 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:109
GAUDI_API std::ostream & histoDump_(const AIDA::IHistogram1D *histo, std::ostream &stream, const std::size_t width=80, const std::size_t height=50, const bool errors=false)
dump the text representation of the histogram
Definition: HistoDump.cpp:584
T endl(T...args)
bool isFailure() const
Definition: StatusCode.h:130
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:50
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IBaseHistogram *  histo,
const std::string title = "" 
)

get the title

Definition at line 109 of file Fill.cpp.

109  {
110  return htitle_( histo, title );
111 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IHistogram *  histo,
const std::string title = "" 
)

get the title

Definition at line 115 of file Fill.cpp.

115  {
116  return htitle_( histo, title );
117 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IHistogram1D *  histo,
const std::string title = "" 
)

get the title

Definition at line 121 of file Fill.cpp.

121  {
122  return htitle_( histo, title );
123 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IHistogram2D *  histo,
const std::string title = "" 
)

get the title

Definition at line 127 of file Fill.cpp.

127  {
128  return htitle_( histo, title );
129 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IHistogram3D *  histo,
const std::string title = "" 
)

get the title

Definition at line 133 of file Fill.cpp.

133  {
134  return htitle_( histo, title );
135 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IProfile *  histo,
const std::string title = "" 
)

get the title

Definition at line 139 of file Fill.cpp.

139  {
140  return htitle_( histo, title );
141 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IProfile1D *  histo,
const std::string title = "" 
)

get the title

Definition at line 145 of file Fill.cpp.

145  {
146  return htitle_( histo, title );
147 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IProfile2D *  histo,
const std::string title = "" 
)

get the title

Definition at line 151 of file Fill.cpp.

151  {
152  return htitle_( histo, title );
153 }
std::string Gaudi::Utils::Histos::path ( const AIDA::IBaseHistogram *  aida)

get the path in THS for AIDA histogram

Author
Vanya BELYAEV ibely.nosp@m.aev@.nosp@m.physi.nosp@m.cs.s.nosp@m.yr.ed.nosp@m.u
Date
2007-08-08

Definition at line 219 of file HistoTableFormat.cpp.

219  {
220  if ( !aida ) { return ""; } // RETURN
221  const auto object = dynamic_cast<const DataObject*>( aida );
222  if ( !object ) { return ""; } // RETURN
223  const auto registry = object->registry();
224  if ( !registry ) { return ""; } // RETURN
225  const auto _path = registry->identifier();
226  const auto n = _path.find( "/stat/" );
227  return ( 0 == n ? std::string( _path, 6 ) : _path ); // RETURN
228 }
STL class.
A DataObject is the base class of any identifiable object on any data store.
Definition: DataObject.h:30
template<class HISTO , class STREAM , class TERMINATOR >
STREAM& Gaudi::Utils::Histos::printList ( HISTO  first,
HISTO  last,
const std::string fmt,
STREAM &  stream,
TERMINATOR  term 
)
inline

print the simple sequence (list-like) of histograms as table

SEQUENCE histos = ... ;
// print a table with three colums path, title and #entries
( histos.begin () ,
histos.end () ,
" | %1$|-40.40s | %2$-30.30s | %3$=7d | " ,
'\n' ) ;
Parameters
firstbegin-iterator for the sequence
lastend-iterator for the sequence
streamthe stream to be used for printout
termthe terminmator for the stream
fmtthe format to be used

Definition at line 277 of file HistoTableFormat.h.

277  {
278  for ( ; first != last; ++first ) { stream << format( *first, fmt ) << term; } // print table rows
279  return stream; // RETURN
280  }
GAUDI_API std::string format(const std::string &val1, const std::string &val2, const std::string &fmt)
helper method to merge the headers for short format table
template<class LIST , class STREAM , class TERMINATOR >
STREAM& Gaudi::Utils::Histos::printList ( const LIST &  histos,
const std::string fmt,
STREAM &  stream,
TERMINATOR  term 
)
inline

print the simple container of histograms as table

using namespace Gaudi::Utils::Histos ;
SEQUENCE histos = ... ;
// print a table with three columns: path, title and #entries
( histos ,
" | %1$|-40.40s | %2$-30.30s | %3$=7d | " ,
'\n' ) ;
Parameters
histosthe sequence of histograms
streamthe stream to be used for printout
termthe terminmator for the stream
fmtthe format to be used

Definition at line 306 of file HistoTableFormat.h.

306  {
307  return printList( histos.begin(), histos.end(), fmt, stream, term );
308  }
STREAM & printList(const LIST &histos, const std::string &fmt, STREAM &stream, TERMINATOR term)
print the simple container of histograms as table
template<class HISTO , class STREAM , class TERMINATOR >
STREAM& Gaudi::Utils::Histos::printMap ( HISTO  begin,
HISTO  end,
const std::string fmt1,
const std::string fmt2,
STREAM &  stream,
TERMINATOR  term 
)
inline

Print the "associative sequence" (e.g.

part of std:map) of histograms as table:

using namespace Gaudi::Utils::Histos ;
( m.begin () ,
m.end () ,
"| %1$-10.10s | %2% " , // short format
Gaudi::Utils::Histos::Formats::histoFormatOnly ,
always() ,
endmsg ) ;

Print only mean and rms:

using namespace Gaudi::Utils::Histos ;
( m.begin () ,
m.end () ,
"| %1$-10.10s | %2% " , // short format
" %8$10.5g+-%10$-10.5g|", // mean+-rms
always() ,
endmsg ) ;
Parameters
begin'begin'-iterator for the mapping sequence
end'end'-iterator for the mapping sequence
fmt1'short' format for the table printout
fmt3format for the printout of the histogram
streamthe stream for printout
termstream terminator

Definition at line 355 of file HistoTableFormat.h.

356  {
357  for ( ; begin != end; ++begin ) {
358  stream << format( begin->second, // the histogram
359  begin->first, // the key
360  fmt1, fmt2 )
361  << term;
362  }
363  return stream;
364  }
GAUDI_API std::string format(const std::string &val1, const std::string &val2, const std::string &fmt)
helper method to merge the headers for short format table
AttribStringParser::Iterator begin(const AttribStringParser &parser)
template<class MAP , class STREAM , class TERMINATOR >
STREAM& Gaudi::Utils::Histos::printMap ( const MAP &  histos,
const std::string fmt1,
const std::string fmt2,
STREAM &  stream,
TERMINATOR  term 
)
inline

Print the "associative sequence" (e.g.

part of std:map) of histograms as table:

using namespace Gaudi::Utils::Histos ;
( m ,
"| %1$-10.10s | %2% " , // short format
Gaudi::Utils::Histos::Formats::histoFormatOnly ,
always() ,
endmsg ) ;

Print only mean and rms:

using namespace Gaudi::Utils::Histos ;
( m ,
"| %1$-10.10s | %2% " , // short format
" %8$10.5g+-%10$-10.5g|", // mean+-rms
always() ,
endmsg ) ;
Parameters
begin'begin'-iterator for the mapping sequence
end'end'-iterator for the mapping sequence
fmt1'short' format for the table printout
fmt3format for the printout of the histogram
streamthe stream for printout
termstream terminator

Definition at line 409 of file HistoTableFormat.h.

410  {
411  return printMap( histos.begin(), histos.end(), fmt1, fmt2, stream, term );
412  }
STREAM & printMap(const MAP &histos, const std::string &fmt1, const std::string &fmt2, STREAM &stream, TERMINATOR term)
Print the "associative sequence" (e.g.
bool Gaudi::Utils::Histos::setAxisLabels ( AIDA::IHistogram1D *  hist,
const std::string xAxis,
const std::string yAxis 
)

Set the axis labels for the given 1D histogram.

Parameters
histPointer to the histogram
xAxisLabel for the x axis
yAxisLabel for the y axis
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 129 of file HistoLabels.cpp.

129  {
130  return setAxisLabels_<TH1D>( hist, xAxis, yAxis );
131  }
bool Gaudi::Utils::Histos::setAxisLabels ( AIDA::IProfile1D *  hist,
const std::string xAxis,
const std::string yAxis 
)

Set the axis labels for the given 1D profile histogram.

Parameters
histPointer to the histogram
xAxisLabel for the x axis
yAxisLabel for the y axis
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 133 of file HistoLabels.cpp.

133  {
134  return setAxisLabels_<TProfile>( hist, xAxis, yAxis );
135  }
bool Gaudi::Utils::Histos::setAxisLabels ( AIDA::IHistogram2D *  hist,
const std::string xAxis,
const std::string yAxis 
)

Set the axis labels for the given 2D histogram.

Parameters
histPointer to the histogram
xAxisLabel for the x axis
yAxisLabel for the y axis
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 137 of file HistoLabels.cpp.

137  {
138  return setAxisLabels_<TH2D>( hist, xAxis, yAxis );
139  }
bool Gaudi::Utils::Histos::setAxisLabels ( AIDA::IProfile2D *  hist,
const std::string xAxis,
const std::string yAxis 
)

Set the axis labels for the given 2D profile histogram.

Parameters
histPointer to the histogram
xAxisLabel for the x axis
yAxisLabel for the y axis
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 141 of file HistoLabels.cpp.

141  {
142  return setAxisLabels_<TProfile2D>( hist, xAxis, yAxis );
143  }
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IHistogram1D *  hist,
const Labels labels 
)

Set the Bin labels for a given 1D histogram.

The labels will be applied in the order they appear in the list, starting at the first bin. If the list of labels is too short, the later bins will be missing a label. If the list is too long, only the first N will be used, where N is the number of bins in the histogram

Parameters
histPointer to the histogram
labelsThe list of labels
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 87 of file HistoLabels.cpp.

87 { return setBinLabels_( hist, labels ); }
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IHistogram1D *  hist,
const BinLabels labels 
)

Set the Bin labels for a given 1D histogram.

Each entry in 'labels' gives the bin number and its associated label

Parameters
histPointer to the histogram
labelsThe list of labels
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 79 of file HistoLabels.cpp.

79  {
80  return setBinLabels_<TH1D>( hist, labels );
81  }
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IProfile1D *  hist,
const Labels labels 
)

Set the Bin labels for a given 1D profile histogram.

The labels will be applied in the order they appear in the list, starting at the first bin. If the list of labels is too short, the later bins will be missing a label. If the list is too long, only the first N will be used, where N is the number of bins in the histogram

Parameters
histPointer to the histogram
labelsThe list of labels
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 89 of file HistoLabels.cpp.

89 { return setBinLabels_( hist, labels ); }
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IProfile1D *  hist,
const BinLabels labels 
)

Set the Bin labels for a given 1D profile histogram.

Each entry in 'labels' gives the bin number and its associated label

Parameters
histPointer to the histogram
labelsThe list of bin numbers and the associated label
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 83 of file HistoLabels.cpp.

83  {
84  return setBinLabels_<TProfile>( hist, labels );
85  }
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IHistogram2D *  hist,
const Labels xlabels,
const Labels ylabels 
)

Set the Bin labels for a given 2D histogram.

The labels will be applied in the order they appear in the lists, starting at the first bin. If the list of labels is too short, the later bins will be missing a label. If the list is too long, only the first N will be used, where N is the number of bins in the histogram

Parameters
histPointer to the histogram
xlabelsThe list of x labels
ylabelsThe list of y labels
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 91 of file HistoLabels.cpp.

91  {
92  if ( !hist ) return false;
93  TH2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
94  if ( !h2d ) return false;
95  BinLabels lx;
96  lx.reserve( xlabels.size() );
97  for ( unsigned int i = 0; i < xlabels.size(); ++i ) { lx.emplace_back( i, xlabels[i] ); }
98  BinLabels ly;
99  ly.reserve( ylabels.size() );
100  for ( unsigned int i = 0; i < ylabels.size(); ++i ) { ly.emplace_back( i, ylabels[i] ); }
101  return ( setBinLabels_( h2d->GetXaxis(), lx ) && setBinLabels_( h2d->GetYaxis(), ly ) );
102  }
std::vector< BinLabel > BinLabels
Typedef for a list of bin numbers and their associated label.
Definition: HistoLabels.h:27
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IHistogram2D *  hist,
const BinLabels xlabels,
const BinLabels ylabels 
)

Set the Bin labels for a given 2D histogram.

Each entry in 'labels' lists gives the bin number and its associated label

Parameters
histPointer to the histogram
xlabelsThe list of x bin numbers and the associated label
ylabelsThe list of y bin numbers and the associated label
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 104 of file HistoLabels.cpp.

104  {
105  TH2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
106  return ( h2d && setBinLabels_( h2d->GetXaxis(), xlabels ) && setBinLabels_( h2d->GetYaxis(), ylabels ) );
107  }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IProfile2D *  hist,
const Labels xlabels,
const Labels ylabels 
)

Set the Bin labels for a given 2D profile histogram.

The labels will be applied in the order they appear in the lists, starting at the first bin. If the list of labels is too short, the later bins will be missing a label. If the list is too long, only the first N will be used, where N is the number of bins in the histogram

Parameters
histPointer to the histogram
xlabelsThe list of x labels
ylabelsThe list of y labels
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 109 of file HistoLabels.cpp.

109  {
110  if ( !hist ) return false;
111  TProfile2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
112  if ( !h2d ) return false;
113  BinLabels lx;
114  lx.reserve( xlabels.size() );
115  for ( unsigned int i = 0; i < xlabels.size(); ++i ) { lx.emplace_back( i, xlabels[i] ); }
116  BinLabels ly;
117  ly.reserve( ylabels.size() );
118  for ( unsigned int i = 0; i < ylabels.size(); ++i ) { ly.emplace_back( i, ylabels[i] ); }
119  return ( setBinLabels_( h2d->GetXaxis(), lx ) && setBinLabels_( h2d->GetYaxis(), ly ) );
120  }
std::vector< BinLabel > BinLabels
Typedef for a list of bin numbers and their associated label.
Definition: HistoLabels.h:27
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
bool Gaudi::Utils::Histos::setBinLabels ( AIDA::IProfile2D *  hist,
const BinLabels xlabels,
const BinLabels ylabels 
)

Set the Bin labels for a given 2D profile histogram.

Each entry in 'labels' lists gives the bin number and its associated label

Parameters
histPointer to the histogram
xlabelsThe list of x bin numbers and the associated label
ylabelsThe list of y bin numbers and the associated label
Returns
Boolean indicating if the labels were successfully applied or not
Return values
TRUELabels were applied OK
FALSELabels were NOT applied

Definition at line 122 of file HistoLabels.cpp.

122  {
123  TProfile2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
124  return ( h2d && setBinLabels_( h2d->GetXaxis(), xlabels ) && setBinLabels_( h2d->GetYaxis(), ylabels ) );
125  }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
AIDA::IBaseHistogram * Gaudi::Utils::Histos::toBase ( AIDA::IHistogram1D *  histo)

Definition at line 155 of file Fill.cpp.

155 { return histo; }
AIDA::IBaseHistogram * Gaudi::Utils::Histos::toBase ( AIDA::IHistogram2D *  histo)

Definition at line 157 of file Fill.cpp.

157 { return histo; }
AIDA::IBaseHistogram * Gaudi::Utils::Histos::toBase ( AIDA::IHistogram3D *  histo)

Definition at line 159 of file Fill.cpp.

159 { return histo; }
AIDA::IBaseHistogram * Gaudi::Utils::Histos::toBase ( AIDA::IProfile1D *  histo)

Definition at line 161 of file Fill.cpp.

161 { return histo; }
AIDA::IBaseHistogram * Gaudi::Utils::Histos::toBase ( AIDA::IProfile2D *  histo)

Definition at line 163 of file Fill.cpp.

163 { return histo; }
std::ostream & Gaudi::Utils::Histos::toXml ( const TH1D &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 64 of file HistoXML.cpp.

64  {
65  return stream << TBufferXML::ConvertToXML( &histo );
66 }
std::ostream & Gaudi::Utils::Histos::toXml ( const TH2D &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 73 of file HistoXML.cpp.

73  {
74  return stream << TBufferXML::ConvertToXML( &histo );
75 }
std::ostream & Gaudi::Utils::Histos::toXml ( const TH3D &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 82 of file HistoXML.cpp.

82  {
83  return stream << TBufferXML::ConvertToXML( &histo );
84 }
std::ostream & Gaudi::Utils::Histos::toXml ( const TProfile &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 118 of file HistoXML.cpp.

118  {
119  return stream << TBufferXML::ConvertToXML( &histo );
120 }
std::ostream & Gaudi::Utils::Histos::toXml ( const TProfile2D &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 127 of file HistoXML.cpp.

127  {
128  return stream << TBufferXML::ConvertToXML( &histo );
129 }
std::ostream & Gaudi::Utils::Histos::toXml ( const AIDA::IHistogram1D &  histo,
std::ostream stream 
)

stream the AIDA histogram into the output stream as XML

Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 136 of file HistoXML.cpp.

136  {
137  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
138  return root ? toXml( *root, stream ) : stream;
139 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:64
std::ostream & Gaudi::Utils::Histos::toXml ( const AIDA::IHistogram2D &  histo,
std::ostream stream 
)

stream the AIDA histogram into the output stream as XML

Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 146 of file HistoXML.cpp.

146  {
147  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
148  return root ? toXml( *root, stream ) : stream;
149 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:64
std::ostream & Gaudi::Utils::Histos::toXml ( const AIDA::IHistogram3D &  histo,
std::ostream stream 
)

stream the AIDA histogram into the output stream as XML

Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 156 of file HistoXML.cpp.

156  {
157  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
158  return root ? toXml( *root, stream ) : stream;
159 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:64
std::ostream & Gaudi::Utils::Histos::toXml ( const AIDA::IProfile1D &  histo,
std::ostream stream 
)

stream the AIDA histogram into the output stream as XML

Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 166 of file HistoXML.cpp.

166  {
167  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
168  return root ? toXml( *root, stream ) : stream;
169 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:64
std::ostream & Gaudi::Utils::Histos::toXml ( const TH1F &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 91 of file HistoXML.cpp.

91  {
92  return stream << TBufferXML::ConvertToXML( &histo );
93 }
std::ostream & Gaudi::Utils::Histos::toXml ( const TH2F &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 100 of file HistoXML.cpp.

100  {
101  return stream << TBufferXML::ConvertToXML( &histo );
102 }
std::ostream & Gaudi::Utils::Histos::toXml ( const TH3F &  histo,
std::ostream stream 
)

stream the ROOT histogram into output stream as XML

See also
TBufferXML
TBufferXML::ConvertToXML
Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 109 of file HistoXML.cpp.

109  {
110  return stream << TBufferXML::ConvertToXML( &histo );
111 }
std::ostream & Gaudi::Utils::Histos::toXml ( const AIDA::IProfile2D &  histo,
std::ostream stream 
)

stream the AIDA histogram into the output stream as XML

Parameters
histo(INPUT) the histogram to be streamed
stream(OUTPUT) the stream

Definition at line 176 of file HistoXML.cpp.

176  {
177  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
178  return root ? toXml( *root, stream ) : stream;
179 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:71
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:64