Gaudi Framework, version v20r4

Generated: 8 Jan 2009

HbookCnv::H2DCnv Class Reference

#include <H2DCnv.h>

Inheritance diagram for HbookCnv::H2DCnv:

Inheritance graph
[legend]
Collaboration diagram for HbookCnv::H2DCnv:

Collaboration graph
[legend]

List of all members.


Detailed Description

Definition at line 30 of file H2DCnv.h.


Public Member Functions

virtual StatusCode createObj (IOpaqueAddress *pAddr, DataObject *&pObj)
 Create the transient representation of an object.
virtual StatusCode updateObj (IOpaqueAddress *pAddr, DataObject *pObj)
 Update the transient object from the other representation.
virtual StatusCode updateRep (IOpaqueAddress *pAddr, DataObject *pObj)
 Update the converted representation of a transient object.
 H2DCnv (ISvcLocator *svc)
 Standard constructor.
virtual ~H2DCnv ()
 Standard destructor.

Static Public Member Functions

static const CLIDclassID ()
 Inquire class type.

Protected Member Functions

virtual StatusCode book (IOpaqueAddress *pAddr, DataObject *pObj)
 Book 2 D histogram.

Friends

class CnvFactory< H2DCnv >

Constructor & Destructor Documentation

HbookCnv::H2DCnv::H2DCnv ( ISvcLocator svc  )  [inline]

Standard constructor.

Definition at line 44 of file H2DCnv.h.

00044                                : HConverter( classID(), svc )   { 
00045       m_prefix = "/stat";
00046       m_deleteAfterSave = true;
00047     }

virtual HbookCnv::H2DCnv::~H2DCnv (  )  [inline, virtual]

Standard destructor.

Definition at line 49 of file H2DCnv.h.

00049                          { 
00050     }


Member Function Documentation

static const CLID& HbookCnv::H2DCnv::classID (  )  [inline, static]

Inquire class type.

Definition at line 34 of file H2DCnv.h.

00034                                    {
00035       return CLID_H2D;
00036     }

StatusCode HbookCnv::H2DCnv::createObj ( IOpaqueAddress pAddr,
DataObject *&  pObj 
) [virtual]

Create the transient representation of an object.

Create 2 D histogram from disk representation.

Reimplemented from Converter.

Definition at line 195 of file H2DCnv.cpp.

00197 {
00198   //MsgStream log(msgSvc(), "HHistoCnv");
00199   IHistogramFactory* pFac = dynamic_cast<IHistogramFactory*>(dataProvider());
00200   if ( 0 != pAddr && pFac ) {
00201     std::string loc = pAddr->registry()->identifier();
00202     long id = pAddr->ipar()[0];
00203     if ( readObject(loc, id).isSuccess() )  {
00204       std::string title;
00205       int nx, ny, nwt, lpaw;
00206       float xmi, xma, ymi, yma;
00207       ::HGIVE(id, title, nx, xmi, xma, ny, ymi, yma, nwt, lpaw);
00208       pAddr->addRef();
00209       refpObj = dynamic_cast<DataObject*>(pFac->createHistogram2D(loc, title,
00210                                                                   nx, xmi, xma,
00211                                                                   ny, ymi, yma));
00212       refpObj->registry()->setAddress(pAddr);
00213       StatusCode sc = updateObj(pAddr, refpObj);
00214       pAddr->release();
00215       return sc;
00216     }
00217   }
00218   return StatusCode::FAILURE;
00219 }

StatusCode HbookCnv::H2DCnv::updateObj ( IOpaqueAddress pAddr,
DataObject pObj 
) [virtual]

Update the transient object from the other representation.

Convert the transient object to the requested representation.

Reimplemented from Converter.

Definition at line 222 of file H2DCnv.cpp.

00223                                                           {
00224   IHistogram2D* h = dynamic_cast<IHistogram2D*>(pObj);
00225 
00226   if ( 0 != h && 0 != pAddr ) {
00227     int    id  = pAddr->ipar()[0]; // Histogram ID
00228     std::string loc = pAddr->registry()->identifier();
00229     if ( !::HEXIST( id ) )  {
00230       readObject(loc, id);
00231     }
00232     if ( ::HEXIST( id ) )  {
00233       const  IAxis& x = h->xAxis();
00234       int    nBinX    = x.bins();         // Number of bins
00235       double xLow     = x.lowerEdge();    // Histogram lower edge
00236       double xHigh    = x.upperEdge();    // Histogram upper edge
00237       double x_uflow  = xLow  - 1.0;      // Underflow coordinate
00238       double x_oflow  = xHigh + 1.0;      // Overflow coordinate
00239 
00240       const  IAxis& y = h->yAxis();
00241       int    nBinY    = y.bins();         // Number of bins in the axis Y
00242       double yLow     = y.lowerEdge();    // Y-lower edge
00243       double yHigh    = y.upperEdge();    // Y-upper edge
00244       double y_uflow  = yLow  - 1.0;      // Underflow coordinate Y
00245       double y_oflow  = yHigh + 1.0;      // Overflow coordinate Y
00246       // float x_mean  = xLow + (xHigh-xLow)/2.0;
00247       // float y_mean  = yLow + (yHigh-yLow)/2.0;
00248       double centreX, centreY, height, err;
00249       int ix, iy;
00250 
00251       h->reset();
00252 
00253       for ( ix = 0; ix < nBinX; ix++ )  {
00254         for ( iy = 0; iy < nBinY; iy++ )  {
00255           centreX = h->binMeanX(ix,iy);
00256           centreY = h->binMeanY(ix,iy);
00257           height  = ::HIJ(id, ix+1, iy+1);
00258           err     = ::HIJE(id, ix+1, iy+1);
00259           // Unfortunately IHistogram2D does not support to set bins
00260           // and bin errors....errors will be SQRT(height)
00261           h->fill(centreX, centreY, height);
00262         }
00263       }
00264 
00265       // Underflow and overflow bins
00266       // W
00267       for ( iy = 0; iy < nBinY; iy++ )  {
00268         centreY = h->binMeanY(IAxis::UNDERFLOW_BIN,iy);
00269         height  = ::HIJ(id, 0, iy+1);
00270         h->fill(x_uflow, centreY, height);
00271       }
00272       // NW
00273       height  = ::HIJ(id, 0, nBinY+1);
00274       h->fill(x_uflow, y_oflow, height);
00275       // N
00276       for ( ix = 0; ix < nBinX; ix++ )  {
00277         double centreX = h->binMeanX(ix,IAxis::OVERFLOW_BIN);
00278         double height  = ::HIJ(id, ix+1, nBinY+1);
00279         h->fill(centreX, y_oflow, height);
00280       }
00281       // NE
00282       height  = ::HIJ(id, nBinX+1, nBinY+1);
00283       h->fill(x_oflow, y_oflow, height);
00284       // E
00285       for ( iy = 0; iy < nBinY; iy++ )  {
00286         double centreY = h->binMeanY(IAxis::OVERFLOW_BIN,iy);
00287         double height  = ::HIJ(id, nBinX+1, iy+1);
00288         h->fill(x_oflow, centreY, height);
00289       }
00290       // SE
00291       height  = ::HIJ(id, nBinX+1, 0);
00292       h->fill(x_oflow, y_uflow, height);
00293       // S
00294       for ( ix = 0; ix < nBinX; ix++ )  {
00295         double centreX = h->binMeanX(ix,IAxis::UNDERFLOW_BIN);
00296         double height  = ::HIJ(id, ix+1, 0);
00297         h->fill(centreX, y_uflow, height);
00298       }
00299       // SW
00300       height  = ::HIJ(id, 0, 0);
00301       h->fill(x_uflow, y_uflow, height);
00302 
00303       ::HDELET(id);
00304       return StatusCode::SUCCESS;
00305     }
00306   }
00307   return StatusCode::FAILURE;
00308 }

StatusCode HbookCnv::H2DCnv::updateRep ( IOpaqueAddress pAddr,
DataObject pObj 
) [virtual]

Update the converted representation of a transient object.

Convert the transient object to the requested representation.

Reimplemented from Converter.

Definition at line 61 of file H2DCnv.cpp.

00062                                                           {
00063   IHistogram2D* h = dynamic_cast<IHistogram2D*>(pObj);
00064   if ( 0 != h ) {
00065     const IAxis& x = h->xAxis();
00066     int    histoID  = pAddr->ipar()[0]; // Histogram ID
00067     int    nBinX    = x.bins();         // Number of bins in the axis X
00068     double xLow     = x.lowerEdge();    // X-lower edge
00069     double xHigh    = x.upperEdge();    // X-upper edge
00070     double x_uflow  = xLow  - 1.0;      // Underflow coordinate X
00071     double x_oflow  = xHigh + 1.0;      // Overflow coordinate X
00072 
00073     const IAxis& y = h->yAxis();
00074     int    nBinY    = y.bins();         // Number of bins in the axis Y
00075     double yLow     = y.lowerEdge();    // Y-lower edge
00076     double yHigh    = y.upperEdge();    // Y-upper edge
00077     double y_uflow  = yLow  - 1.0;      // Underflow coordinate Y
00078     double y_oflow  = yHigh + 1.0;      // Overflow coordinate Y
00079 
00080     int    entries  = 0;                // Number of entries in a bin
00081     double height   = 0.;               // Height if a bin
00082     double centreX  = 0.;               // Centre of a bin on the axis X
00083     double centreY  = 0.;               // Centre of a bin on the axis Y
00084     double weight   = 0.;               // height / entries
00085     int    ix       = 0;                // Bin index on the axis X
00086     int    iy       = 0;                // Bin index on the axis Y
00087     int    j        = 0;                // Fill index
00088 
00089     double* errs     = new double[nBinX*nBinY];  // Array of bin errors
00090     int    ind      = 0;                       // Index of bin errors
00091 
00092     double x_mean  = xLow + (xHigh-xLow)/2.0;
00093     double y_mean  = yLow + (yHigh-yLow)/2.0;
00094 
00095     ::HRESET(histoID, " ", 1);
00096 
00097     // Fill the HBOOK 2D histogram with fixed binning
00098     for ( ix = 0; ix < nBinX; ix++ )  {
00099       for ( iy = 0; iy < nBinY; iy++, ind++ )  {
00100         errs[ind] = h->binError(ix,iy);
00101         height = h->binHeight(ix,iy);
00102         if( 0 != height ) {
00103           centreX = h->binMeanX( ix, iy );
00104           centreY = h->binMeanY( ix, iy );
00105           entries = h->binEntries(ix,iy);
00106           weight  = height / entries;
00107           for( j = 0; j < entries; j++ ) {
00108             ::HFILL(histoID, (float)centreX, (float)centreY, (float)weight);
00109           }
00110         }
00111       }
00112     }
00113     ::HPAKE(histoID,errs);
00114     delete [] errs;
00115 
00116     // Underflow and overflow bins
00117     // W
00118     entries          = h->binEntriesX( IAxis::UNDERFLOW_BIN ) - 
00119                        h->binEntries ( IAxis::UNDERFLOW_BIN,  IAxis::OVERFLOW_BIN ) - 
00120                        h->binEntries ( IAxis::UNDERFLOW_BIN, IAxis::UNDERFLOW_BIN );
00121     height           = h->binHeightX ( IAxis::UNDERFLOW_BIN ) - 
00122                        h->binHeight  ( IAxis::UNDERFLOW_BIN, IAxis::OVERFLOW_BIN ) -
00123                        h->binHeight  ( IAxis::UNDERFLOW_BIN, IAxis::UNDERFLOW_BIN );
00124     weight  = height / entries;
00125     for( j = 0; j < entries; j++ ) {
00126       ::HFILL(histoID, (float)x_uflow, (float)y_mean, (float)weight);
00127     }
00128     // NW
00129     entries          = h->binEntries ( IAxis::UNDERFLOW_BIN, IAxis::OVERFLOW_BIN  );
00130     height           = h->binHeight  ( IAxis::UNDERFLOW_BIN, IAxis::OVERFLOW_BIN  );
00131     weight  = height / entries;
00132     for( j = 0; j < entries; j++ ) {
00133       ::HFILL(histoID, (float)x_uflow, (float)y_oflow, (float)weight);
00134     }
00135     // N
00136     entries          = h->binEntriesY( IAxis::OVERFLOW_BIN  ) - 
00137                        h->binEntries ( IAxis::UNDERFLOW_BIN, IAxis::OVERFLOW_BIN  ) -
00138                        h->binEntries ( IAxis::OVERFLOW_BIN,  IAxis::OVERFLOW_BIN  );
00139     height           = h->binHeightY ( IAxis::OVERFLOW_BIN  ) -
00140                        h->binHeight  ( IAxis::UNDERFLOW_BIN, IAxis::OVERFLOW_BIN  ) -
00141                        h->binHeight  ( IAxis::OVERFLOW_BIN,  IAxis::OVERFLOW_BIN  );
00142     weight  = height / entries;
00143     for( j = 0; j < entries; j++ ) {
00144       ::HFILL(histoID, (float)x_mean, (float)y_oflow, (float)weight);
00145     }
00146     // NE
00147     entries          = h->binEntries ( IAxis::OVERFLOW_BIN,  IAxis::OVERFLOW_BIN  );
00148     height           = h->binHeight  ( IAxis::OVERFLOW_BIN,  IAxis::OVERFLOW_BIN  );
00149     weight  = height / entries;
00150     for( j = 0; j < entries; j++ ) {
00151       ::HFILL(histoID, (float)x_oflow, (float)y_oflow, (float)weight);
00152     }
00153     // E
00154     entries          = h->binEntriesX( IAxis::OVERFLOW_BIN  ) - 
00155                        h->binEntries ( IAxis::OVERFLOW_BIN, IAxis::UNDERFLOW_BIN  ) -
00156                        h->binEntries ( IAxis::OVERFLOW_BIN, IAxis::OVERFLOW_BIN  );
00157     height           = h->binHeightX ( IAxis::OVERFLOW_BIN  ) -
00158                        h->binHeight  ( IAxis::OVERFLOW_BIN, IAxis::UNDERFLOW_BIN  ) -
00159                        h->binHeight  ( IAxis::OVERFLOW_BIN, IAxis::OVERFLOW_BIN  );
00160     weight  = height / entries;
00161     for( j = 0; j < entries; j++ ) {
00162       ::HFILL(histoID, (float)x_oflow, (float)y_mean, (float)weight);
00163     }
00164     // SE
00165     entries          = h->binEntries ( IAxis::OVERFLOW_BIN,  IAxis::UNDERFLOW_BIN );
00166     height           = h->binHeight  ( IAxis::OVERFLOW_BIN,  IAxis::UNDERFLOW_BIN );
00167     weight  = height / entries;
00168     for( j = 0; j < entries; j++ ) {
00169       ::HFILL(histoID, (float)x_oflow, (float)y_uflow, (float)weight);
00170     }
00171     // S
00172     entries          = h->binEntriesY( IAxis::UNDERFLOW_BIN  ) - 
00173                        h->binEntries ( IAxis::OVERFLOW_BIN,  IAxis::UNDERFLOW_BIN  ) -
00174                        h->binEntries ( IAxis::UNDERFLOW_BIN, IAxis::UNDERFLOW_BIN  );
00175     height           = h->binHeightY ( IAxis::UNDERFLOW_BIN  ) -
00176                        h->binHeight  ( IAxis::OVERFLOW_BIN,  IAxis::UNDERFLOW_BIN  ) -
00177                        h->binHeight  ( IAxis::UNDERFLOW_BIN, IAxis::UNDERFLOW_BIN  );
00178     weight  = height / entries;
00179     for( j = 0; j < entries; j++ ) {
00180       ::HFILL(histoID, (float)x_mean, (float)y_uflow, (float)weight);
00181     }
00182     // SW
00183     entries          = h->binEntries ( IAxis::UNDERFLOW_BIN, IAxis::UNDERFLOW_BIN );
00184     height           = h->binHeight  ( IAxis::UNDERFLOW_BIN, IAxis::UNDERFLOW_BIN );
00185     weight  = height / entries;
00186     for( j = 0; j < entries; j++ ) {
00187       ::HFILL(histoID, (float)x_uflow, (float)y_uflow, (float)weight);
00188     }
00189     return StatusCode::SUCCESS;
00190   }
00191   return StatusCode::FAILURE;
00192 }

StatusCode HbookCnv::H2DCnv::book ( IOpaqueAddress pAddr,
DataObject pObj 
) [protected, virtual]

Book 2 D histogram.

Reimplemented from HbookCnv::HConverter.

Definition at line 36 of file H2DCnv.cpp.

00036                                                                      {
00037   IHistogram2D* h = dynamic_cast<IHistogram2D*>(pObj);
00038   if ( 0 != h ) {
00039     int idh = pAdd->ipar()[0];
00040     const IAxis& x = h->xAxis();
00041     const IAxis& y = h->yAxis();
00042     ::HBOOK2( idh, h->title(),
00043               x.bins(), float(x.lowerEdge()), float(x.upperEdge()),
00044               y.bins(), float(y.lowerEdge()), float(y.upperEdge()));
00045     ::HIDOPT(idh, "STAT");
00046 
00047     if ( !x.isFixedBinning() || !y.isFixedBinning() ) {
00048       MsgStream log( msgSvc(), "HbookCnv::H2DCnv" );
00049       log << MSG::WARNING 
00050           << "2D variable bin size histograms not supported in HBOOK" << endreq
00051           << "  -> Persistent 2D histogram " << idh << " '" + h->title() + "' INACCURATE!"
00052           << endreq;
00053     }
00054 
00055     return StatusCode::SUCCESS;
00056   }
00057   return StatusCode::FAILURE;
00058 }


Friends And Related Function Documentation

friend class CnvFactory< H2DCnv > [friend]

Definition at line 31 of file H2DCnv.h.


The documentation for this class was generated from the following files:

Generated at Thu Jan 8 17:54:06 2009 for Gaudi Framework, version v20r4 by Doxygen version 1.5.6 written by Dimitri van Heesch, © 1997-2004