00001
00002 #ifndef HISTOGRAMSVC_GENHISTO2D_H
00003 #define HISTOGRAMSVC_GENHISTO2D_H 1
00004
00005
00006
00007 #include <string>
00008 #include "GaudiKernel/xtoa.h"
00009 #include "GaudiKernel/Kernel.h"
00010
00011 #include "GaudiKernel/StatusCode.h"
00012 #include "GaudiKernel/DataObject.h"
00013 #include "GaudiKernel/IHistogram2D.h"
00014 #include "Axis.h"
00015
00016 #include "HTL/Histograms.h"
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049 template <class TYPE>
00050 class GenHisto2D : virtual public IHistogram2D, public DataObject {
00051
00052 protected:
00053
00057
00058 GenHisto2D(
00059 const std::string& title,
00060 const std::string& histoID,
00061 int nBinsX,
00062 double lowX,
00063 double highX,
00064 int nBinsY,
00065 double lowY,
00066 double highY,
00067 End_Point_Convention epc )
00068 {
00069 std::string temp = histoID;
00070 temp += '|';
00071 temp += title;
00072 m_histogram = new TYPE( temp.c_str(), nBinsX, lowX, highX,
00073 nBinsY, lowY, highY, epc );
00074
00075 m_xAxis = new Axis( m_histogram, 0, nBinsX, lowX, highX );
00076
00077 m_yAxis = new Axis( m_histogram, 1, nBinsY, lowY, highY );
00078 }
00079
00080 GenHisto2D(
00081 const std::string& title,
00082 int histoID,
00083 int nBinsX,
00084 double lowX,
00085 double highX,
00086 int nBinsY,
00087 double lowY,
00088 double highY,
00089 End_Point_Convention epc )
00090 {
00091
00092 char buffer[32];
00093 _itoa(histoID,buffer,10);
00094 std::string temp = buffer;
00095 temp += '|';
00096 temp += title;
00097 m_histogram = new TYPE( temp.c_str(), nBinsX, lowX, highX,
00098 nBinsY, lowY, highY, epc );
00099
00100 m_xAxis = new Axis( m_histogram, 0, nBinsX, lowX, highX );
00101
00102 m_yAxis = new Axis( m_histogram, 1, nBinsY, lowY, highY );
00103 }
00104
00105
00109
00110 GenHisto2D(
00111 const std::string& title,
00112 const std::string& histoID,
00113 std::vector<double> edgesX,
00114 std::vector<double> edgesY,
00115 End_Point_Convention epc )
00116 {
00117 std::string temp = histoID;
00118 temp += '|';
00119 temp += title;
00120 m_histogram = new TYPE( temp.c_str(), edgesX, edgesY );
00121
00122 m_xAxis = new Axis( m_histogram, 0, edgesX );
00123
00124 m_yAxis = new Axis( m_histogram, 1, edgesY );
00125 }
00126
00127 GenHisto2D(
00128 const std::string& title,
00129 int histoID,
00130 std::vector<double> edgesX,
00131 std::vector<double> edgesY,
00132 End_Point_Convention epc )
00133 {
00134
00135 char buffer[32];
00136 _itoa(histoID,buffer,10);
00137 std::string temp = buffer;
00138 temp += '|';
00139 temp += title;
00140 m_histogram = new TYPE( temp.c_str(), edgesX, edgesY );
00141
00142 m_xAxis = new Axis( m_histogram, 0, edgesX );
00143
00144 m_yAxis = new Axis( m_histogram, 1, edgesY );
00145 }
00146
00147
00151
00152 virtual ~GenHisto2D() {
00153 delete m_histogram;
00154 delete m_xAxis;
00155 delete m_yAxis;
00156 }
00157
00158
00159 public:
00160
00161
00165
00166
00170
00174 virtual std::string title() const {
00175 return m_histogram->name();
00176 }
00177 virtual void setTitle( std::string value ) {
00178 m_histogram->set_name( value.c_str() );
00179 }
00180
00182
00183
00185 virtual int dimensions() const {
00186 return m_histogram->dim();
00187 }
00188
00190 virtual void reset() {
00191 m_histogram->reset();
00192 }
00193
00194
00198
00200 virtual int entries() const {
00201 return (int) HStat::in_range_entries_count( *m_histogram );
00202 }
00203
00206 virtual int allEntries() const {
00207 return (int) HStat::entries_count( *m_histogram );
00208 }
00209
00211 virtual int extraEntries() const {
00212 return (int) HStat::extra_entries_count( *m_histogram );
00213 }
00214
00217 virtual double equivalentBinEntries() const {
00218 return HStat::nequival( *m_histogram );
00219 }
00220
00221
00225
00227 virtual double sumBinHeights() const {
00228 return HInfo::in_range_value( *m_histogram );
00229 }
00230
00232 virtual double sumAllBinHeights() const {
00233 return sumBinHeights() + sumExtraBinHeights();
00234 }
00235
00237 virtual double sumExtraBinHeights() const {
00238 return binHeightX( UNDERFLOW_BIN )
00239 + binHeight ( UNDERFLOW_BIN, OVERFLOW_BIN )
00240 + binHeightY( OVERFLOW_BIN )
00241 + binHeight ( OVERFLOW_BIN, OVERFLOW_BIN )
00242 + binHeightX( OVERFLOW_BIN )
00243 + binHeight ( OVERFLOW_BIN, UNDERFLOW_BIN )
00244 + binHeightY( UNDERFLOW_BIN )
00245 + binHeight ( UNDERFLOW_BIN, UNDERFLOW_BIN );
00246 }
00247
00248
00252
00253
00257
00259 virtual void fill( double x, double y, double weight = 1 ) {
00260 m_histogram->fill( x, y, weight );
00261 }
00262
00263
00277
00279 virtual int binEntries( int indexX, int indexY ) const {
00280 int i = checkIndexX(indexX);
00281 int j = checkIndexY(indexY);
00282 return htlBin(i,j)->count();
00283 }
00284
00287 virtual int binEntriesX( int indexX ) const {
00288 int i = checkIndexX(indexX);
00289 if( UNDERFLOW_BIN == i ) {
00290 return m_histogram->extra_bin( H_UNDERFLOW, H_IN_RANGE ).count();
00291 }
00292 else if( OVERFLOW_BIN == i ) {
00293 return m_histogram->extra_bin( H_OVERFLOW , H_IN_RANGE ).count();
00294 }
00295
00296 I_Histo::I_Bin_Location loc(2);
00297 loc[0] = i;
00298 int sum = 0;
00299 for( int j=0; j<m_yAxis->bins(); j++ ) {
00300 loc[1] = j;
00301 sum += m_histogram->i_bin( loc ).count();
00302 }
00303 return sum;
00304 }
00305
00308 virtual int binEntriesY( int indexY ) const {
00309 int j = checkIndexY(indexY);
00310 if( UNDERFLOW_BIN == j ) {
00311 return m_histogram->extra_bin( H_IN_RANGE , H_UNDERFLOW ).count();
00312 }
00313 else if( OVERFLOW_BIN == j ) {
00314 return m_histogram->extra_bin( H_IN_RANGE , H_OVERFLOW ).count();
00315 }
00316
00317 I_Histo::I_Bin_Location loc(2);
00318 loc[1] = j;
00319 int sum = 0;
00320 for( int i=0; i<m_xAxis->bins(); i++ ) {
00321 loc[0] = i;
00322 sum += m_histogram->i_bin( loc ).count();
00323 }
00324 return sum;
00325 }
00326
00328 virtual double binHeight( int indexX, int indexY ) const {
00329 int i = checkIndexX(indexX);
00330 int j = checkIndexY(indexY);
00331 return htlBin(i,j)->value();
00332 }
00333
00336 virtual double binHeightX( int indexX ) const {
00337 int i = checkIndexX(indexX);
00338 if( UNDERFLOW_BIN == i ) {
00339 return m_histogram->extra_bin( H_UNDERFLOW, H_IN_RANGE ).value();
00340 }
00341 else if( OVERFLOW_BIN == i ) {
00342 return m_histogram->extra_bin( H_OVERFLOW , H_IN_RANGE ).value();
00343 }
00344
00345 I_Histo::I_Bin_Location loc(2);
00346 loc[0] = i;
00347 double sum = 0.;
00348 for( int j=0; j<m_yAxis->bins(); j++ ) {
00349 loc[1] = j;
00350 sum += m_histogram->i_bin( loc ).value();
00351 }
00352 return sum;
00353 }
00354
00357 virtual double binHeightY( int indexY ) const {
00358 int j = checkIndexY(indexY);
00359 if( UNDERFLOW_BIN == j ) {
00360 return m_histogram->extra_bin( H_IN_RANGE , H_UNDERFLOW ).value();
00361 }
00362 else if( OVERFLOW_BIN == j ) {
00363 return m_histogram->extra_bin( H_IN_RANGE , H_OVERFLOW ).value();
00364 }
00365
00366 I_Histo::I_Bin_Location loc(2);
00367 loc[1] = j;
00368 double sum = 0.;
00369 for( int i=0; i<m_xAxis->bins(); i++ ) {
00370 loc[0] = i;
00371 sum += m_histogram->i_bin( loc ).value();
00372 }
00373 return sum;
00374 }
00375
00377 virtual double binError( int indexX, int indexY ) const {
00378 int i = checkIndexX(indexX);
00379 int j = checkIndexY(indexY);
00380 return htlBin(i,j)->error();
00381 }
00382
00383 protected:
00384 virtual I_Bin* htlBin( int i, int j ) const {
00385 if( UNDERFLOW_BIN != i && OVERFLOW_BIN != i
00386 && UNDERFLOW_BIN != j && OVERFLOW_BIN != j ) {
00387 I_Histo::I_Bin_Location loc(2); loc[0] = i; loc[1] = j;
00388 return &m_histogram->i_bin( loc );
00389 }
00390 else if( UNDERFLOW_BIN == i && UNDERFLOW_BIN != j && OVERFLOW_BIN != j ) {
00391 return &m_histogram->extra_bin( H_UNDERFLOW, H_IN_RANGE );
00392 }
00393 else if( UNDERFLOW_BIN == i && OVERFLOW_BIN == j ) {
00394 return &m_histogram->extra_bin( H_UNDERFLOW, H_OVERFLOW );
00395 }
00396 else if( UNDERFLOW_BIN != i && OVERFLOW_BIN != i && OVERFLOW_BIN == j ) {
00397 return &m_histogram->extra_bin( H_IN_RANGE , H_OVERFLOW );
00398 }
00399 else if( OVERFLOW_BIN == i && OVERFLOW_BIN == j ) {
00400 return &m_histogram->extra_bin( H_OVERFLOW , H_OVERFLOW );
00401 }
00402 else if( OVERFLOW_BIN == i && UNDERFLOW_BIN != j && OVERFLOW_BIN != j ) {
00403 return &m_histogram->extra_bin( H_OVERFLOW , H_IN_RANGE );
00404 }
00405 else if( OVERFLOW_BIN == i && UNDERFLOW_BIN == j ) {
00406 return &m_histogram->extra_bin( H_OVERFLOW , H_UNDERFLOW );
00407 }
00408 else if( UNDERFLOW_BIN != i && OVERFLOW_BIN != i && UNDERFLOW_BIN == j ) {
00409 return &m_histogram->extra_bin( H_IN_RANGE , H_UNDERFLOW );
00410 }
00411
00412 return &m_histogram->extra_bin( H_UNDERFLOW, H_UNDERFLOW );
00413 }
00414
00415 public:
00416
00417
00421
00424 virtual double meanX() const {
00425 std::cerr << "Sorry, not yet implemented" << std::endl;
00426 return 0.;
00427 }
00430 virtual double meanY() const {
00431 std::cerr << "Sorry, not yet implemented" << std::endl;
00432 return 0.;
00433 }
00434
00437 virtual double rmsX() const {
00438 std::cerr << "Sorry, not yet implemented" << std::endl;
00439 return 0.;
00440 }
00443 virtual double rmsY() const {
00444 std::cerr << "Sorry, not yet implemented" << std::endl;
00445 return 0.;
00446 }
00447
00448
00452
00454 virtual double minBinHeight() const {
00455 return HInfo::in_range_min_value( *m_histogram );
00456 }
00457
00459 virtual int minBinX() const {
00460 I_Histo::I_Bin_Location minLoc(2);
00461 minLoc = htlMinBin();
00462 return minLoc[0];
00463 }
00464
00466 virtual int minBinY() const {
00467 I_Histo::I_Bin_Location minLoc(2);
00468 minLoc = htlMinBin();
00469 return minLoc[1];
00470 }
00471
00473 virtual double maxBinHeight() const {
00474 return HInfo::in_range_max_value( *m_histogram );
00475 }
00476
00478 virtual int maxBinX() const {
00479 I_Histo::I_Bin_Location maxLoc(2);
00480 maxLoc = htlMinBin();
00481 return maxLoc[0];
00482 }
00483
00485 virtual int maxBinY() const {
00486 I_Histo::I_Bin_Location maxLoc(2);
00487 maxLoc = htlMinBin();
00488 return maxLoc[1];
00489 }
00490
00491 protected:
00492
00494 virtual I_Histo::I_Bin_Location htlMinBin() const {
00495 I_Histo::I_Bin_Location minLoc(2);
00496 I_Histo::I_Bin_Location loc(2);
00497 loc[0] = 0;
00498 loc[1] = 0;
00499 double val = m_histogram->i_bin( loc ).value();
00500 for( int i=1; i<m_xAxis->bins(); i++ ) {
00501 for( int j=1; j<m_yAxis->bins(); j++ ) {
00502 loc[0] = i;
00503 loc[1] = j;
00504 if( m_histogram->i_bin( loc ).value() < val ) {
00505 minLoc[0] = i;
00506 minLoc[1] = j;
00507 }
00508 }
00509 }
00510 return minLoc;
00511 }
00512
00514 virtual I_Histo::I_Bin_Location htlMaxBin() const {
00515 I_Histo::I_Bin_Location maxLoc(2);
00516 I_Histo::I_Bin_Location loc(2);
00517 loc[0] = 0;
00518 loc[1] = 0;
00519 double val = m_histogram->i_bin( loc ).value();
00520 for( int i=1; i<m_xAxis->bins(); i++ ) {
00521 for( int j=1; j<m_yAxis->bins(); j++ ) {
00522 loc[0] = i;
00523 loc[1] = j;
00524 if( m_histogram->i_bin( loc ).value() > val ) {
00525 maxLoc[0] = i;
00526 maxLoc[1] = j;
00527 }
00528 }
00529 }
00530 return maxLoc;
00531 }
00532
00533
00534 public:
00535
00539
00541 virtual IAxis* xAxis() const {
00542 return m_xAxis;
00543 }
00544
00546 virtual IAxis* yAxis() const {
00547 return m_yAxis;
00548 }
00549
00550
00554
00556 virtual int coordToIndexX( double coordX ) const {
00557 return xAxis()->coordToIndex( coordX );
00558 }
00559
00561 virtual int coordToIndexY( double coordY ) const {
00562 return yAxis()->coordToIndex( coordY );
00563 }
00564
00565
00569
00572 virtual IHistogram1D* projectionX() const {
00573 HSlicer s;
00574
00575 std::string newTitle = title();
00576 int sep = newTitle.find('|');
00577 if ( 0 < sep ) {
00578
00579 std::string temp( newTitle, sep+1, newTitle.length() );
00580 newTitle = temp;
00581 }
00582
00583 newTitle = newTitle + " ( Projection on the axis X )";
00584
00585 IHistogram1D* proj =
00586 new H1DVar( newTitle, *s.xProject( *m_histogram ), xAxis()->edges() );
00587 return proj;
00588 }
00589
00592 virtual IHistogram1D* projectionY() const {
00593 HSlicer s;
00594
00595 std::string newTitle = title();
00596 int sep = newTitle.find('|');
00597 if ( 0 < sep ) {
00598
00599 std::string temp( newTitle, sep+1, newTitle.length() );
00600 newTitle = temp;
00601 }
00602
00603 newTitle = newTitle + " ( Projection on the axis Y )";
00604
00605 IHistogram1D* proj =
00606 new H1DVar( newTitle, *s.yProject( *m_histogram ), yAxis()->edges() );
00607 return proj;
00608 }
00609
00610
00614
00617 virtual IHistogram1D* sliceX( int indexY ) const {
00618 HSlicer s;
00619
00620 char buffer[32];
00621 _itoa( indexY, buffer, 10 );
00622 std::string indexAsString = buffer;
00623
00624 std::string newTitle = title();
00625 int sep = newTitle.find('|');
00626 if ( 0 < sep ) {
00627
00628 std::string temp( newTitle, sep+1, newTitle.length() );
00629 newTitle = temp;
00630 }
00631
00632 newTitle = newTitle + " ( Slice parallel with the axis X, at the index "
00633 + indexAsString
00634 + " )";
00635
00636 IHistogram1D* proj =
00637 new H1DVar( newTitle,
00638 *s.xBand( *m_histogram, indexY ),
00639 xAxis()->edges() );
00640 return proj;
00641 }
00642
00645 virtual IHistogram1D* sliceY( int indexX ) const {
00646 HSlicer s;
00647
00648 char buffer[32];
00649 _itoa( indexX, buffer, 10 );
00650 std::string indexAsString = buffer;
00651
00652 std::string newTitle = title();
00653 int sep = newTitle.find('|');
00654 if ( 0 < sep ) {
00655
00656 std::string temp( newTitle, sep+1, newTitle.length() );
00657 newTitle = temp;
00658 }
00659
00660 newTitle = newTitle + " ( Slice parallel with the axis Y, at the index "
00661 + indexAsString
00662 + " )";
00663
00664 IHistogram1D* proj =
00665 new H1DVar( newTitle,
00666 *s.yBand( *m_histogram, indexX ),
00667 yAxis()->edges() );
00668 return proj;
00669 }
00670
00672 virtual IHistogram1D* sliceX( int , int ) const {
00673 std::cerr << "Sorry, not yet implemented" << std::endl;
00674 return (IHistogram1D*) 0;
00675 }
00676
00678 virtual IHistogram1D* sliceY( int , int ) const {
00679 std::cerr << "Sorry, not yet implemented" << std::endl;
00680 return (IHistogram1D*) 0;
00681 }
00682
00683
00687
00690
00692 virtual std::ostream& print( std::ostream& s ) const {
00693 HPrinter hp( s );
00694 hp.print( *m_histogram );
00695 return s;
00696 }
00697
00700
00702 virtual std::ostream& write( std::ostream& s ) const {
00703 s << "2-dimensional histogram" << std::endl;
00704 for( int i = 0; i < m_xAxis->bins(); i++ ) {
00705 for( int j = 0; j < m_yAxis->bins(); j++ ) {
00706 s << m_xAxis->binCentre(i) << " "
00707 << m_yAxis->binCentre(j) << " "
00708 << binHeight(i,j) << " "
00709 << binError(i,j) << std::endl;
00710 }
00711 }
00712 s << std::endl;
00713 return s;
00714 }
00715
00717 virtual int write( const char* file_name ) const {
00718 HistoTable2D ht( file_name );
00719 int status = ht.write( *m_histogram );
00720 return status;
00721 }
00722
00723
00724 protected:
00725
00727 virtual int calcIndex( double coord, int part ) const {
00728 Index an_index;
00729 Extra_Index an_extra_index;
00730 m_histogram->i_partition(part).map_point( coord, an_index, an_extra_index );
00731 if (an_extra_index == H_IN_RANGE) {
00732 return an_index;
00733 }
00734 return 0;
00735 }
00736
00740
00743 virtual int checkIndexX( int indexX ) const {
00744 return m_xAxis->checkIndex( indexX );
00745 }
00746
00749 virtual int checkIndexY( int indexY ) const {
00750 return m_yAxis->checkIndex( indexY );
00751 }
00752
00753
00757
00759 TYPE* m_histogram;
00760
00762 Axis* m_xAxis;
00763
00765 Axis* m_yAxis;
00766
00767 };
00768
00769
00770 #endif // HISTOGRAMSVC_GENHISTO2D_H