The Gaudi Framework  v29r0 (ff2e7097)
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 historgam 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...
 
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 29 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 31 of file HistoLabels.h.

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

Typedef for a list of labels.

Definition at line 27 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.

36 {
37  if ( histo ) {
38  histo->fill( value, weight );
39  }
40 }
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 52 of file Fill.cpp.

54 {
55  if ( histo ) {
56  histo->fill( valueX, valueY, weight );
57  }
58 }
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 71 of file Fill.cpp.

73 {
74  if ( histo ) {
75  histo->fill( valueX, valueY, valueZ, weight );
76  }
77 }
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 89 of file Fill.cpp.

91 {
92  if ( histo ) {
93  histo->fill( valueX, valueY, weight );
94  }
95 }
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 108 of file Fill.cpp.

110 {
111  if ( histo ) {
112  histo->fill( valueX, valueY, valueZ, weight );
113  }
114 }
std::string Gaudi::Utils::Histos::format ( const AIDA::IHistogram1D *  histo,
const std::string fmt 
)

Make the string representation of the historgam 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.

235 {
236  if ( 0 == histo ) {
237  return "<NULL>";
238  }
239  using namespace Gaudi::Utils;
240  using namespace boost::io;
241  boost::format _fmt( fmt );
242  // allow various number of arguments
243  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
244 
245  _fmt % ( "\"" + path( histo ) + "\"" ) // 1) histogram path
246  % ( "\"" + histo->title() + "\"" ) // 2) title
247  % histo->allEntries() // 3) # entries
248  % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) // 4) # underflow
249  % histo->binEntries( AIDA::IAxis::OVERFLOW_BIN ) // 5) # overflow
250  % histo->equivalentBinEntries() // 6) equivalent entries
251  % histo->sumBinHeights() // 7) integral
252  % histo->mean() // 8) mean value
253  % HistoStats::meanErr( histo ) // 9) error in mean
254  % histo->rms() // 10) rms
255  % HistoStats::rmsErr( histo ) // 11) error in rms
256  % HistoStats::skewness( histo ) // 12) skewness
257  % HistoStats::skewnessErr( histo ) // 13) error in skewness
258  % HistoStats::kurtosis( histo ) // 14) kurtosis
259  % HistoStats::kurtosisErr( histo ) // 15) error in kurtosis
260  //
261  % histo->sumAllBinHeights() // 16) full integral (in and out range)
262  % HistoStats::sumAllBinHeightErr( histo ) // 17) error on 16
263  % HistoStats::sumBinHeightErr( histo ) // 18) error on 7
264  //
265  % ( 100 * HistoStats::underflowEntriesFrac( histo ) ) // 19) fraction of underflow entries
266  % ( 100 * HistoStats::underflowEntriesFracErr( histo ) ) // 20) error on 19
267  % ( 100 * HistoStats::overflowEntriesFrac( histo ) ) // 21) fraction of overflow entries
268  % ( 100 * HistoStats::overflowEntriesFracErr( histo ) ) // 22) error on 21
269  //
270  % HistoStats::underflowIntegralFrac( histo ) // 23) fraction of underflow integral
271  % HistoStats::underflowIntegralFracErr( histo ) // 24) error on 23
272  % HistoStats::overflowIntegralFrac( histo ) // 25) fraction of overflow intergal
273  % HistoStats::overflowIntegralFracErr( histo ); // 26) error on 25
274  //
275  return _fmt.str();
276 }
static double underflowIntegralFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of underflow integral
Definition: HistoStats.cpp:462
static double overflowIntegralFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of overflow intergal
Definition: HistoStats.cpp:433
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
static double underflowEntriesFrac(const AIDA::IHistogram1D *histo)
the fraction of underflow entries (useful for shape comparison)
Definition: HistoStats.cpp:329
static double rmsErr(const AIDA::IHistogram1D *histo)
get an error in the rms value
Definition: HistoStats.cpp:250
static double underflowIntegralFrac(const AIDA::IHistogram1D *histo)
the fraction of underflow integral (useful for shape comparison)
Definition: HistoStats.cpp:370
static double meanErr(const AIDA::IHistogram1D *histo)
get an error in the mean value
Definition: HistoStats.cpp:235
static double overflowIntegralFrac(const AIDA::IHistogram1D *histo)
the fraction of overflow intergal (useful for shape comparison)
Definition: HistoStats.cpp:351
static double overflowEntriesFrac(const AIDA::IHistogram1D *histo)
the fraction of overflow entries (useful for shape comparison)
Definition: HistoStats.cpp:307
static double sumBinHeightErr(const AIDA::IHistogram1D *histo)
get an error in the sum bin height ("in-range integral")
Definition: HistoStats.cpp:267
static double overflowEntriesFracErr(const AIDA::IHistogram1D *histo)
error on fraction of overflow entries (useful for shape comparison)
Definition: HistoStats.cpp:389
static double sumAllBinHeightErr(const AIDA::IHistogram1D *histo)
get an error in the sum of all bin height ("integral")
Definition: HistoStats.cpp:288
static double kurtosisErr(const AIDA::IHistogram1D *histo)
get the error in kurtosis for the histogram
Definition: HistoStats.cpp:206
static double skewnessErr(const AIDA::IHistogram1D *histo)
get the error in skewness for the histogram
Definition: HistoStats.cpp:179
static double underflowEntriesFracErr(const AIDA::IHistogram1D *histo)
the error on fraction of underflow entries (useful for shape comparison)
Definition: HistoStats.cpp:411
static double kurtosis(const AIDA::IHistogram1D *histo)
get the kurtosis for the histogram
Definition: HistoStats.cpp:194
static double skewness(const AIDA::IHistogram1D *histo)
get the skewness for the histogram
Definition: HistoStats.cpp:167
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 282 of file HistoTableFormat.cpp.

284 {
285  using namespace boost::io;
286  boost::format _fmt( fmt1 );
287  // allow various number of arguments
288  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
289  //
290  _fmt % ID // format ID
291  % format( histo, fmt2 ); // format the histogram
292  //
293  return _fmt.str();
294 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
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 298 of file HistoTableFormat.cpp.

299 {
300  using namespace boost::io;
301  boost::format _fmt( fmt );
302  // allow various number of arguments
303  _fmt.exceptions( all_error_bits ^ ( too_many_args_bit | too_few_args_bit ) );
304  //
305  _fmt % val1 // format ID
306  % val2; // format the histogram
307  //
308  return _fmt.str();
309 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
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 202 of file HistoXML.cpp.

203 {
204  //
205  result.Reset(); // RESET old histogram
206  //
207 
208  auto histo = _Xml<TH1D>( input );
209  if ( !histo ) {
210  return StatusCode::FAILURE;
211  } // RETURN
212  //
213  histo->Copy( result );
214  //
215  return StatusCode::SUCCESS;
216 }
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 224 of file HistoXML.cpp.

225 {
226  //
227  result.Reset(); // RESET old histogram
228  //
229  auto histo = _Xml<TH2D>( input );
230  if ( !histo ) {
231  return StatusCode::FAILURE;
232  } // RETURN
233  //
234  histo->Copy( result );
235  //
236  return StatusCode::SUCCESS;
237 }
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 245 of file HistoXML.cpp.

246 {
247  //
248  result.Reset(); // RESET old histogram
249  //
250  auto histo = _Xml<TH3D>( input );
251  if ( !histo ) {
252  return StatusCode::FAILURE;
253  } // RETURN
254  //
255  histo->Copy( result );
256  //
257  return StatusCode::SUCCESS;
258 }
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 329 of file HistoXML.cpp.

330 {
331  //
332  result.Reset(); // RESET old histogram
333  //
334  auto histo = _Xml<TProfile>( input );
335  if ( !histo ) {
336  return StatusCode::FAILURE;
337  } // RETURN
338  //
339  histo->Copy( result );
340  //
341  return StatusCode::SUCCESS;
342 }
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 350 of file HistoXML.cpp.

351 {
352  //
353  result.Reset(); // RESET old histogram
354  //
355  auto histo = _Xml<TProfile2D>( input );
356  if ( !histo ) {
357  return StatusCode::FAILURE;
358  } // RETURN
359  //
360  histo->Copy( result );
361  //
362  return StatusCode::SUCCESS;
363 }
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 266 of file HistoXML.cpp.

267 {
268  //
269  result.Reset(); // RESET old histogram
270  //
271  auto histo = _Xml<TH1F>( input );
272  if ( !histo ) {
273  return StatusCode::FAILURE;
274  } // RETURN
275  //
276  histo->Copy( result );
277  //
278  return StatusCode::SUCCESS;
279 }
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 287 of file HistoXML.cpp.

288 {
289  //
290  result.Reset(); // RESET old histogram
291  //
292  auto histo = _Xml<TH2F>( input );
293  if ( !histo ) {
294  return StatusCode::FAILURE;
295  } // RETURN
296  //
297  histo->Copy( result );
298  //
299  return StatusCode::SUCCESS;
300 }
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 308 of file HistoXML.cpp.

309 {
310  //
311  result.Reset(); // RESET old histogram
312  //
313  auto histo = _Xml<TH3F>( input );
314  if ( !histo ) {
315  return StatusCode::FAILURE;
316  } // RETURN
317  //
318  histo->Copy( result );
319  //
320  return StatusCode::SUCCESS;
321 }
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 372 of file HistoXML.cpp.

373 {
374  if ( result ) {
375  return fromXml( *result, input );
376  }
377  //
378  auto histo = _Xml<TH1D>( input );
379  if ( !histo ) {
380  return StatusCode::FAILURE;
381  } // RETURN
382  //
383  result = histo.release(); // ASSIGN
384  //
385  return StatusCode::SUCCESS;
386 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
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 394 of file HistoXML.cpp.

395 {
396  if ( result ) {
397  return fromXml( *result, input );
398  }
399  //
400  auto histo = _Xml<TH2D>( input );
401  if ( !histo ) {
402  return StatusCode::FAILURE;
403  } // RETURN
404  //
405  result = histo.release(); // ASSIGN
406  //
407  return StatusCode::SUCCESS;
408 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
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 416 of file HistoXML.cpp.

417 {
418  if ( result ) {
419  return fromXml( *result, input );
420  }
421  //
422  auto histo = _Xml<TH3D>( input );
423  if ( !histo ) {
424  return StatusCode::FAILURE;
425  } // RETURN
426  //
427  result = histo.release(); // ASSIGN
428  //
429  return StatusCode::SUCCESS;
430 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
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 439 of file HistoXML.cpp.

440 {
441  if ( result ) {
442  return fromXml( *result, input );
443  }
444  //
445  auto histo = _Xml<TProfile>( input );
446  if ( !histo ) {
447  return StatusCode::FAILURE;
448  } // RETURN
449  //
450  result = histo.release(); // ASSIGN
451  //
452  return StatusCode::SUCCESS;
453 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
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 461 of file HistoXML.cpp.

462 {
463  if ( result ) {
464  return fromXml( *result, input );
465  }
466  //
467  auto histo = _Xml<TProfile2D>( input );
468  if ( !histo ) {
469  return StatusCode::FAILURE;
470  } // RETURN
471  //
472  result = histo.release(); // ASSIGN
473  //
474  return StatusCode::SUCCESS;
475 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
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 483 of file HistoXML.cpp.

484 {
485  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
486  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
487 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 495 of file HistoXML.cpp.

496 {
497  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
498  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
499 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 507 of file HistoXML.cpp.

508 {
509  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
510  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
511 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 519 of file HistoXML.cpp.

520 {
521  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
522  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
523 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 531 of file HistoXML.cpp.

532 {
533  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &result );
534  return root ? fromXml( *root, input ) : StatusCode::FAILURE; // RETURN
535 }
GAUDI_API StatusCode fromXml(TH1D &result, const std::string &input)
parse the histogram from standard ROOT XML
Definition: HistoXML.cpp:202
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 703 of file HistoDump.cpp.

705 {
706  std::ostringstream stream;
707  histoDump_( histo, stream, width, height, errors );
708  return stream.str();
709 }
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:639
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 777 of file HistoDump.cpp.

779 {
780  std::ostringstream stream;
781  histoDump_( histo, stream, width, height, spread );
782  return stream.str();
783 }
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:639
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 903 of file HistoDump.cpp.

904 {
905  std::ostringstream stream;
906  histoDump_( histo, stream, width, height );
907  return stream.str();
908 }
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:639
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 886 of file HistoDump.cpp.

888 {
889  std::ostringstream stream;
890  histoDump_( histo, stream, width, height, errors );
891  return stream.str();
892 }
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:639
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 639 of file HistoDump.cpp.

641 {
642  stream << std::endl;
643  if ( !histo ) {
644  return stream;
645  } // RETURN
646  Histo hist;
647  StatusCode sc = _getHisto( histo, hist );
648  if ( sc.isFailure() ) {
649  return stream;
650  } // RETURN
651  //
652  stream << boost::format( " Histo TES : \"%s\"" ) % path( histo ) << std::endl
653  << boost::format( " Histo Title : \"%s\"" ) % histo->title() << std::endl
654  << std::endl;
655  //
656  stream << boost::format( " Mean : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::mean( histo ) %
658  << std::endl
659  << boost::format( " Rms : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::rms( histo ) %
661  << std::endl
662  << boost::format( " Skewness : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::skewness( histo ) %
664  << std::endl
665  << boost::format( " Kurtosis : %11.5g +- %-10.4g " ) % Gaudi::Utils::HistoStats::kurtosis( histo ) %
667  << std::endl
668  << std::endl;
669  //
670  stream << boost::format( " Entries :\n | %=9s | %=9s | %=9s | %9s | %=11s | %=11s | %=11s |" ) % "All" %
671  "In Range" % "Underflow" % "Overflow" % "#Equivalent" % "Integral" % "Total"
672  << std::endl
673  << boost::format( " | %=9d | %=9d | %=9d | %=9d | %=11.5g | %=11.5g | %=11.5g |" ) % histo->allEntries() %
674  histo->entries() % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) %
675  histo->binEntries( AIDA::IAxis::OVERFLOW_BIN ) % histo->equivalentBinEntries() %
676  histo->sumBinHeights() % histo->sumAllBinHeights()
677  << std::endl
678  << std::endl;
679  //
680  const AIDA::IAnnotation& a = histo->annotation();
681  if ( 0 != a.size() ) {
682  stream << " Annotation" << std::endl;
683  for ( int i = 0; i < a.size(); ++i ) {
684  stream << boost::format( " | %-25.25s : %-45.45s | " ) % a.key( i ) % a.value( i ) << std::endl;
685  }
686  stream << std::endl;
687  }
688  //
689  return dumpText( hist, width, height, errors, stream );
690 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
static double rmsErr(const AIDA::IHistogram1D *histo)
get an error in the rms value
Definition: HistoStats.cpp:250
T endl(T...args)
static double meanErr(const AIDA::IHistogram1D *histo)
get an error in the mean value
Definition: HistoStats.cpp:235
bool isFailure() const
Test for a status code of FAILURE.
Definition: StatusCode.h:86
static double rms(const AIDA::IHistogram1D *histo)
get the rms value for the histogram (just for completeness)
Definition: HistoStats.cpp:246
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
static double kurtosisErr(const AIDA::IHistogram1D *histo)
get the error in kurtosis for the histogram
Definition: HistoStats.cpp:206
static double mean(const AIDA::IHistogram1D *histo)
get the mean value for the histogram (just for completeness)
Definition: HistoStats.cpp:231
static double skewnessErr(const AIDA::IHistogram1D *histo)
get the error in skewness for the histogram
Definition: HistoStats.cpp:179
static double kurtosis(const AIDA::IHistogram1D *histo)
get the kurtosis for the histogram
Definition: HistoStats.cpp:194
static double skewness(const AIDA::IHistogram1D *histo)
get the skewness for the histogram
Definition: HistoStats.cpp:167
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 722 of file HistoDump.cpp.

724 {
725  stream << std::endl;
726  if ( !histo ) {
727  return stream;
728  } // RETURN
729  Histo hist;
730  StatusCode sc = _getHisto( histo, hist, spread );
731  if ( sc.isFailure() ) {
732  return stream;
733  } // RETURN
734  //
735  stream << boost::format( " Histo TES : \"%s\"" ) % path( histo ) << std::endl
736  << boost::format( " Histo Title : \"%s\"" ) % histo->title() << std::endl
737  << std::endl;
738  //
739  stream << boost::format( " Mean : %11.5g " ) % histo->mean() << std::endl
740  << boost::format( " Rms : %11.5g " ) % histo->rms() << std::endl
741  << std::endl;
742  //
743  stream << boost::format( " Entries :\n | %=9s | %=9s | %=9s | %9s | %=11s | %=11s |" ) % "All" % "In Range" %
744  "Underflow" % "Overflow"
745  // % "#Equivalent"
746  % "Integral" % "Total"
747  << std::endl
748  << boost::format( " | %=9d | %=9d | %=9d | %=9d | %=11.5g | %=11.5g |" ) % histo->allEntries() %
749  histo->entries() % histo->binEntries( AIDA::IAxis::UNDERFLOW_BIN ) %
750  histo->binEntries( AIDA::IAxis::OVERFLOW_BIN )
751  // % histo -> equivalentBinEntries ()
752  % histo->sumBinHeights() % histo->sumAllBinHeights()
753  << std::endl
754  << std::endl;
755  //
756  const AIDA::IAnnotation& a = histo->annotation();
757  if ( 0 != a.size() ) {
758  stream << " Annotation" << std::endl;
759  for ( int i = 0; i < a.size(); ++i ) {
760  stream << boost::format( " | %-25.25s : %-45.45s | " ) % a.key( i ) % a.value( i ) << std::endl;
761  }
762  stream << std::endl;
763  }
764  //
765  return dumpText( hist, width, height, true, stream );
766 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
T endl(T...args)
bool isFailure() const
Test for a status code of FAILURE.
Definition: StatusCode.h:86
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
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 845 of file HistoDump.cpp.

847 {
848  stream << std::endl;
849  if ( !histo ) {
850  return stream;
851  } // RETURN
852  Histo hist;
853  StatusCode sc = _getHisto( histo, hist, true );
854  if ( sc.isFailure() ) {
855  return stream;
856  } // RETURN
857  //
858  stream << boost::format( " Profile Name : \"%s\"" ) % histo->GetName() << std::endl
859  << boost::format( " Profile Title : \"%s\"" ) % histo->GetTitle() << std::endl
860  << std::endl;
861  //
862  stream << boost::format( " Mean : %11.5g " ) % histo->GetMean() << std::endl
863  << boost::format( " Rms : %11.5g " ) % histo->GetRMS() << std::endl
864  << std::endl;
865  //
866  stream << boost::format( " Entries :\n | %=11s | %=11s | %=11s | %=11s |" ) % "All" % "Underflow" % "Overflow" %
867  "Integral"
868  << std::endl
869  << boost::format( " | %=11.5g | %=11.5g | %=11.5g | %=11.5g |" ) % histo->GetEntries() %
870  histo->GetBinContent( 0 ) % histo->GetBinContent( histo->GetNbinsX() + 1 ) % histo->Integral()
871  << std::endl
872  << std::endl;
873  //
874  return dumpText( hist, width, height, true, stream );
875 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
T endl(T...args)
bool isFailure() const
Test for a status code of FAILURE.
Definition: StatusCode.h:86
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
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 795 of file HistoDump.cpp.

797 {
798  const TProfile* profile = dynamic_cast<const TProfile*>( histo );
799  if ( profile ) {
800  return histoDump_( profile, stream, width, height );
801  }
802  //
803  stream << std::endl;
804  if ( !histo ) {
805  return stream;
806  } // RETURN
807  Histo hist;
808  StatusCode sc = _getHisto( histo, hist );
809  if ( sc.isFailure() ) {
810  return stream;
811  } // RETURN
812  //
813  stream << boost::format( " Histo Name : \"%s\"" ) % histo->GetName() << std::endl
814  << boost::format( " Histo Title : \"%s\"" ) % histo->GetTitle() << std::endl
815  << std::endl;
816  //
817  stream << boost::format( " Mean : %11.5g +- %-10.4g " ) % histo->GetMean() % histo->GetMeanError() << std::endl
818  << boost::format( " Rms : %11.5g +- %-10.4g " ) % histo->GetRMS() % histo->GetRMSError() << std::endl
819  << boost::format( " Skewness : %11.5g " ) % histo->GetSkewness() << std::endl
820  << boost::format( " Kurtosis : %11.5g " ) % histo->GetKurtosis() << std::endl
821  << std::endl;
822  //
823  stream << boost::format( " Entries :\n | %=11s | %=11s | %=11s | %=11s | %=11s |" ) % "All" % "Underflow" %
824  "Overflow" % "#Equivalent" % "Integral"
825  << std::endl
826  << boost::format( " | %=11.5g | %=11.5g | %=11.5g | %=11.5g | %=11.5g |" ) % histo->GetEntries() %
827  histo->GetBinContent( 0 ) % histo->GetBinContent( histo->GetNbinsX() + 1 ) %
828  histo->GetEffectiveEntries() % histo->Integral()
829  << std::endl
830  << std::endl;
831  //
832  return dumpText( hist, width, height, errors, stream );
833 }
GAUDI_API std::string format(const char *,...)
MsgStream format utility "a la sprintf(...)".
Definition: MsgStream.cpp:120
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:639
T endl(T...args)
bool isFailure() const
Test for a status code of FAILURE.
Definition: StatusCode.h:86
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:28
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IBaseHistogram *  histo,
const std::string title = "" 
)

get the title

Definition at line 126 of file Fill.cpp.

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

get the title

Definition at line 133 of file Fill.cpp.

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

get the title

Definition at line 140 of file Fill.cpp.

141 {
142  return htitle_( histo, title );
143 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IHistogram2D *  histo,
const std::string title = "" 
)

get the title

Definition at line 147 of file Fill.cpp.

148 {
149  return htitle_( histo, title );
150 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IHistogram3D *  histo,
const std::string title = "" 
)

get the title

Definition at line 154 of file Fill.cpp.

155 {
156  return htitle_( histo, title );
157 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IProfile *  histo,
const std::string title = "" 
)

get the title

Definition at line 161 of file Fill.cpp.

162 {
163  return htitle_( histo, title );
164 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IProfile1D *  histo,
const std::string title = "" 
)

get the title

Definition at line 168 of file Fill.cpp.

169 {
170  return htitle_( histo, title );
171 }
std::string Gaudi::Utils::Histos::htitle ( const AIDA::IProfile2D *  histo,
const std::string title = "" 
)

get the title

Definition at line 175 of file Fill.cpp.

176 {
177  return htitle_( histo, title );
178 }
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 209 of file HistoTableFormat.cpp.

210 {
211  if ( 0 == aida ) {
212  return "";
213  } // RETURN
214  const DataObject* object = dynamic_cast<const DataObject*>( aida );
215  if ( 0 == object ) {
216  return "";
217  } // RETURN
218  IRegistry* registry = object->registry();
219  if ( 0 == registry ) {
220  return "";
221  } // RETURN
222  std::string _path = registry->identifier();
223  std::string::size_type n = _path.find( "/stat/" );
224  if ( 0 == n ) {
225  return std::string( _path, 6 );
226  } // RETURN
227  return _path; // RETURN
228 }
STL class.
virtual const id_type & identifier() const =0
Full identifier (or key)
The IRegistry represents the entry door to the environment any data object residing in a transient da...
Definition: IRegistry.h:22
T find(T...args)
A DataObject is the base class of any identifiable object on any data store.
Definition: DataObject.h:29
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 225 of file HistoTableFormat.h.

226  {
227  for ( ; first != last; ++first ) {
228  stream << format( *first, fmt ) << term;
229  } // print table rows
230  return stream; // RETURN
231  }
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 257 of file HistoTableFormat.h.

258  {
259  return printList( histos.begin(), histos.end(), fmt, stream, term );
260  }
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 307 of file HistoTableFormat.h.

309  {
310  for ( ; begin != end; ++begin ) {
311  stream << format( begin->second, // the histogram
312  begin->first, // the key
313  fmt1, fmt2 )
314  << term;
315  }
316  return stream;
317  }
auto begin(reverse_wrapper< T > &w)
Definition: reverse.h:58
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
auto end(reverse_wrapper< T > &w)
Definition: reverse.h:64
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 362 of file HistoTableFormat.h.

364  {
365  return printMap( histos.begin(), histos.end(), fmt1, fmt2, stream, term );
366  }
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 151 of file HistoLabels.cpp.

152  {
153  return setAxisLabels_<TH1D>( hist, xAxis, yAxis );
154  }
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 156 of file HistoLabels.cpp.

157  {
158  return setAxisLabels_<TProfile>( hist, xAxis, yAxis );
159  }
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 161 of file HistoLabels.cpp.

162  {
163  return setAxisLabels_<TH2D>( hist, xAxis, yAxis );
164  }
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 166 of file HistoLabels.cpp.

167  {
168  return setAxisLabels_<TProfile2D>( hist, xAxis, yAxis );
169  }
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 97 of file HistoLabels.cpp.

97 { 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 87 of file HistoLabels.cpp.

88  {
89  return setBinLabels_<TH1D>( hist, labels );
90  }
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 99 of file HistoLabels.cpp.

99 { 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 92 of file HistoLabels.cpp.

93  {
94  return setBinLabels_<TProfile>( hist, labels );
95  }
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 101 of file HistoLabels.cpp.

102  {
103  if ( !hist ) return false;
104  TH2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
105  if ( !h2d ) return false;
106  BinLabels lx;
107  lx.reserve( xlabels.size() );
108  for ( unsigned int i = 0; i < xlabels.size(); ++i ) {
109  lx.emplace_back( i, xlabels[i] );
110  }
111  BinLabels ly;
112  ly.reserve( ylabels.size() );
113  for ( unsigned int i = 0; i < ylabels.size(); ++i ) {
114  ly.emplace_back( i, ylabels[i] );
115  }
116  return ( setBinLabels_( h2d->GetXaxis(), lx ) && setBinLabels_( h2d->GetYaxis(), ly ) );
117  }
std::vector< BinLabel > BinLabels
Typedef for a list of bin numbers and their associated label.
Definition: HistoLabels.h:31
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 119 of file HistoLabels.cpp.

120  {
121  TH2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
122  return ( h2d && setBinLabels_( h2d->GetXaxis(), xlabels ) && setBinLabels_( h2d->GetYaxis(), ylabels ) );
123  }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 125 of file HistoLabels.cpp.

126  {
127  if ( !hist ) return false;
128  TProfile2D* h2d = Gaudi::Utils::Aida2ROOT::aida2root( hist );
129  if ( !h2d ) return false;
130  BinLabels lx;
131  lx.reserve( xlabels.size() );
132  for ( unsigned int i = 0; i < xlabels.size(); ++i ) {
133  lx.emplace_back( i, xlabels[i] );
134  }
135  BinLabels ly;
136  ly.reserve( ylabels.size() );
137  for ( unsigned int i = 0; i < ylabels.size(); ++i ) {
138  ly.emplace_back( i, ylabels[i] );
139  }
140  return ( setBinLabels_( h2d->GetXaxis(), lx ) && setBinLabels_( h2d->GetYaxis(), ly ) );
141  }
std::vector< BinLabel > BinLabels
Typedef for a list of bin numbers and their associated label.
Definition: HistoLabels.h:31
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
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 143 of file HistoLabels.cpp.

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

Definition at line 180 of file Fill.cpp.

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

Definition at line 182 of file Fill.cpp.

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

Definition at line 184 of file Fill.cpp.

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

Definition at line 186 of file Fill.cpp.

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

Definition at line 188 of file Fill.cpp.

188 { 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 66 of file HistoXML.cpp.

67 {
68  return stream << TBufferXML::ConvertToXML( &histo );
69 }
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 76 of file HistoXML.cpp.

77 {
78  return stream << TBufferXML::ConvertToXML( &histo );
79 }
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 86 of file HistoXML.cpp.

87 {
88  return stream << TBufferXML::ConvertToXML( &histo );
89 }
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 126 of file HistoXML.cpp.

127 {
128  return stream << TBufferXML::ConvertToXML( &histo );
129 }
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 136 of file HistoXML.cpp.

137 {
138  return stream << TBufferXML::ConvertToXML( &histo );
139 }
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 146 of file HistoXML.cpp.

147 {
148  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
149  return root ? toXml( *root, stream ) : stream;
150 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:66
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 157 of file HistoXML.cpp.

158 {
159  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
160  return root ? toXml( *root, stream ) : stream;
161 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:66
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 168 of file HistoXML.cpp.

169 {
170  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
171  return root ? toXml( *root, stream ) : stream;
172 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:66
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 179 of file HistoXML.cpp.

180 {
181  auto root = Gaudi::Utils::Aida2ROOT::aida2root( &histo );
182  return root ? toXml( *root, stream ) : stream;
183 }
static TH1D * aida2root(AIDA::IHistogram1D *aida)
get the underlying pointer for 1D-histogram
Definition: Aida2ROOT.cpp:73
GAUDI_API std::ostream & toXml(const TH1D &histo, std::ostream &stream)
stream the ROOT histogram into output stream as XML
Definition: HistoXML.cpp:66
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 96 of file HistoXML.cpp.

97 {
98  return stream << TBufferXML::ConvertToXML( &histo );
99 }
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 106 of file HistoXML.cpp.

107 {
108  return stream << TBufferXML::ConvertToXML( &histo );
109 }
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 116 of file HistoXML.cpp.

117 {
118  return stream << TBufferXML::ConvertToXML( &histo );
119 }
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 190 of file HistoXML.cpp.

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