TTable.cxx

Go to the documentation of this file.
00001 // @(#)root/table:$Id: TTable.cxx 36215 2010-10-09 09:01:54Z brun $
00002 // Author: Valery Fine(fine@bnl.gov)   03/07/98
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 
00013 ////////////////////////////////////////////////////////////////////////////
00014 //                                                                        //
00015 // TTable                                                                 //
00016 //                                                                        //
00017 // Wraps the array of the plain C-structures (one C-structure per element)//
00018 //                                                                        //
00019 // class TTable provides the automatic schema evolution for               //
00020 // the derived "table" classes saved with ROOT format.                    //
00021 //                                                                        //
00022 // "Automatic Schema evolution" provides:                                 //
00023 //   -  skipping data-member if it is not present for the current         //
00024 //      implementation of the "table" but was present at the time the     //
00025 //      table was written;                                                //
00026 //   -  assign a default value ZERO for the brand-new data-members,       //
00027 //      those were not in the structure when the object was written but   //
00028 //      present now;                                                      //
00029 //   -  trace propely any change in the order of the data-members         //
00030 //                                                                        //
00031 // To enjoy this class one has to derive his/her own custom class:        //
00032 //                                                                        //
00033 // St_dst_track_Table.h:                                                  //
00034 // ---------------------                                                  //
00035 //  #ifndef STAF_St_dst_track_Table                                       //
00036 //  #define STAF_St_dst_track_Table                                       //
00037 //                                                                        //
00038 //  #include "TTable.h"                                                   //
00039 //                                                                        //
00040 // // #include "dst_track.h"  the C-structure defintion may be kept       //
00041 //                            separately                                  //
00042 //  typedef struct dst_track_st {                                         //
00043 //     float r0;             /* radius at start (cm). See also comments*/
00044 //     float phi0;           /* azimuthal angle at start (deg)         */
00045 //     float z0;             /* z-coord. at start (cm)                 */
00046 //     float psi;            /* azimuthal angle of pT vector (deg)     */
00047 //     float tanl;           /* tan(dip) =pz/pt at start               */
00048 //     float invpt;          /* 1/pt at start (GeV/c)^(-1)             */
00049 //     float curvature;      /* Track curvature (1/cm)                 */
00050 //     float covar[15];      /* full covariance matrix                 */
00051 //     float chisq[2];       /* Chi-square per degree of freedom       */
00052 //     float x_first[3];     /* coord. of first measured point (cm)    */
00053 //     float x_last[3];      /* coord. of last measured point (cm)     */
00054 //     float length;         /* from first to last point (cm)          */
00055 //     float impact;         /* primary vertex (cm)                    */
00056 //     unsigned long map[2]; /* extrap. info. (see preceeding comments)*/
00057 //     int id;               /* Primary key (see comments)             */
00058 //     int iflag;            /* bitmask quality info. (see comments)   */
00059 //     int det_id;           /* Detector id information                */
00060 //     int method;           /* Track finding/fitting method, packed   */
00061 //     int pid;              /* Geant particle ID for assumed mass     */
00062 //     int n_point;          /* SVT, TPC, FTPC component #s are packed */
00063 //     int n_max_point;      /* SVT, TPC, FTPC component #s are packed */
00064 //     int n_fit_point;      /* SVT, TPC, FTPC component #s are packed */
00065 //     int icharge;          /* Particle charge in units of |e|        */
00066 //     int id_start_vertex;  /* final fit and primary track candidates */
00067 //  } DST_TRACK_ST;                                                       //
00068 //                                                                        //
00069 //  class St_dst_track : public TTable                                    //
00070 //  {                                                                     //
00071 //   public:                                                              //
00072 //     ClassDefTable(St_dst_track,dst_track_st)                           //
00073 //     ClassDef(St_dst_track,2) //C++ wrapper for <dst_track> StAF table  //
00074 //  };                                                                    //
00075 //  #endif                                                                //
00076 // ---------------------                                                  //
00077 //                                                                        //
00078 //  where the CPP macro defines several convinient methods for the        //
00079 //  "table" class (see: $ROOTSYS/table/inc/Ttypes.h for details:          //
00080 //                                                                        //
00081 //  #define ClassDefTable(className,structName)
00082 //    protected:
00083 //       static TTableDescriptor *fgColDescriptors;
00084 //       virtual TTableDescriptor *GetDescriptorPointer() const { return fgColDescriptors;}
00085 //       virtual void SetDescriptorPointer(TTableDescriptor *list) const { fgColDescriptors = list;}
00086 //    public:
00087 //      typedef structName* iterator;
00088 //      className() : TTable(_QUOTE_(className),sizeof(structName))    {SetType(_QUOTE_(structName));}
00089 //      className(const char *name) : TTable(name,sizeof(structName)) {SetType(_QUOTE_(structName));}
00090 //      className(Int_t n) : TTable(_QUOTE_(className),n,sizeof(structName)) {SetType(_QUOTE_(structName));}
00091 //      className(const char *name,Int_t n) : TTable(name,n,sizeof(structName)) {SetType(_QUOTE_(structName));}
00092 //      structName *GetTable(Int_t i=0) const { return ((structName *)GetArray())+i;}
00093 //      structName &operator[](Int_t i){ assert(i>=0 && i < GetNRows()); return *GetTable(i); }
00094 //      const structName &operator[](Int_t i) const { assert(i>=0 && i < GetNRows()); return *((const structName *)(GetTable(i))); }
00095 //      structName *begin() const  {                      return GetNRows()? GetTable(0):0;}
00096 //      structName *end()   const  {Int_t i = GetNRows(); return          i? GetTable(i):0;}
00097 //                                                                        //
00098 //  The class implementation file may 2 lines and look as follows:        //
00099 //  (for the example above):                                              //
00100 //                                                                        //
00101 //  St_dst_track_Table.cxx:                                               //
00102 //  -----------------------                                               //
00103 //       #include "St_dst_track_Table.h"                                  //
00104 //       TableClassImpl(St_dst_track, dst_track_st)                       //
00105 //  -----------------------                                               //
00106 //  LinkDef.h                                                             //
00107 //  -----------------------                                               //
00108 //  To provide ROOT I/O for this class TWO CINT dictonary entries         //
00109 //  should be defined with your custom LinkDef.h file                     //
00110 //     1. First entry (as usually) for the class derived from TTable      //
00111 //        for example:                                                    //
00112 // #pragma C++ class St_dst_track                                         //
00113 //     2. Second entry for the C-structure wrapped into the class.        //
00114 //         Since C-structuire is not derived from TObject it must be      //
00115 //         properly defined as "foreign" ROOT class                       //
00116 //    #pragma C++ class dst_track_st+;                                    //
00117 //  -----------------------                                               //
00118 // meta-variables i$ and n$ introduced                                    //
00119 // where "i$" stands for the current row index                            //
00120 //       "n$" stands for the total number of rows                         //
00121 // meta-variable can be used along the normal                             //
00122 // table column names in the expressions (see for example                 //
00123 // method TTable::Draw                                                    //
00124 //                                                                        //
00125 ////////////////////////////////////////////////////////////////////////////
00126 
00127 #include <assert.h>
00128 
00129 #ifdef WIN32
00130 # include <float.h>
00131 #endif
00132 
00133 //#if ROOT_VERSION_CODE >= ROOT_VERSION(3,03,5)
00134 #include "Riosfwd.h"
00135 #include "Riostream.h"
00136 //#include <iomanip.h>
00137 
00138 // #endif
00139 
00140 #include "TROOT.h"
00141 #include "TBaseClass.h"
00142 #include "TSystem.h"
00143 #include "TBuffer.h"
00144 #include "TMath.h"
00145 #include "TClass.h"
00146 #include "TBrowser.h"
00147 #include "TString.h"
00148 #include "Api.h"
00149 #include "TDataSetIter.h"
00150 #include "TTable.h"
00151 #include "TTableDescriptor.h"
00152 #include "TColumnView.h"
00153 
00154 #include "TGaxis.h"
00155 #include "TH1.h"
00156 #include "TH2.h"
00157 #include "TProfile.h"
00158 #include "TVirtualPad.h"
00159 #include "TEventList.h"
00160 #include "TPolyMarker.h"
00161 #include "TView.h"
00162 #include "TGaxis.h"
00163 #include "TPolyMarker3D.h"
00164 
00165 #include "THLimitsFinder.h"
00166 
00167 #include "TTableMap.h"
00168 
00169 static TH1 *gCurrentTableHist = 0;
00170 
00171 static const char *gDtorName = "dtor";
00172 static   Int_t         gNbins[4] = {100,100,100,100};     //Number of bins per dimension
00173 static   Float_t       gVmin[4]  = {0,0,0,0};             //Minima of varexp columns
00174 static   Float_t       gVmax[4]  = {20,20,20,20};         //Maxima of varexp columns
00175 
00176 const char *TTable::fgTypeName[] = {
00177    "NAN", "float", "int", "long", "short", "double",
00178    "unsigned int", "unsigned long","unsigned short",
00179    "unsigned char", "char", "Ptr_t"
00180 };
00181 
00182 //______________________________________________________________________________
00183 static void ArrayLayout(UInt_t *layout,const UInt_t *size, Int_t dim)
00184 {
00185   //
00186   // ArrayLayout - calculates the array layout recursively
00187   //
00188   // Input:
00189   // -----
00190   // dim   - dimension of the targeted array
00191   // size  - the max index for each dimension
00192   //
00193   // Output:
00194   // ------
00195   // layout - the "start index" for each dimension of an array
00196   //
00197 
00198    if (dim && layout && size) {
00199       if (++layout[dim-1] >= size[dim-1]) {
00200          layout[dim-1] = 0;
00201          dim--;
00202          ArrayLayout(layout,size, dim);
00203       }
00204    }
00205 }
00206 
00207 ClassImp(TTable)
00208 
00209 //______________________________________________________________________________
00210 TTableDescriptor *TTable::GetTableDescriptors() const {
00211  // protected: create a new TTableDescriptor descriptor for this table
00212    assert(0);
00213    return new TTableDescriptor(this);
00214 }
00215 
00216 //______________________________________________________________________________
00217 void TTable::AsString(void *buf, EColumnType type, Int_t width,ostream &out) const
00218 {
00219   //
00220   // AsString represents the value provided via "void *b" with type defined
00221   //          by "name"
00222   //
00223   //   void *buf  - the pointer to the value to be printed out.
00224   //        type  - the basic data type for the value above
00225   //       width  - the number of psotion to be used to print the value out
00226   //
00227    int prevPrec = out.precision();
00228    const std::ios_base::fmtflags prevFmt = out.flags();
00229    
00230    switch (type) {
00231       case kFloat:
00232          out << dec  << setw(width) << setprecision(width-3) << *(float *)buf;
00233          break;
00234       case kInt:
00235          out << dec  <<  setw(width) << *(int *)buf;
00236          break;
00237       case kLong:
00238          out << dec  << setw(width) << *(long *)buf;
00239          break;
00240       case kShort:
00241          out << dec  << setw(width) << *(short *)buf;
00242          break;
00243       case kDouble:
00244          out << dec  << setw(width) << setprecision(width-3) << *(double *)buf;
00245          break;
00246       case kUInt:
00247          out << dec  << setw(width) << *(unsigned int *)buf;
00248          break;
00249       case kULong:
00250          out << dec  << setw(width) << *(unsigned long *)buf;
00251          break;
00252       case kUShort:
00253          out  << setw(width) << "0x" << hex << *(unsigned short *)buf;
00254          break;
00255       case kUChar:
00256          out  << setw(width) << "0x" << hex << int(*(unsigned char *)buf);
00257          break;
00258       case kChar:
00259          out << setw(width) << *(char *)buf;
00260          break;
00261       case kBool:
00262          out << setw(width) << *(bool *)buf;
00263          break;
00264       case kPtr:
00265          out << "->" << setw(width) << *(void **)buf;
00266          break;
00267       default:
00268          out << "\"NaN\"";
00269          break;
00270    };
00271    out.precision(prevPrec);
00272    out.setf(prevFmt);
00273 }
00274 
00275 //________________________________________________________________________
00276 const char *TTable::GetTypeName(TTable::EColumnType type)
00277 {  
00278    //return table type name
00279    return  fgTypeName[type]; 
00280 }
00281 
00282 //________________________________________________________________________
00283 TTable::EColumnType TTable::GetTypeId(const char *typeName)
00284 {
00285    // return the Id of the C basic type by given name
00286    // return kNAN if the name provided fits no knwn basic name.
00287    //
00288    Int_t allTypes = sizeof(fgTypeName)/sizeof(const char *);
00289    for (int i = 0; i < allTypes; i++)
00290    if (!strcmp(fgTypeName[i],typeName)) return EColumnType(i);
00291    return kNAN;
00292 }
00293 
00294 //______________________________________________________________________________
00295 const void *TTable::At(Int_t i) const
00296 {
00297    // Returns a pointer to the i-th row of the table
00298    if (!BoundsOk("TTable::At", i)) {
00299        Warning("TTable::At","%s.%s",GetName(),GetType());
00300       i = 0;
00301    }
00302    return (const void *)(fTable+i*fSize);
00303 }
00304 
00305 //______________________________________________________________________________
00306 Int_t TTable::CopyRows(const TTable *srcTable, Long_t srcRow, Long_t dstRow, Long_t nRows, Bool_t expand)
00307 {
00308  // CopyRows copies nRows from starting from the srcRow of srcTable
00309  // to the dstRow in this table upto nRows or by the end of this table.
00310  //
00311  // This table if automaticaly increased if expand = kTRUE.
00312  // The old values of this table rows are to be destroyed and
00313  // replaced with the new ones.
00314  //
00315  // PARAMETERS:
00316  //   srcTable - a pointer to the table "donor"
00317  //   srcRow   - the index of the first row of the table donor to copy from
00318  //   dstRow   - the index of the first row of this table to copy to
00319  //   nRows    - the total number of rows to be copied. This table will be expanded
00320  //              as needed if expand = kTRUE (it is kFALSE "by default")
00321  //          = 0 to copy ALL remain rows from the srcTable.
00322  //   expand   - flag whether this table should reallocated if needed.
00323  //
00324  // RETURN:
00325  //          the number of the rows been copied
00326 
00327    assert(!TestBit(kIsNotOwn));
00328    if (!(srcTable && srcTable->GetNRows()) || srcRow > srcTable->GetNRows()-1   )   return 0;
00329    if (strcmp(GetType(),srcTable->GetType())) {
00330       // check this table current capacity
00331       if (!nRows) nRows = srcTable->GetNRows();
00332       Long_t tSize = GetTableSize();
00333       Long_t extraRows = (tSize - dstRow) - nRows;
00334       if (extraRows < 0) {
00335          if (expand) {
00336             ReAllocate(tSize - extraRows);
00337             extraRows = 0;
00338          }
00339          nRows += extraRows;
00340       }
00341       if (dstRow+nRows > GetNRows()) SetNRows(dstRow+nRows);
00342       ::memmove((*this)[dstRow],(*srcTable)[srcRow],(size_t)GetRowSize()*nRows);
00343       return nRows;
00344    } else
00345       Error("CopyRows",
00346            "This table is <%s> but the src table has a wrong type <%s>",GetType()
00347            ,srcTable->GetType());
00348    return 0;
00349 }
00350 //______________________________________________________________________________
00351 void TTable::DeleteRows(Long_t indx, UInt_t nRows)
00352 {
00353   // Delete one or several rows from the table
00354   //
00355   //  Int_t indx  - index of the first row to be deleted
00356   //  Int_t nRows - the total number of rows to be deleted
00357   //              = 1 "by default
00358    if (CopyRows(this, indx+nRows, indx, GetNRows()-indx-nRows))
00359       SetUsedRows(GetNRows() - nRows);
00360 }
00361 //______________________________________________________________________________
00362 TH1  *TTable::Draw(TCut varexp, TCut selection, Option_t *option, Int_t nentries, Int_t firstentry)
00363 {
00364 //*-*-*-*-*-*-*-*-*-*-*Draw expression varexp for specified entries-*-*-*-*-*
00365 //*-*                  ===========================================
00366 //
00367 //   This function accepts TCut objects as arguments.
00368 //   Useful to use the string operator +
00369 //         example:
00370 //            table.Draw("x",cut1+cut2+cut3);
00371 //
00372 //   TCutG object with "CUTG" name can be created via the graphics editor.
00373 //
00374 
00375    return TTable::Draw(varexp.GetTitle(), selection.GetTitle(), option, nentries, firstentry);
00376 }
00377 
00378 //______________________________________________________________________________
00379 TH1 *TTable::Draw(const char *varexp00, const char *selection, Option_t *option,Int_t nentries, Int_t firstentry)
00380 {
00381 //*-*-*-*-*-*-*-*-*-*-*Draw expression varexp for specified entries-*-*-*-*-*
00382 //*-*                  ===========================================
00383 //
00384 //  varexp is an expression of the general form e1:e2:e3
00385 //    where e1,etc is a C++ expression referencing a combination of the TTable columns
00386 //          One can use two extra meta variable "i$" and "n$" along with the table
00387 //          column names.
00388 //          i$ is to involve the current row number
00389 //          n$ refers the total num,ber of rows of this table provided by TTable::GetNRows()
00390 //
00391 //  Example:
00392 //     varexp = x     simplest case: draw a 1-Dim distribution of column named x
00393 //            = sqrt(x)            : draw distribution of sqrt(x)
00394 //            = x*y/z
00395 //            = y:sqrt(x) 2-Dim dsitribution of y versus sqrt(x)
00396 //            = i$:sqrt(x) 2-Dim dsitribution of i versus sqrt(x[i])
00397 //            = phep[0]:sqrt(phep[3]) 2-Dim dsitribution of phep[0] versus sqrt(phep[3])
00398 //
00399 //  Note that the variables e1, e2 or e3 may contain a boolean expression as well.
00400 //  example, if e1= x*(y<0), the value histogrammed will be x if y<0
00401 //  and will be 0 otherwise.
00402 //
00403 //  selection is a C++ expression with a combination of the columns.
00404 //  The value corresponding to the selection expression is used as a weight
00405 //  to fill the histogram.
00406 //  If the expression includes only boolean operations, the result
00407 //  is 0 or 1. If the result is 0, the histogram is not filled.
00408 //  In general, the expression may be of the form:
00409 //
00410 //      value*(boolean expression)
00411 //
00412 //  if boolean expression is true, the histogram is filled with
00413 //  a weight = value.
00414 //  Examples:
00415 //      selection1 = "x<y && sqrt(z)>3.2 && 6 < i$ && i$ < n$"
00416 //      selection2 = "(x+y)*(sqrt(z)>3.2"
00417 //      selection3 = "signal*(log(signal)>1.2)"
00418 //  selection1 returns a weigth = 0 or 1
00419 //  selection2 returns a weight = x+y if sqrt(z)>3.2
00420 //             returns a weight = 0 otherwise.
00421 //  selection3 returns a weight = signal if log(signal)>1.2
00422 //
00423 //  option is the drawing option
00424 //      see TH1::Draw for the list of all drawing options.
00425 //      If option contains the string "goff", no graphics is generated.
00426 //
00427 //  nentries is the number of entries to process (default is all)
00428 //  first is the first entry to process (default is 0)
00429 //
00430 //     Saving the result of Draw to an histogram
00431 //     =========================================
00432 //  By default the temporary histogram created is called htemp.
00433 //  If varexp0 contains >>hnew (following the variable(s) name(s),
00434 //  the new histogram created is called hnew and it is kept in the current
00435 //  directory.
00436 //  Example:
00437 //    tree.Draw("sqrt(x)>>hsqrt","y>0")
00438 //    will draw sqrt(x) and save the histogram as "hsqrt" in the current
00439 //    directory.
00440 //
00441 //  By default, the specified histogram is reset.
00442 //  To continue to append data to an existing histogram, use "+" in front
00443 //  of the histogram name;
00444 //    table.Draw("sqrt(x)>>+hsqrt","y>0")
00445 //      will not reset hsqrt, but will continue filling.
00446 //
00447 //     Making a Profile histogram
00448 //     ==========================
00449 //  In case of a 2-Dim expression, one can generate a TProfile histogram
00450 //  instead of a TH2F histogram by specyfying option=prof or option=profs.
00451 //  The option=prof is automatically selected in case of y:x>>pf
00452 //  where pf is an existing TProfile histogram.
00453 //
00454 //     Saving the result of Draw to a TEventList
00455 //     =========================================
00456 //  TTable::Draw can be used to fill a TEventList object (list of entry numbers)
00457 //  instead of histogramming one variable.
00458 //  If varexp0 has the form >>elist , a TEventList object named "elist"
00459 //  is created in the current directory. elist will contain the list
00460 //  of entry numbers satisfying the current selection.
00461 //  Example:
00462 //    tree.Draw(">>yplus","y>0")
00463 //    will create a TEventList object named "yplus" in the current directory.
00464 //    In an interactive session, one can type (after TTable::Draw)
00465 //       yplus.Print("all")
00466 //    to print the list of entry numbers in the list.
00467 //
00468 //  By default, the specified entry list is reset.
00469 //  To continue to append data to an existing list, use "+" in front
00470 //  of the list name;
00471 //    table.Draw(">>+yplus","y>0")
00472 //      will not reset yplus, but will enter the selected entries at the end
00473 //      of the existing list.
00474 //
00475 
00476    if (GetNRows() == 0 || varexp00 == 0 || varexp00[0]==0) return 0;
00477    TString  opt;
00478 //   char *hdefault = (char *)"htemp";
00479    const char *hdefault = "htemp";
00480    Int_t i,j,action;
00481    Int_t hkeep = 0;
00482    opt = option;
00483    opt.ToLower();
00484    char *varexp0 = StrDup(varexp00);
00485    char *hname = strstr(varexp0,">>");
00486    TH1 *oldh1 = 0;
00487    TEventList *elist = 0;
00488    Bool_t profile = kFALSE;
00489 
00490    gCurrentTableHist = 0;
00491    if (hname) {
00492       *hname  = 0;
00493       hname += 2;
00494       hkeep  = 1;
00495       i = strcspn(varexp0,">>");
00496       Bool_t hnameplus = kFALSE;
00497       while (*hname == ' ') hname++;
00498       if (*hname == '+') {
00499          hnameplus = kTRUE;
00500          hname++;
00501          while (*hname == ' ') hname++;
00502          j = strlen(hname)-1;
00503          while (j) {
00504             if (hname[j] != ' ') break;
00505             hname[j] = 0;
00506             j--;
00507          }
00508       }
00509       if (i) {
00510          oldh1 = (TH1*)gDirectory->Get(hname);
00511          if (oldh1 && !hnameplus) oldh1->Reset();
00512       } else {
00513          elist = (TEventList*)gDirectory->Get(hname);
00514          if (!elist) {
00515             elist = new TEventList(hname,selection,1000,0);
00516          }
00517          if (elist && !hnameplus) elist->Reset();
00518       }
00519    }
00520    if (!hname || *hname==0) {
00521       hkeep  = 0;
00522       if (gDirectory) {
00523          oldh1 = (TH1*)gDirectory->Get(hdefault);
00524          if (oldh1 ) { oldh1->Delete(); oldh1 = 0;}
00525       }
00526    }
00527 
00528    // Look for colons
00529    const Char_t *expressions[] ={varexp0,0,0,0,selection};
00530    Int_t maxExpressions = sizeof(expressions)/sizeof(Char_t *);
00531    Char_t *nextColon    = varexp0;
00532    Int_t colIndex       = 1;
00533    while ((nextColon = strchr(nextColon,':')) && ( colIndex < maxExpressions - 1 ) ) {
00534       *nextColon = 0;
00535       nextColon++;
00536       expressions[colIndex] = nextColon;
00537       colIndex++;
00538    }
00539 
00540    expressions[colIndex] = selection;
00541 
00542 
00543 //--------------------------------------------------
00544    Printf(" Draw %s for <%s>\n", varexp00, selection);
00545    Char_t *exprFileName = MakeExpression(expressions,colIndex+1);
00546    if (!exprFileName) {
00547       delete [] varexp0;
00548       return 0;
00549    }
00550 
00551 //--------------------------------------------------
00552 //   if (!fVar1 && !elist) return 0;
00553 
00554 //*-*- In case oldh1 exists, check dimensionality
00555    Int_t dimension = colIndex;
00556 
00557    TString title = expressions[0];
00558    for (i=1;i<colIndex;i++) {
00559       title += ":";
00560       title += expressions[i];
00561    }
00562    Int_t nsel = strlen(selection);
00563    if (nsel > 1) {
00564       if (nsel < 80-title.Length()) {
00565          title += "{";
00566          title += selection;
00567          title += "}";
00568       } else
00569          title += "{...}";
00570    }
00571 
00572    const Char_t *htitle = title.Data();
00573 
00574    if (oldh1) {
00575       Int_t mustdelete = 0;
00576       if (oldh1->InheritsFrom(TProfile::Class())) profile = kTRUE;
00577       if (opt.Contains("prof")) {
00578          if (!profile) mustdelete = 1;
00579       } else {
00580          if (oldh1->GetDimension() != dimension) mustdelete = 1;
00581       }
00582       if (mustdelete) {
00583          Warning("Draw","Deleting old histogram with different dimensions");
00584          delete oldh1; oldh1 = 0;
00585       }
00586    }
00587 //*-*- Create a default canvas if none exists
00588    if (!gPad && !opt.Contains("goff") && dimension > 0) {
00589       gROOT->MakeDefCanvas();
00590    }
00591 //*-*- 1-D distribution
00592    if (dimension == 1) {
00593       action = 1;
00594       if (!oldh1) {
00595          gNbins[0] = 100;
00596          if (gPad && opt.Contains("same")) {
00597             TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
00598             if (oldhtemp) {
00599                gNbins[0] = oldhtemp->GetXaxis()->GetNbins();
00600                gVmin[0]  = oldhtemp->GetXaxis()->GetXmin();
00601                gVmax[0]  = oldhtemp->GetXaxis()->GetXmax();
00602             } else {
00603                gVmin[0]  = gPad->GetUxmin();
00604                gVmax[0]  = gPad->GetUxmax();
00605             }
00606          } else {
00607             action = -1;
00608          }
00609       }
00610       TH1F *h1;
00611       if (oldh1) {
00612          h1 = (TH1F*)oldh1;
00613          gNbins[0] = h1->GetXaxis()->GetNbins();  // for proofserv
00614       } else {
00615          h1 = new TH1F(hname,htitle,gNbins[0],gVmin[0],gVmax[0]);
00616          if (!hkeep) {
00617             h1->SetBit(kCanDelete);
00618             h1->SetDirectory(0);
00619          }
00620          if (opt.Length() && opt[0] == 'e') h1->Sumw2();
00621       }
00622 
00623       EntryLoop(exprFileName,action, h1, nentries, firstentry, option);
00624 
00625 //    if (!fDraw && !opt.Contains("goff")) h1->Draw(option);
00626       if (!opt.Contains("goff")) h1->Draw(option);
00627 
00628 //*-*- 2-D distribution
00629    } else if (dimension == 2) {
00630       action = 2;
00631       if (!opt.Contains("same") && gPad)  gPad->Clear();
00632       if (!oldh1 || !opt.Contains("same")) {
00633          gNbins[0] = 40;
00634          gNbins[1] = 40;
00635          if (opt.Contains("prof")) gNbins[1] = 100;
00636          if (opt.Contains("same")) {
00637             TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
00638             if (oldhtemp) {
00639                gNbins[1] = oldhtemp->GetXaxis()->GetNbins();
00640                gVmin[1]  = oldhtemp->GetXaxis()->GetXmin();
00641                gVmax[1]  = oldhtemp->GetXaxis()->GetXmax();
00642                gNbins[0] = oldhtemp->GetYaxis()->GetNbins();
00643                gVmin[0]  = oldhtemp->GetYaxis()->GetXmin();
00644                gVmax[0]  = oldhtemp->GetYaxis()->GetXmax();
00645             } else {
00646                gNbins[1] = 40;
00647                gVmin[1]  = gPad->GetUxmin();
00648                gVmax[1]  = gPad->GetUxmax();
00649                gNbins[0] = 40;
00650                gVmin[0]  = gPad->GetUymin();
00651                gVmax[0]  = gPad->GetUymax();
00652             }
00653          } else {
00654             action = -2;
00655          }
00656       }
00657       if (profile || opt.Contains("prof")) {
00658          TProfile *hp;
00659          if (oldh1) {
00660             action = 4;
00661             hp = (TProfile*)oldh1;
00662          } else {
00663             if (action < 0) action = -4;
00664             if (opt.Contains("profs"))
00665                hp = new TProfile(hname,htitle,gNbins[1],gVmin[1], gVmax[1],"s");
00666             else
00667                hp = new TProfile(hname,htitle,gNbins[1],gVmin[1], gVmax[1],"");
00668             if (!hkeep) {
00669                hp->SetBit(kCanDelete);
00670                hp->SetDirectory(0);
00671             }
00672          }
00673 
00674          EntryLoop(exprFileName,action,hp,nentries, firstentry, option);
00675 
00676          if (!opt.Contains("goff")) hp->Draw(option);
00677       } else {
00678          TH2F *h2;
00679          if (oldh1) {
00680             h2 = (TH2F*)oldh1;
00681          } else {
00682             h2 = new TH2F(hname,htitle,gNbins[1],gVmin[1], gVmax[1], gNbins[0], gVmin[0], gVmax[0]);
00683             if (!hkeep) {
00684                const Int_t kNoStats = BIT(9);
00685                h2->SetBit(kCanDelete);
00686                h2->SetBit(kNoStats);
00687                h2->SetDirectory(0);
00688             }
00689          }
00690          Int_t noscat = strlen(option);
00691          if (opt.Contains("same")) noscat -= 4;
00692          if (noscat) {
00693             EntryLoop(exprFileName,action,h2,nentries, firstentry, option);
00694  //           if (!fDraw && !opt.Contains("goff")) h2->Draw(option);
00695             if (!opt.Contains("goff")) h2->Draw(option);
00696          } else {
00697             action = 12;
00698             if (!oldh1 && !opt.Contains("same")) action = -12;
00699             EntryLoop(exprFileName,action,h2,nentries, firstentry, option);
00700 //            if (oldh1 && !fDraw && !opt.Contains("goff")) h2->Draw(option);
00701             if (oldh1 && !opt.Contains("goff")) h2->Draw(option);
00702          }
00703       }
00704 
00705 //*-*- 3-D distribution
00706    } else if (dimension == 3) {
00707       action = 13;
00708       if (!opt.Contains("same")) action = -13;
00709       EntryLoop(exprFileName,action,0,nentries, firstentry, option);
00710 
00711 //*-* an Event List
00712    //} else if (elist) {
00713    //   action = 5;
00714 //      Int_t oldEstimate = fEstimate;
00715 //      SetEstimate(1);
00716    //   EntryLoop(exprFileName,action,elist,nentries, firstentry, option);
00717 //      SetEstimate(oldEstimate);
00718    }
00719    delete [] exprFileName;
00720    delete [] varexp0;
00721    return gCurrentTableHist;
00722 }
00723 
00724 //______________________________________________________________________________
00725 static void FindGoodLimits(Int_t nbins, Int_t &newbins, Float_t &xmin, Float_t &xmax)
00726 {
00727 //*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
00728 //*-*              ==========================
00729 //*-*  This mathod is a straight copy of void TTree::FindGoodLimits method
00730 //*-*
00731 
00732    Double_t binlow,binhigh,binwidth=0;
00733    Int_t n;
00734    Double_t dx = 0.1*(xmax-xmin);
00735    Double_t umin = xmin - dx;
00736    Double_t umax = xmax + dx;
00737    if (umin < 0 && xmin >= 0) umin = 0;
00738    if (umax > 0 && xmax <= 0) umax = 0;
00739 
00740 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,03,5)
00741    THLimitsFinder::Optimize(umin,umax,nbins,binlow,binhigh,n,binwidth,"");
00742 #else
00743    TGaxis::Optimize(umin,umax,nbins,binlow,binhigh,n,binwidth,"");
00744 #endif
00745 
00746    if (binwidth <= 0 || binwidth > 1.e+39) {
00747       xmin = -1;
00748       xmax = 1;
00749    } else {
00750       xmin    = binlow;
00751       xmax    = binhigh;
00752    }
00753 
00754    newbins = nbins;
00755 }
00756 
00757 //______________________________________________________________________________
00758 Bool_t TTable::EntryLoop(const Char_t *exprFileName,Int_t &action, TObject *obj
00759                           ,Int_t nentries, Int_t firstentry, Option_t *option)
00760 {
00761  //
00762  // EntryLoop creates a CINT bytecode to evaluate the given expressions for
00763  // all table rows in loop and fill the appropriated histograms.
00764  //
00765  // Solution for Byte code
00766  // From: Masaharu Goto <MXJ02154@nifty.ne.jp>
00767  // To: <fine@bnl.gov>
00768  // Cc: <rootdev@hpsalo.cern.ch>
00769  // Sent: 13-th august 1999 year  23:01
00770  //
00771  //  action =  1  Fill 1-D histogram obj
00772  //         =  2  Fill 2-D histogram obj
00773  //         =  3  Fill 3-D histogram obj
00774  //         =  4  Fill Profile histogram obj
00775  //         =  5  Fill a TEventlist
00776  //         = 11  Estimate Limits
00777  //         = 12  Fill 2-D PolyMarker obj
00778  //         = 13  Fill 3-D PolyMarker obj
00779  //  action < 0   Evaluate Limits for case abs(action)
00780  //
00781  //  Load file
00782    Double_t rmin[3],rmax[3];
00783    switch(G__loadfile((Char_t *)exprFileName)) {
00784       case G__LOADFILE_SUCCESS:
00785       case G__LOADFILE_DUPLICATE:
00786          break;
00787       default:
00788          Error("EntryLoop","Error: loading file %s",exprFileName);
00789          G__unloadfile((Char_t *)exprFileName);
00790          return kFALSE; // can not load file
00791    }
00792 
00793    // Float_t  Selection(Float_t *results[], void *address[], int& i$, int n$)
00794    //   where  i$ - meta variable to set current row index
00795    //          n$ - meta variable to set the total table size
00796    const Char_t *funcName = "SelectionQWERTY";
00797 #define BYTECODE
00798 #ifdef BYTECODE
00799    const Char_t *argtypes = "Float_t *,float **, int&, int& ";
00800    Long_t offset;
00801    G__ClassInfo globals;
00802    G__MethodInfo func = globals.GetMethod(funcName,argtypes,&offset);
00803 
00804    // Compile bytecode
00805    struct G__bytecodefunc *pbc = func.GetBytecode();
00806    if(!pbc) {
00807       Error("EntryLoop","Bytecode compilation %s",funcName);
00808       G__unloadfile((Char_t *)exprFileName);
00809       return kFALSE; // can not get bytecode
00810    }
00811 #endif
00812    // Prepare callfunc object
00813    int i;
00814    int nRows =  GetNRows();
00815    TTableDescriptor    *tabsDsc   = GetRowDescriptors();
00816    tableDescriptor_st  *descTable = tabsDsc->GetTable();
00817    Float_t  results[]    = {1,1,1,1,1};
00818    Char_t **addressArray = (Char_t **)new ULong_t[tabsDsc->GetNRows()];
00819    Char_t *thisTable     = (Char_t *)GetArray();
00820 #ifdef BYTECODE
00821    G__CallFunc callfunc;
00822    callfunc.SetBytecode(pbc);
00823 
00824    callfunc.SetArg((long)(&results[0]));   // give 'Float_t *results[5]'     as 1st argument
00825    callfunc.SetArg((long)(addressArray));  // give 'void    *addressArray[]' as 2nd argument
00826    callfunc.SetArg((long)(&i));            // give 'int& i$'                 as 3nd argument
00827    callfunc.SetArg((long)(&nRows));        // give 'int& n$= nRows           as 4th argument
00828 #else
00829    char buf[200];
00830    sprintf(buf,"%s((Float_t*)(%ld),(void**)(%ld),*(int*)(%ld),*(int*)(%ld))"
00831              ,funcName
00832              ,(long int)results,(long int)addressArray,(long int)(&i),(long int)(&nRows));
00833 #endif
00834 
00835    // Call bytecode in loop
00836 
00837 #ifdef BYTECODE
00838 #  define CALLMETHOD callfunc.Exec(0);
00839 #else
00840 #  define CALLMETHOD G__calc(buf);
00841 #endif
00842 
00843 #define TAKEACTION_BEGIN                                                                    \
00844             descTable = tabsDsc->GetTable();                                                \
00845             for (i=0; i < tabsDsc->GetNRows(); i++,descTable++ )                            \
00846                addressArray[i] = addressEntry + descTable->fOffset;                         \
00847             for(i=firstentry;i<lastEntry;i++) {                                             \
00848             CALLMETHOD
00849 
00850 #define TAKEACTION_END  for (int j=0; j < tabsDsc->GetNRows(); j++ ) addressArray[j] += rSize;}
00851 
00852 
00853    if (firstentry < nRows ) {
00854       Long_t rSize         = GetRowSize();
00855       Char_t *addressEntry = thisTable + rSize*firstentry;
00856       Int_t lastEntry = TMath::Min(UInt_t(firstentry+nentries),UInt_t(nRows));
00857       if (action < 0) {
00858          gVmin[0] = gVmin[1] = gVmin[2] = 1e30;
00859          gVmax[0] = gVmax[1] = gVmax[2] = -gVmin[0];
00860       }
00861       Int_t nchans = 0;
00862       switch ( action ) {
00863          case -1: {
00864             TAKEACTION_BEGIN
00865             if (results[1]) {
00866                if (gVmin[0] > results[0]) gVmin[0] = results[0];
00867                if (gVmax[0] < results[0]) gVmax[0] = results[0];
00868             }
00869             TAKEACTION_END
00870 
00871             nchans = gNbins[0];
00872             if (gVmin[0] >= gVmax[0]) { gVmin[0] -= 1; gVmax[0] += 1;}
00873             FindGoodLimits(nchans,gNbins[0],gVmin[0],gVmax[0]);
00874             ((TH1 *)obj)->SetBins(gNbins[0],gVmin[0],gVmax[0]);
00875          }
00876          // Intentional fall though
00877          case  1:
00878             TAKEACTION_BEGIN
00879             if (results[1]) ((TH1 *)obj)->Fill(Axis_t(results[0]),Stat_t(results[1]));
00880             TAKEACTION_END
00881             gCurrentTableHist = ((TH1 *)obj);
00882             break;
00883          case  -2:
00884             TAKEACTION_BEGIN
00885             if (results[2]) {
00886                if (gVmin[0] > results[1]) gVmin[0] = results[1];
00887                if (gVmax[0] < results[1]) gVmax[0] = results[1];
00888                if (gVmin[1] > results[0]) gVmin[1] = results[0];
00889                if (gVmax[1] < results[0]) gVmax[1] = results[0];
00890             }
00891             TAKEACTION_END
00892             nchans = gNbins[0];
00893             if (gVmin[0] >= gVmax[0]) { gVmin[0] -= 1; gVmax[0] += 1;}
00894             FindGoodLimits(nchans,gNbins[0],gVmin[0],gVmax[0]);
00895             if (gVmin[1] >= gVmax[1]) { gVmin[1] -= 1; gVmax[1] += 1;}
00896             FindGoodLimits(nchans,gNbins[1],gVmin[1],gVmax[1]);
00897             ((TH1*)obj)->SetBins(gNbins[1],gVmin[1],gVmax[1],gNbins[0],gVmin[0],gVmax[0]);
00898             // Intentional fall though
00899          case   2:
00900             if (obj->IsA() == TH2F::Class()) {
00901                TAKEACTION_BEGIN
00902                if (results[2]) ((TH2F*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
00903                TAKEACTION_END
00904             } else if (obj->IsA() == TH2S::Class()) {
00905                TAKEACTION_BEGIN
00906                if (results[2]) ((TH2S*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
00907                TAKEACTION_END
00908             } else if (obj->IsA() == TH2C::Class()) {
00909                TAKEACTION_BEGIN
00910                if (results[2]) ((TH2C*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
00911                TAKEACTION_END
00912             } else if (obj->IsA() == TH2D::Class()) {
00913                TAKEACTION_BEGIN
00914                if (results[2]) ((TH2D*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
00915                TAKEACTION_END
00916             }
00917             gCurrentTableHist =  ((TH1 *)obj);
00918             break;
00919          case -4:
00920             TAKEACTION_BEGIN
00921             if (results[2]) {
00922                if (gVmin[0] > results[1]) gVmin[0] = results[1];
00923                if (gVmax[0] < results[1]) gVmax[0] = results[1];
00924                if (gVmin[1] > results[0]) gVmin[1] = results[0];
00925                if (gVmax[1] < results[0]) gVmax[1] = results[0];
00926             }
00927             TAKEACTION_END
00928             nchans = gNbins[1];
00929             if (gVmin[1] >= gVmax[1]) { gVmin[1] -= 1; gVmax[1] += 1;}
00930             FindGoodLimits(nchans,gNbins[1],gVmin[1],gVmax[1]);
00931             ((TProfile*)obj)->SetBins(gNbins[1],gVmin[1],gVmax[1]);
00932             // Intentional fall though
00933          case  4:
00934             TAKEACTION_BEGIN
00935             if (results[2]) ((TProfile*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
00936             TAKEACTION_END
00937             break;
00938          case -12:
00939             TAKEACTION_BEGIN
00940             if (results[2]) {
00941                if (gVmin[0] > results[1]) gVmin[0] = results[1];
00942                if (gVmax[0] < results[1]) gVmax[0] = results[1];
00943                if (gVmin[1] > results[0]) gVmin[1] = results[0];
00944                if (gVmax[1] < results[0]) gVmax[1] = results[0];
00945             }
00946             TAKEACTION_END
00947             nchans = gNbins[0];
00948             if (gVmin[0] >= gVmax[0]) { gVmin[0] -= 1; gVmax[0] += 1;}
00949             FindGoodLimits(nchans,gNbins[0],gVmin[0],gVmax[0]);
00950             if (gVmin[1] >= gVmax[1]) { gVmin[1] -= 1; gVmax[1] += 1;}
00951             FindGoodLimits(nchans,gNbins[1],gVmin[1],gVmax[1]);
00952             ((TH2F*)obj)->SetBins(gNbins[1],gVmin[1],gVmax[1],gNbins[0],gVmin[0],gVmax[0]);
00953             // Intentional fall though
00954          case  12: {
00955             if (!strstr(option,"same") && !strstr(option,"goff")) {
00956                ((TH2F*)obj)->DrawCopy(option);
00957                gPad->Update();
00958             }
00959 //            pm->SetMarkerStyle(GetMarkerStyle());
00960 //            pm->SetMarkerColor(GetMarkerColor());
00961 //            pm->SetMarkerSize(GetMarkerSize());
00962             Float_t *x = new Float_t[lastEntry-firstentry]; // pm->GetX();
00963             Float_t *y = new Float_t[lastEntry-firstentry]; // pm->GetY();
00964             Float_t u, v;
00965             Float_t umin = gPad->GetUxmin();
00966             Float_t umax = gPad->GetUxmax();
00967             Float_t vmin = gPad->GetUymin();
00968             Float_t vmax = gPad->GetUymax();
00969             Int_t pointIndex = 0;
00970             TAKEACTION_BEGIN
00971             if (results[2]) {
00972                u = gPad->XtoPad(results[0]);
00973                v = gPad->YtoPad(results[1]);
00974                if (u < umin) u = umin;
00975                if (u > umax) u = umax;
00976                if (v < vmin) v = vmin;
00977                if (v > vmax) v = vmax;
00978                x[pointIndex] = u;
00979                y[pointIndex] = v;
00980                pointIndex++;
00981             }
00982             TAKEACTION_END
00983             if (pointIndex && !strstr(option,"goff")) {
00984                TPolyMarker *pm = new TPolyMarker(pointIndex,x,y);
00985                pm->Draw();
00986                pm->SetBit(kCanDelete);
00987             }
00988             if (!((TH2F*)obj)->TestBit(kCanDelete))
00989                if (pointIndex)
00990                   for(i=0;i<pointIndex;i++) ((TH2F*)obj)->Fill(x[i], y[i]);
00991             delete [] x; delete [] y;
00992             gCurrentTableHist = ((TH1*)obj);
00993             }
00994             break;
00995          case -13:
00996             TAKEACTION_BEGIN
00997             if (results[3]) {
00998                if (gVmin[0] > results[2]) gVmin[0] = results[2];
00999                if (gVmax[0] < results[2]) gVmax[0] = results[2];
01000                if (gVmin[1] > results[1]) gVmin[1] = results[1];
01001                if (gVmax[1] < results[1]) gVmax[1] = results[1];
01002                if (gVmin[2] > results[0]) gVmin[2] = results[0];
01003                if (gVmax[2] < results[0]) gVmax[2] = results[0];
01004             }
01005             TAKEACTION_END
01006             rmin[0] = gVmin[2]; rmin[1] = gVmin[1]; rmin[2] = gVmin[0];
01007             rmax[0] = gVmax[2]; rmax[1] = gVmax[1]; rmax[2] = gVmax[0];
01008             gPad->Clear();
01009             gPad->Range(-1,-1,1,1);
01010             TView::CreateView(1,rmin,rmax);
01011             // Intentional fall though
01012          case 13: {
01013             TPolyMarker3D *pm3d = new TPolyMarker3D(lastEntry-firstentry);
01014             pm3d->SetBit(kCanDelete);
01015 //            pm3d->SetMarkerStyle(GetMarkerStyle());
01016 //            pm3d->SetMarkerColor(GetMarkerColor());
01017 //            pm3d->SetMarkerSize(GetMarkerSize());
01018             TAKEACTION_BEGIN
01019             if (results[3]) pm3d->SetNextPoint(results[0],results[1],results[2]);
01020             TAKEACTION_END
01021             pm3d->Draw();
01022          }
01023          break;
01024          default:
01025             Error("EntryLoop","unknown action \"%d\" for table <%s>", action, GetName());
01026             break;
01027       };
01028    }
01029    G__unloadfile((Char_t *)exprFileName);
01030    delete [] addressArray;
01031    return kTRUE;
01032 }
01033 
01034 //______________________________________________________________________________
01035 TTable::TTable(const char *name, Int_t size) : TDataSet(name),
01036          fSize(size),fN(0), fTable(0),fMaxIndex(0)
01037 {
01038    // Default TTable ctor.
01039    if (size == 0) Warning("TTable(0)","Wrong table format");
01040 }
01041 
01042 //______________________________________________________________________________
01043 TTable::TTable(const char *name, Int_t n,Int_t size) : TDataSet(name),
01044         fSize(size),fN(0),fTable(0),fMaxIndex(0)
01045 {
01046    // Create TTable object and set array size to n longs.
01047    if (n > 0) Set(n);
01048 }
01049 
01050 //______________________________________________________________________________
01051 TTable::TTable(const char *name, Int_t n, Char_t *table,Int_t size) : TDataSet(name),
01052          fSize(size),fN(0),fTable(0),fMaxIndex(0)
01053 {
01054    // Create TTable object and initialize it with values of array.
01055 
01056    Set(n, table);
01057 }
01058 
01059 //______________________________________________________________________________
01060 TTable::TTable(const char *name, const char *type, Int_t n, Char_t *array, Int_t size)
01061          : TDataSet(name),fSize(size),fTable(0),fMaxIndex(0)
01062 {
01063    // Create TTable object and initialize it with values of array.
01064 
01065    fTable = array;
01066    SetType(type);
01067    SetfN(n);
01068 }
01069 
01070 //______________________________________________________________________________
01071 TTable::TTable(const TTable &table):TDataSet(table)
01072 {
01073    // Copy constructor.
01074    fTable    = 0;
01075    SetUsedRows(table.GetNRows());
01076    fSize     = table.GetRowSize();
01077    SetfN(table.fN);
01078    Set(table.fN, table.fTable);
01079 }
01080 
01081 //______________________________________________________________________________
01082 TTable &TTable::operator=(const TTable &rhs)
01083 {
01084    // TTable assignment operator.
01085    // This operator REALLOCATEs this table to fit the number of
01086    // the USED rows of the source table if any
01087 
01088    if (strcmp(GetType(),rhs.GetType()) == 0) {
01089       if (this != &rhs && rhs.GetNRows() >0 ){
01090          Set(rhs.GetNRows(), rhs.fTable);
01091          SetUsedRows(rhs.GetNRows());
01092       }
01093    } else
01094       Error("operator=","Can not copy <%s> table into <%s> table", rhs.GetType(),GetType());
01095    return *this;
01096 }
01097 
01098 //______________________________________________________________________________
01099 TTable::~TTable()
01100 {
01101    // Delete TTable object.
01102    Delete();
01103 }
01104 
01105 //______________________________________________________________________________
01106 void TTable::Adopt(Int_t n, void *arr)
01107 {
01108    // Adopt array arr into TTable, i.e. don't copy arr but use it directly
01109    // in TTable. User may not delete arr, TTable dtor will do it.
01110 
01111    Clear();
01112 
01113    SetfN(n); SetUsedRows(n);
01114    fTable = (char *)arr;
01115 }
01116 
01117 //______________________________________________________________________________
01118 Int_t TTable::AddAt(const void *row)
01119 {
01120   // Add        the "row" at the GetNRows() position, and
01121   // reallocate the table if neccesary,               and
01122   // return     the row index the "row" has occupied.
01123   //
01124   // row == 0 see method TTable::AddAt(const void *row, Int_t i)
01125 
01126    Int_t gap = GetTableSize() - GetNRows();
01127    // do we need to add an extra space?
01128    if (gap < 1) ReAllocate(GetTableSize() + TMath::Max(1,Int_t(0.3*GetTableSize())));
01129    Int_t indx = GetNRows();
01130    AddAt(row,indx);
01131    return indx;
01132 }
01133 //______________________________________________________________________________
01134 void TTable::AddAt(const void *row, Int_t i)
01135 {
01136    // Add    one element ("row") of structure at position "i".
01137    // Check  for out of bounds.
01138    //
01139    //        If the row == 0 the "i" cell is still occupied and
01140    // filled with the pattern "ff"
01141 
01142    if (!BoundsOk("TTable::AddAt", i))
01143       i = 0;
01144    if (row) memcpy(fTable+i*fSize,row,fSize);
01145    else memset(fTable+i*fSize,127,fSize);
01146    SetUsedRows(TMath::Max((Int_t)i+1,Int_t(fMaxIndex)));
01147 }
01148 
01149 //______________________________________________________________________________
01150 void TTable::CopyStruct(Char_t *dest, const Char_t *src)
01151 {
01152  // Copy the C-structure src into the new location
01153  // the length of the strucutre is defined by this class descriptor
01154    ::memcpy(dest,src,fSize*fN);
01155 }
01156 //______________________________________________________________________________
01157 void TTable::CopySet(TTable &array)
01158 {
01159    //to be documented
01160    array.Set(fN);
01161    CopyStruct(array.fTable,fTable);
01162 }
01163 //______________________________________________________________________________
01164 const Char_t *TTable::GetColumnComment(Int_t columnIndex) const {
01165    // Get a comment from the table descriptor
01166    TDataSetIter nextComment(GetRowDescriptors()->MakeCommentField(kFALSE));
01167    TDataSet *nxc = 0;
01168    for (int i=0; i<= columnIndex; i++) nxc = nextComment();
01169    return nxc ? nxc->GetTitle() : 0;
01170 }
01171 //______________________________________________________________________________
01172 Long_t TTable::AppendRows(const void *row, UInt_t nRows)
01173 {
01174    // Append nRows row of the array "row" to the table
01175    // return
01176    //    - the new table size (# of table rows)
01177    //    - 0 if the object doesn't own the internal array and can not expand it
01178    Long_t size = 0;
01179    if (!TestBit(kIsNotOwn) && row && nRows ) {
01180       Int_t indx = GetNRows();
01181       ReAllocate(nRows);
01182       // Copy (insert) the extra staff in
01183       ::memmove(fTable+indx*fSize,row,fSize*nRows);
01184       size = GetSize();
01185    }
01186    return TestBit(kIsNotOwn) ? 0 : GetSize();
01187 }
01188 //______________________________________________________________________________
01189 Long_t TTable::InsertRows(const void *row, Long_t indx, UInt_t nRows)
01190 {
01191   // void InsertRows(cons void *row, Long_t indx, UInt_t nRows)
01192   //
01193   // Insert one or several rows into the table at "indx" position
01194   // The rest table stuff is shifted down
01195   //
01196   //  cons void    - a pointer to the array of rows to be inserted
01197   //  Long_t indx =  The position these rows will be inserted to
01198   //  Int_t nRows  - the total number of rows to be inserted
01199   //                 = 1 "by default
01200   //  return:
01201   //  The number of the rows has been shifted to accomodate
01202   //  the new rows.
01203   //
01204    Long_t nShifted = 0;
01205    if (nRows > 0) {
01206       // Shift the table down
01207       nShifted = CopyRows(this, indx, indx+nRows, GetNRows()+nRows);
01208       // Copy (insert) the extra staff in
01209       ::memmove(fTable+indx*fSize,row,fSize*nRows);
01210    }
01211    return nShifted;
01212 
01213 }
01214 //______________________________________________________________________________
01215 void *TTable::ReAllocate()
01216 {
01217   // Reallocate this table leaving only (used rows)+1 allocated
01218   // GetTableSize() = GetNRows() + 1
01219   // returns a pointer to the first row of the reallocated table
01220   // Note:
01221   // The table is reallocated if it is an owner of the internal array
01222 
01223    ReAlloc(GetNRows()+1);
01224    return (void *)fTable;
01225 }
01226 //______________________________________________________________________________
01227 void *TTable::ReAllocate(Int_t newsize)
01228 {
01229   // Reallocate this table leaving only <newsize> allocated
01230   // GetTableSize() = newsize;
01231   // returns a pointer to the first row of the reallocated table
01232   // Note:
01233   // The table is reallocated if it is an owner of the internal array
01234 
01235    if (newsize > fN) ReAlloc(newsize);
01236    return (void *)fTable;
01237 }
01238 
01239 //______________________________________________________________________________
01240 void TTable::ReAlloc(Int_t newsize)
01241 {
01242   // The table is reallocated if it is an owner of the internal array
01243    if (!TestBit(kIsNotOwn) && newsize > 0) {
01244       void *arr = 0;
01245       Int_t sleepCounter = 0;
01246       while (!(arr =  realloc(fTable,fSize*newsize))) {
01247          sleepCounter++;
01248          Warning("ReAlloc",
01249               "Not enough memory to Reallocate %d bytes for table <%s::%s>. Please cancel some jobs",
01250               newsize, GetType(),GetName());
01251          gSystem->Sleep(1000*600);
01252          if (sleepCounter > 30) {
01253             Error("ReAlloc","I can not wait anymore. Good bye");
01254             assert(0);
01255          }
01256       }
01257       SetfN(newsize);
01258       fTable = (char *)arr;
01259    }
01260 }
01261 
01262 //______________________________________________________________________________
01263 Char_t *TTable::Create()
01264 {
01265    // Allocate a space for the new table, if any
01266    // Sleep for a while if space is not available and try again
01267    if (!fTable) {
01268       void *ptr = 0;
01269       Int_t sleepCounter = 0;
01270       while (!(ptr = malloc(fSize*fN))) {
01271          sleepCounter++;
01272          Warning("Create",
01273             "Not enough memory to allocate %d rows for table <%s::%s>. Please cancel some jobs",
01274             fN, GetType(),GetName());
01275          gSystem->Sleep(1000*600);
01276          if (sleepCounter > 30){
01277             Error("Create","I can not wait anymore. Good bye");
01278             assert(0);
01279          }
01280       }
01281       fTable = (Char_t *)ptr;
01282       // make sure all link-columns are zero
01283       memset(fTable,0,fSize*fN);
01284    }
01285    return fTable;
01286 }
01287 
01288 //______________________________________________________________________________
01289 void TTable::Browse(TBrowser *b){
01290    // Wrap each table coulumn with TColumnView object to browse.
01291    if (!b) return;
01292    TDataSet::Browse(b);
01293    Int_t nrows = TMath::Min(Int_t(GetNRows()),6);
01294    if (nrows == 0) nrows = 1;
01295    Print(0,nrows);
01296    // Add the table columns to the browser
01297    UInt_t nCol = GetNumberOfColumns();
01298    for (UInt_t i = 0;i<nCol;i++){
01299       TColumnView *view = 0;
01300       UInt_t nDim = GetDimensions(i);
01301       const Char_t *colName = GetColumnName(i);
01302       if (!nDim) { // scalar
01303          // This will cause a small memory leak
01304          // unless TBrowser recognizes kCanDelete bit
01305          if( GetColumnType(i)== kPtr) {
01306             UInt_t offset = GetOffset(i);
01307             TTableMap *m = *(TTableMap **)(((char *)GetArray())+offset);
01308             if (m) {
01309                TString nameMap = "*";
01310                nameMap += m->Table()->GetName();
01311                b->Add(m,nameMap.Data());
01312             }
01313          } else {
01314             view = new TColumnView(GetColumnName(i),this);
01315             view->SetBit(kCanDelete);
01316             b->Add(view,view->GetName());
01317          }
01318       } else {     // array
01319          const UInt_t *indx = GetIndexArray(i);
01320          UInt_t totalSize = 1;
01321          UInt_t k;
01322          for (k=0;k<nDim; k++) totalSize *= indx[k];
01323          for (k=0;k<totalSize;k++) {
01324             TString buffer;
01325             buffer.Form("%s[%d]",colName,k);
01326             view = new TColumnView(buffer,this);
01327             view->SetBit(kCanDelete);
01328             b->Add(view,view->GetName());
01329          }
01330       }
01331    }
01332 }
01333 
01334 //______________________________________________________________________________
01335 void TTable::Clear(Option_t *opt)
01336 {
01337    // Deletes the internal array of this class
01338    // if this object does own its internal table
01339 
01340    if (!fTable) return;
01341    Bool_t dtor = kFALSE;
01342    dtor = opt && (strcmp(opt,gDtorName)==0);
01343    if (!opt || !opt[0] || dtor ) {
01344       if (! TestBit(kIsNotOwn)) {
01345          if (!dtor) ResetMap();
01346          free(fTable);
01347       }
01348       fTable    = 0;
01349       fMaxIndex = 0;
01350       SetfN(0);
01351       return;
01352    }
01353 }
01354 
01355 //______________________________________________________________________________
01356 void TTable::Delete(Option_t *opt)
01357 {
01358    //
01359    // Delete the internal array and free the memory it occupied
01360    // if this object did own this array
01361    //
01362    // Then perform TDataSet::Delete(opt)
01363    Clear(gDtorName);
01364    TDataSet::Delete(opt);
01365 }
01366 
01367 //______________________________________________________________________________
01368 TClass  *TTable::GetRowClass() const
01369 {
01370    //to be documented
01371    TClass *cl = 0;
01372    TTableDescriptor *dsc = GetRowDescriptors();
01373    if (dsc) cl = dsc->RowClass();
01374    else Error("GetRowClass()","Table descriptor of <%s::%s> table lost",
01375              GetName(),GetType());
01376    return cl;
01377 }
01378 
01379 //______________________________________________________________________________
01380 Long_t TTable::GetNRows() const {
01381 // Returns the number of the used rows for the wrapped table
01382    return fMaxIndex;
01383 }
01384 
01385 //______________________________________________________________________________
01386 Long_t TTable::GetRowSize() const {
01387 // Returns the size (in bytes) of one table row
01388    return fSize;
01389 }
01390 
01391 //______________________________________________________________________________
01392 Long_t TTable::GetTableSize() const {
01393 // Returns the number of the allocated rows
01394    return fN;
01395 }
01396 
01397 //______________________________________________________________________________
01398 void TTable::Fit(const char *formula ,const char *varexp, const char *selection,Option_t *option ,Option_t *goption,Int_t nentries, Int_t firstentry)
01399 {
01400 //*-*-*-*-*-*-*-*-*Fit a projected item(s) from a TTable*-*-*-*-*-*-*-*-*-*
01401 //*-*              =======================================
01402 //
01403 //  formula is a TF1 expression.
01404 //
01405 //  See TTable::Draw for explanations of the other parameters.
01406 //
01407 //  By default the temporary histogram created is called htemp.
01408 //  If varexp contains >>hnew , the new histogram created is called hnew
01409 //  and it is kept in the current directory.
01410 //  Example:
01411 //    table.Fit(pol4,"sqrt(x)>>hsqrt","y>0")
01412 //    will fit sqrt(x) and save the histogram as "hsqrt" in the current
01413 //    directory.
01414 //
01415 
01416    TString opt(option);
01417    opt += "goff";
01418 
01419    Draw(varexp,selection,opt,nentries,firstentry);
01420 
01421    TH1 *hfit = gCurrentTableHist;
01422    if (hfit) {
01423       Printf("hname=%s, formula=%s, option=%s, goption=%s\n",hfit->GetName(),formula,option,goption);
01424       // remove bit temporary
01425       Bool_t canDeleteBit = hfit->TestBit(kCanDelete);
01426       if (canDeleteBit)  hfit->ResetBit(kCanDelete);
01427       hfit->Fit(formula,option,goption);
01428       if (TestBit(canDeleteBit))   hfit->SetBit(kCanDelete);
01429    }
01430    else      Printf("ERROR hfit=0\n");
01431 }
01432 
01433 //______________________________________________________________________________
01434 const Char_t *TTable::GetType() const
01435 {
01436 //Returns the type of the wrapped C-structure kept as the TNamed title
01437    return GetTitle();
01438 }
01439 
01440 //______________________________________________________________________________
01441 Bool_t TTable::IsFolder() const {
01442    // return Folder flag to be used by TBrowse object
01443    // The table is a folder if
01444    //  - it has sub-dataset
01445    //    or
01446    //  - GetNRows > 0
01447    return kTRUE; // to provide the "fake" folder bit to workaround TKey::Browse()
01448 
01449 #if 0
01450    // this became useless due TKey::Browse new implementation
01451    return
01452    (fList && fList->Last() ? kTRUE : kFALSE)
01453    ||
01454      (GetNRows() > 0);
01455 #endif
01456 }
01457 
01458 //______________________________________________________________________________
01459 Int_t TTable::NaN()
01460 {
01461 //
01462 // return the total number of the NaN for float/double cells of this table
01463 // Thanks Victor Perevoztchikov
01464 //
01465 
01466    EColumnType code;
01467    char const *cell,*colname,*table;
01468    double word;
01469    int icol,irow,colsize,wordsize,nwords,iword,nerr,offset;
01470 
01471    TTableDescriptor *rowDes = GetRowDescriptors();
01472    assert(rowDes!=0);
01473    table = (const char*)GetArray();
01474 
01475    int ncols = rowDes->GetNumberOfColumns();
01476 
01477    int lrow  = GetRowSize();
01478    int nrows = GetNRows  ();
01479    nerr =0;
01480    for (icol=0; icol < ncols; icol++) {// loop over cols
01481       code = rowDes->GetColumnType(icol);
01482       if (code!=kFloat && code!=kDouble) continue;
01483 
01484       offset   = rowDes->GetOffset    (icol);
01485       colsize  = rowDes->GetColumnSize(icol);
01486       wordsize = rowDes->GetTypeSize  (icol);
01487       nwords = colsize/wordsize;
01488       for (irow=0; irow < nrows; irow++) { //loop over rows
01489          cell = table + offset + irow*lrow;
01490          for (iword=0;iword<nwords; iword++,cell+=wordsize) { //words in col
01491             word = (code==kDouble) ? *(double*)cell : *(float*)cell;
01492             if (TMath::Finite(word))     continue;
01493 //              ERROR FOUND
01494             nerr++; colname = rowDes->GetColumnName(icol);
01495             Warning("NaN"," Table %s.%s.%d\n",GetName(),colname,irow);
01496          }
01497       }
01498    }
01499    return nerr;
01500 }
01501 
01502 //______________________________________________________________________________
01503 TTable *TTable::New(const Char_t *name, const Char_t *type, void *array, UInt_t size)
01504 {
01505   // This static method creates a new TTable object if provided
01506 
01507    TTable *table = 0;
01508    if (type && name) {
01509       TString tableType(type);
01510       TString t = tableType.Strip();
01511 
01512       TString classname("St_");
01513       classname += t;
01514       TClass *cl = TClass::GetClass(classname);
01515       if (cl) {
01516          table = (TTable *)cl->New();
01517          if (table) {
01518             table->SetTablePointer(array);
01519             table->SetName(name);
01520             table->SetfN(size);
01521             table->SetUsedRows(size);
01522          }
01523       }
01524    }
01525    return table;
01526 }
01527 //______________________________________________________________________________
01528 Bool_t TTable::OutOfBoundsError(const char *where, Int_t i) const
01529 {
01530    // Generate an out-of-bounds error. Always returns false.
01531    Error(where, "index %d out of bounds (size: %d, this: 0x%lx)", i, fN, (Long_t)this);
01532    return kFALSE;
01533 }
01534 //______________________________________________________________________________
01535 Char_t *TTable::Print(Char_t *strbuf,Int_t lenbuf) const
01536 {
01537    // Create IDL table defintion (to be used for XDF I/O)
01538    Int_t iOut = 0;
01539 
01540    TTableDescriptor *dscT = GetRowDescriptors();
01541    if (!dscT ) {
01542       Error("Print"," No dictionary entry for <%s> structure", GetTitle());
01543       if (lenbuf>0) iOut += snprintf(strbuf,lenbuf," *** Errror ***");
01544       return strbuf;
01545    }
01546 
01547    TROOT::IndentLevel();
01548    if (lenbuf>0) {
01549       // cut of the "_st" suffix
01550       Char_t *typenam =  new Char_t [strlen(dscT->GetName())+1];
01551       strlcpy(typenam,dscT->GetName(),strlen(dscT->GetName())+1);
01552       // look for the last "_"
01553       Char_t *last = strrchr(typenam,'_');
01554       // Check whether it is "_st"
01555       Char_t *eon = 0;
01556       if (last) eon = strstr(last,"_st");
01557       // Cut it off if any
01558       if (eon) *eon = '\0';
01559       iOut += snprintf(strbuf+iOut,lenbuf-iOut,"struct %s {",typenam);
01560       delete [] typenam;
01561    } else {
01562       cout << "struct " << dscT->GetName() << " {" << endl;
01563    }
01564 
01565    TTableDescriptor::iterator dsc  = dscT->begin();
01566    TTableDescriptor::iterator dscE = dscT->end();
01567    TDataSetIter nextComment(dscT->MakeCommentField(kFALSE));
01568    for (;dsc != dscE; dsc++) {
01569       TROOT::IndentLevel();
01570       TString name = GetTypeName(EColumnType((*dsc).fType));
01571       if (lenbuf>0) {
01572          // convert C type names to CORBA type names
01573          name.ReplaceAll("unsigned char","octet");
01574          name.ReplaceAll("int","long");
01575          iOut += snprintf(strbuf+iOut,lenbuf-iOut," %s %s",name.Data(),(*dsc).fColumnName);
01576       } else
01577          cout << '\t'<< name.Data() << '\t'<< (*dsc).fColumnName;
01578 
01579       Int_t indx;
01580       Int_t dim = (*dsc).fDimensions;
01581       for  (indx = 0; indx < dim; indx++) {
01582          if (lenbuf>0)
01583             iOut += snprintf(strbuf+iOut,lenbuf-iOut,"[%d]",(*dsc).fIndexArray[indx]);
01584          else
01585             cout <<  "[" << dec << (*dsc).fIndexArray[indx]<<"]";
01586       }
01587       // print comment if any
01588       TDataSet *nxc = nextComment();
01589       if (lenbuf>0)
01590          iOut += snprintf(strbuf+iOut,lenbuf-iOut, ";");
01591       else {
01592          const char *title = nxc ? nxc->GetTitle() : " ";
01593          cout << ";\t//" << title << endl;
01594       }
01595    } /* dsc */
01596 
01597    TROOT::IndentLevel();
01598    if (lenbuf>0)
01599       iOut += snprintf(strbuf+iOut,lenbuf-iOut, "}");
01600    else
01601       cout << "}" << endl;
01602    return strbuf;
01603 }
01604 
01605 //______________________________________________________________________________
01606 const Char_t *TTable::PrintHeader() const
01607 {
01608   // Print general table inforamtion
01609    cout << endl << " ---------------------------------------------------------------------------------------" << endl
01610         <<  " " << Path()
01611                 <<"  Allocated rows: "<<fN
01612                 <<"\t Used rows: "<<fMaxIndex
01613                 <<"\t Row size: "      << fSize << " bytes"
01614         <<endl;
01615    return 0;
01616 }
01617 
01618 //______________________________________________________________________________
01619 const Char_t *TTable::Print(Int_t row, Int_t rownumber, const Char_t *, const Char_t *) const
01620 {
01621   //const Char_t *TTable::Print(Int_t row, Int_t rownumber, const Char_t *colfirst, const Char_t *collast) const
01622   //
01623   //  Print the contents of internal table per COLUMN.
01624   //
01625   //  row       - the index of the first row to print (counting from ZERO)
01626   //  rownumber - the total number of rows to print out (=10 by default)
01627   //
01628   //  (No use !) Char_t *colfirst, *collast - the names of the first/last
01629   //                                          to print out (not implemented yet)
01630   //
01631   //--------------------------------------------------------------
01632    // Check bounds and adjust it
01633    Int_t const width = 8;
01634    Int_t rowStep = 10; // The maximun values to print per line
01635    Int_t rowNumber = rownumber;
01636    if (row  > Int_t(GetSize()) || GetSize() == UInt_t(0))  {
01637       PrintHeader();
01638       cout  << " ======================================================================================" << endl
01639            << "   There are " << GetSize() << " allocated rows for this table only"                     << endl
01640              << " ======================================================================================" << endl;
01641       return 0;
01642    }
01643    if (rowNumber > Int_t(GetSize()-row)) rowNumber = GetSize()-row;
01644    if (!rowNumber) return 0;
01645    rowStep = TMath::Min(rowStep,rowNumber);
01646 
01647    Int_t cdate = 0;
01648    Int_t ctime = 0;
01649    UInt_t *cdatime = 0;
01650    Bool_t isdate = kFALSE;
01651 
01652    TTableDescriptor *dscT = GetRowDescriptors();
01653    if (!dscT ) return 0;
01654 
01655    //  3. Loop by "rowStep x lines"
01656 
01657    const Char_t  *startRow = (const Char_t *)GetArray() + row*GetRowSize();
01658    Int_t rowCount = rowNumber;
01659    Int_t thisLoopLenth = 0;
01660    const Char_t  *nextRow = 0;
01661    while (rowCount) {
01662       PrintHeader();
01663       if  (GetNRows() == 0) {// to Print empty table header
01664          cout  << " ======================================================================================" << endl
01665                << "   There is NO filled row in this table"                                                 << endl
01666                << " ======================================================================================" << endl;
01667          return 0;
01668       }
01669       cout << " Table: " << dscT->GetName()<< "\t";
01670       for (Int_t j = row+rowNumber-rowCount; j<row+rowNumber-rowCount+rowStep && j < row+rowNumber ;j++) {
01671          Int_t hW = width-2;
01672          if (j>=10) hW -= (int)TMath::Log10(float(j))-1;
01673          cout  << setw(hW) << "["<<j<<"]";
01674          cout  << " :" ;
01675       }
01676       cout << endl
01677       <<       " ======================================================================================" << endl;
01678       TTableDescriptor::iterator member = dscT->begin();
01679       TTableDescriptor::iterator   dscE = dscT->end();
01680       TDataSetIter nextComment(dscT->MakeCommentField(kFALSE));
01681 
01682       for (; member != dscE; member++){
01683          TString membertype = GetTypeName(EColumnType((*member).fType));
01684          isdate = kFALSE;
01685          if (strcmp((*member).fColumnName,"fDatime") == 0 && membertype == "UInt_t")
01686                                                                                    isdate = kTRUE;
01687          cout << membertype.Data();
01688 
01689          // Add the dimensions to "array" members
01690          Int_t dim = (*member).fDimensions;
01691          Int_t indx = 0;
01692          UInt_t *arrayLayout = 0;
01693          if (dim) {
01694             arrayLayout = new UInt_t[dim];
01695             memset(arrayLayout,0,dim*sizeof(Int_t));
01696          }
01697          Int_t arrayLength  = 1;
01698          while (indx < dim ){ // Take in account the room this index will occupy
01699             arrayLength *= (*member).fIndexArray[indx];
01700             indx++;
01701          }
01702          // Encode data value or pointer value
01703          Int_t offset = (*member).fOffset;
01704          Int_t thisStepRows;
01705          thisLoopLenth = TMath::Min(rowCount,rowStep);
01706          Int_t indexOffset;
01707          Bool_t breakLoop = kFALSE;
01708 
01709          for (indexOffset=0; indexOffset < arrayLength && !breakLoop; indexOffset++) {
01710             nextRow = startRow;
01711 
01712             if (!indexOffset) cout << "\t" << (*member).fColumnName;
01713             else              cout << "\t" << setw(strlen((*member).fColumnName)) << " ";
01714 
01715             if (dim) {
01716                for (Int_t i=0;i<dim;i++) cout << "["<<dec<<arrayLayout[i]<<"]";
01717                ArrayLayout(arrayLayout,(*member).fIndexArray,dim);
01718             }
01719             cout << "\t";
01720             if ( strlen((*member).fColumnName)+3*dim < 8) cout << "\t";
01721 
01722             for (thisStepRows = 0;thisStepRows < thisLoopLenth; thisStepRows++,nextRow += GetRowSize()) {
01723                const char *pointer = nextRow + offset  + indexOffset*(*member).fTypeSize;
01724                if (isdate) {
01725                   cdatime = (UInt_t*)pointer;
01726                   TDatime::GetDateTime(cdatime[0],cdate,ctime);
01727                   cout << cdate << "/" << ctime;
01728                } else if ((*member).fType == kChar && dim == 1) {
01729                   char charbuffer[11];
01730                   strlcpy(charbuffer,pointer,TMath::Min(10,arrayLength)+1);
01731                   charbuffer[10] = 0;
01732                   cout << "\"" << charbuffer;
01733                   if (arrayLength > 10)
01734                      cout << " . . . ";
01735                   cout << "\"";
01736                   breakLoop = kTRUE;
01737                } else {
01738                   AsString((void *)pointer,EColumnType((*member).fType),width,cout);
01739                   cout << " :";
01740                }
01741             }
01742             // Encode  the column's comment
01743             if (indexOffset==0) {
01744                TDataSet *nxc = nextComment();
01745                cout << " " << (const char *)(nxc ? nxc->GetTitle() : "no comment");
01746             }
01747             cout << endl;
01748          }
01749          if (arrayLayout) delete [] arrayLayout;
01750       }
01751       rowCount -= thisLoopLenth;
01752       startRow  = nextRow;
01753    }
01754    cout << "---------------------------------------------------------------------------------------" << endl;
01755    return 0;
01756 }
01757 //______________________________________________________________________________
01758 void TTable::PrintContents(Option_t *) const
01759 {
01760    //to be documented
01761    TDataSet::PrintContents();
01762    TROOT::IndentLevel();
01763    Printf("\tclass %s: public TTable\t --> Allocated rows: %d\t Used rows: %d\t Row size: %d bytes\n",
01764          IsA()->GetName(),int(fN),int(fMaxIndex),int(fSize));
01765 
01766 }
01767 
01768 //______________________________________________________________________________
01769 void TTable::Project(const char *hname, const char *varexp, const char *selection, Option_t *option,Int_t nentries, Int_t firstentry)
01770 {
01771 //*-*-*-*-*-*-*-*-*Make a projection of a TTable using selections*-*-*-*-*-*-*
01772 //*-*              =============================================
01773 //
01774 //   Depending on the value of varexp (described in Draw) a 1-D,2-D,etc
01775 //   projection of the TTable will be filled in histogram hname.
01776 //   Note that the dimension of hname must match with the dimension of varexp.
01777 //
01778 
01779    TString var;
01780    var.Form("%s>>%s",varexp,hname);
01781 
01782    TString opt(option);
01783    opt += "goff";
01784 
01785    Draw(var,selection,opt,nentries,firstentry);
01786 
01787    delete [] var;
01788    delete [] opt;
01789 }
01790 
01791 //______________________________________________________________________________
01792 Int_t TTable::Purge(Option_t *opt)
01793 {
01794    // Shrink the table to free the unused but still allocated rows
01795    ReAllocate();
01796    return TDataSet::Purge(opt);
01797 }
01798 
01799 //______________________________________________________________________________
01800 void TTable::SavePrimitive(ostream &out, Option_t * /*= ""*/)
01801 {
01802 //   Save a primitive as a C++ statement(s) on output stream "out".
01803    UInt_t arrayLayout[10],arraySize[10];
01804    const unsigned char *pointer=0,*startRow=0;
01805    int i,rowCount;unsigned char ic;
01806 
01807    out << "TDataSet *CreateTable() { " << endl;
01808 
01809    Int_t rowNumber =  GetNRows();
01810    TTableDescriptor *dscT = GetRowDescriptors();
01811 
01812 //                      Is anything Wrong??
01813    if (!rowNumber || !dscT ) {//
01814       out << "// The output table was bad-defined!" << endl
01815           << " fprintf(stderr, \"Bad table found. Please remove me\\n\");" << endl
01816           << " return 0; } "    << endl;
01817       return;
01818    }
01819 
01820    startRow = (const UChar_t *)GetArray();
01821    assert(startRow!=0);
01822 
01823    const Char_t *rowId = "row";
01824    const Char_t *tableId = "tableSet";
01825 
01826 //                      Generate the header
01827 
01828    const char *className = IsA()->GetName();
01829 
01830    out << "// -----------------------------------------------------------------" << endl;
01831    out << "// "   << Path()
01832        << " Allocated rows: "<< rowNumber
01833        <<"  Used rows: "<<      rowNumber
01834        <<"  Row size: " << fSize << " bytes"                 << endl;
01835    out << "// "  << " Table: " << dscT->GetName()<<"[0]--> "
01836        << dscT->GetName()<<"["<<rowNumber-1 <<"]"            << endl;
01837    out << "// ====================================================================" << endl;
01838    out << "// ------  Test whether this table share library was loaded ------"      << endl;
01839    out << "  if (!TClass::GetClass(\"" << className << "\")) return 0;"    << endl;
01840    out <<    dscT->GetName() << " " << rowId << ";" << endl
01841        <<  className << " *" << tableId << " = new "
01842        <<  className
01843        << "(\""<<GetName()<<"\"," << GetNRows() << ");" << endl
01844        << "//" <<endl ;
01845 
01846 //              Row loop
01847    TDataSetIter nextComment(dscT->MakeCommentField(kFALSE));
01848    for (rowCount=0;rowCount<rowNumber; rowCount++,startRow += fSize, nextComment.Reset()) {     //row loop
01849       out << "memset(" << "&" << rowId << ",0," << tableId << "->GetRowSize()" << ");" << endl ;
01850 
01851 //              Member loop
01852    TTableDescriptor::iterator member  = dscT->begin();
01853    TTableDescriptor::iterator   dscE  = dscT->end();
01854    for (; member != dscE; member++) {  //LOOP over members
01855       TString memberType = GetTypeName(EColumnType((*member).fType));
01856       TString memberName((*member).fColumnName);
01857 
01858        // Encode  the column's comment
01859       TDataSet *nxc = nextComment();
01860       TString memberTitle(nxc ? nxc->GetTitle() : "no comment");
01861 
01862       Int_t offset = (*member).fOffset;
01863       int mayBeName = 0;
01864       if (memberName.Index("name",0,TString::kIgnoreCase)>=0) mayBeName=1999;
01865       if (memberName.Index("file",0,TString::kIgnoreCase)>=0) mayBeName=1999;
01866       int typeSize = (*member).fTypeSize;
01867 
01868 //              Add the dimensions to "array" members
01869       Int_t dim = (*member).fDimensions;
01870       if (dim) memset(arrayLayout,0,dim*sizeof(Int_t));
01871       Int_t arrayLength  = 1;
01872       for (int indx=0;indx < dim ;indx++){
01873          arraySize[indx] =  (*member).fIndexArray[indx];;
01874          arrayLength *= arraySize[indx];
01875       }
01876 
01877 //                      Special case, character array
01878       int charLen = (memberType.CompareTo("char")==0);
01879       if (charLen) {    //Char case
01880          charLen=arrayLength;
01881          pointer = startRow + offset;
01882 //                      Actual size of char array
01883          if (mayBeName) {
01884             charLen = strlen((const char*)pointer)+1;
01885             if (charLen>arrayLength) charLen = arrayLength;
01886          } else {
01887             for(;charLen && !pointer[charLen-1];charLen--){;}
01888             if (!charLen) charLen=1;
01889          }
01890 
01891          out << " memcpy(&" << rowId << "." << (const char*)memberName;
01892          out << ",\"";
01893          for (int ii=0; ii<charLen;ii++) {
01894             ic = pointer[ii];
01895             if (ic && (isalnum(ic)
01896              || strchr("!#$%&()*+-,./:;<>=?@{}[]_|~",ic))) {//printable
01897                out << ic;
01898             } else {                                      //nonprintable
01899                out << "\\x" << setw(2) << setfill('0') << hex << (unsigned)ic ;
01900                out << setw(1) << setfill(' ') << dec;
01901             }
01902          }
01903          out << "\"," << dec << charLen << ");";
01904          out << "// " << (const char*)memberTitle << endl;
01905          continue;
01906       } //EndIf of char case
01907 
01908 //                      Normal member
01909       Int_t indexOffset;
01910       for (indexOffset=0; indexOffset < arrayLength ; indexOffset++) {//array loop
01911          out << setw(3) << " " ;
01912          out << " " << rowId << "." << (const char*)memberName;
01913 
01914          if (dim) {
01915             for (i=0;i<dim;i++) {out << "["<<dec<<arrayLayout[i]<<"]";}
01916             ArrayLayout(arrayLayout,arraySize,dim);}
01917 
01918 //                      Generate "="
01919             out << "\t = ";
01920 
01921             pointer = startRow + offset  + indexOffset*typeSize;
01922 
01923             AsString((void *)pointer,EColumnType((*member).fType),10,out);
01924 
01925 //                      Encode data member title
01926             if (indexOffset==0)  out << "; // " << (const char*)memberTitle;
01927             out << ";" << endl;
01928          }//end array loop
01929       }//end of member loop
01930 
01931       out << tableId << "->AddAt(&" << rowId <<");" << endl;
01932 
01933    }//end of row loop
01934    out << "// ----------------- end of code ---------------" << endl
01935        << " return (TDataSet *)tableSet;" << endl
01936        << "}"  << endl;
01937    return;
01938 }
01939 
01940 //______________________________________________________________________________
01941 void TTable::Set(Int_t n)
01942 {
01943    // Set array size of TTable object to n longs. If n<0 leave array unchanged.
01944    if (n < 0) return;
01945    if (fN != n)  Clear();
01946    SetfN(n);
01947    if (fN == 0) return;
01948    Create();
01949    if (TTable::GetNRows()) Reset();
01950 }
01951 //______________________________________________________________________________
01952 void TTable::SetTablePointer(void *table)
01953 {
01954    //to be documented
01955    if (fTable) free(fTable);
01956    fTable = (Char_t *)table;
01957 }
01958 
01959 //______________________________________________________________________________
01960 void TTable::SetType(const char *const type)
01961 {
01962    //to be documented
01963    SetTitle(type);
01964 }
01965 
01966 //______________________________________________________________________________
01967 static Char_t *GetExpressionFileName()
01968 {
01969    // Create a name of the file in the temporary directory if any
01970    const Char_t *tempDirs =  gSystem->Getenv("TEMP");
01971    if (!tempDirs)  tempDirs =  gSystem->Getenv("TMP");
01972    if (!tempDirs) tempDirs = "/tmp";
01973    if (gSystem->AccessPathName(tempDirs)) tempDirs = ".";
01974    if (gSystem->AccessPathName(tempDirs)) return 0;
01975    TString fileName;
01976    fileName.Form("Selection.C.%d.tmp",gSystem->GetPid());
01977    return  gSystem->ConcatFileName(tempDirs,fileName.Data());
01978 }
01979 
01980 //______________________________________________________________________________
01981 Char_t *TTable::MakeExpression(const Char_t *expressions[],Int_t nExpressions)
01982 {
01983   // Create CINT macro to evaluate the user-provided expresssion
01984   // Expression may contains:
01985   //   -  the table columen names
01986   //   - 2 meta names: i$ - the current column index,
01987   //                   n$ - the total table size provided by TTable::GetNRows() method
01988   //
01989   // return the name of temporary file with the current expressions
01990   //
01991    const Char_t *typeNames[] = {"NAN","float", "int",  "long",  "short",         "double"
01992                                 ,"unsigned int","unsigned long", "unsigned short","unsigned char"
01993                                 ,"char", "TTableMap &"};
01994    const char *resID     = "results";
01995    const char *addressID = "address";
01996    Char_t *fileName = GetExpressionFileName();
01997    if (!fileName) {
01998       Error("MakeExpression","Can not create a temporary file");
01999       return 0;
02000    }
02001 
02002    ofstream str;
02003    str.open(fileName);
02004    if (str.bad() ) {
02005       Error("MakeExpression","Can not open the temporary file <%s>",fileName);
02006       delete [] fileName;
02007       return 0;
02008    }
02009 
02010    TTableDescriptor *dsc = GetRowDescriptors();
02011    const tableDescriptor_st *descTable  = dsc->GetTable();
02012    // Create function
02013    str << "void SelectionQWERTY(float *"<<resID<<", float **"<<addressID<< ", int& i$, int& n$ )"   << endl;
02014    str << "{"                                                        << endl;
02015    int i = 0;
02016    for (i=0; i < dsc->GetNRows(); i++,descTable++ ) {
02017       // Take the column name
02018       const Char_t *columnName = descTable->fColumnName;
02019       const Char_t *type = 0;
02020       // First check whether we do need this column
02021       for (Int_t exCount = 0; exCount < nExpressions; exCount++) {
02022          if (expressions[exCount] && expressions[exCount][0] && strstr(expressions[exCount],columnName)) goto LETSTRY;
02023       }
02024       continue;
02025 LETSTRY:
02026       Bool_t isScalar = !(descTable->fDimensions);
02027       Bool_t isFloat = descTable->fType == kFloat;
02028       type = typeNames[descTable->fType];
02029       str << type << " ";
02030       if (!isScalar)  str << "*";
02031 
02032       str << columnName << " = " ;
02033       if (isScalar)   str << "*(";
02034       if (!isFloat)   str << "(" << type << "*)";
02035       str << addressID << "[" << i << "]";
02036       if (isScalar)   str << ")" ;
02037       str << ";" << endl;
02038    }
02039    // Create expressions
02040    for (i=0; i < nExpressions; i++ ) {
02041       if (expressions[i] && expressions[i][0])
02042          str << " "<<resID<<"["<<i<<"]=(float)(" << expressions[i] << ");"  << endl;
02043 //      if (i == nExpressions-1 && i !=0 )
02044 //          str  << "  if ("<<resID<<"["<<i<<"] == 0){ return; }" << endl;
02045    };
02046    str << "}" << endl;
02047    str.close();
02048    // Create byte code and check syntax
02049    if (str.good()) return fileName;
02050    delete [] fileName;
02051    return 0;
02052 }
02053 
02054 //______________________________________________________________________________
02055 void TTable::Reset(Int_t c)
02056 {
02057    // Fill the entire table with byte "c" ;
02058    ///     c=0 "be default"
02059    if (fTable) {
02060       ResetMap(kTRUE);
02061       ::memset(fTable,c,fSize*fN);
02062       if (c) ResetMap(kFALSE);
02063    }
02064 }
02065 
02066 //______________________________________________________________________________
02067 void TTable::ResetMap(Bool_t wipe)
02068 {
02069    // Clean all filled columns with the pointers to TTableMap
02070    // if any
02071    //  wipe = kTRUE - delete all object the Map's point to
02072    //         kFALSE - zero pointer, do not call "delete" though
02073    piterator links     = pbegin();
02074    piterator lastLinks = pend();
02075    for (;links != lastLinks;links++) {
02076       TTableMap **mp = (TTableMap **)(*links);
02077       if (wipe) delete *mp;
02078       *mp = 0;
02079    }
02080 }
02081 //______________________________________________________________________________
02082 void TTable::Set(Int_t n, Char_t *array)
02083 {
02084    // Set array size of TTable object to n longs and copy array.
02085    // If n<0 leave array unchanged.
02086 
02087    if (n < 0) return;
02088    if (fN < n) Clear();
02089 
02090    SetfN(n);
02091 
02092    if (fN == 0) return;
02093    Create();
02094    CopyStruct(fTable,array);
02095    fMaxIndex = n;
02096 }
02097 
02098 //_______________________________________________________________________
02099 void TTable::StreamerTable(TBuffer &b,Version_t version)
02100 {
02101    // Stream an object of class TTable.
02102    if (b.IsReading()) {
02103       TDataSet::Streamer(b);
02104       b >> fN;
02105       StreamerHeader(b,version);
02106       //   Create a table to fit nok rows
02107       Set(fMaxIndex);
02108    } else {
02109       TDataSet::Streamer(b);
02110       b << fN;
02111       StreamerHeader(b,version);
02112    }
02113 }
02114 
02115 //_______________________________________________________________________
02116 void TTable::StreamerHeader(TBuffer &b, Version_t version)
02117 {
02118    // Read "table parameters first"
02119    if (b.IsReading()) {
02120       Long_t rbytes;
02121       if (version) { }   // version to remove compiler warning
02122 #ifdef __STAR__
02123       if (version < 3) {
02124          // skip obsolete  STAR fields (for the sake of the backward compatibility)
02125          //   char name[20];   /* table name */
02126          //   char type[20];   /* table type */
02127          //   long maxlen;     /* # rows allocated */
02128          long len = b.Length() + (20+4) + (20+4) + 4;
02129          b.SetBufferOffset(len);
02130       }
02131 #endif
02132       b >> fMaxIndex;         // fTableHeader->nok;          /* # rows filled */
02133       b >> rbytes;            /* number of bytes per row */
02134       if (GetRowSize() == -1) fSize = rbytes;
02135       if (rbytes - GetRowSize()) {
02136          Warning("StreamerHeader","Schema evolution warning: row size mismatch: expected %ld, read %ld bytes\n",GetRowSize(),rbytes);
02137       }
02138 
02139 #ifdef __STAR__
02140       if (version < 3) {
02141          // skip obsolete  STAR fields (for the sake of the backward compatibility)
02142          //    long dsl_pointer;  /* swizzled (DS_DATASET_T*) */
02143          //    long data_pointer; /* swizzled (char*) */
02144          long len = b.Length() + (4) + (4);
02145          b.SetBufferOffset(len);
02146       }
02147 #endif
02148    } else {
02149       b << fMaxIndex;         //fTableHeader->nok;          /* # rows filled */
02150       b << fSize;             //  fTableHeader->rbytes;     /* number of bytes per row */
02151    }
02152 }
02153 //_______________________________________________________________________
02154 Int_t TTable::SetfN(Long_t len)
02155 {
02156    //to be documented
02157    fN = len;
02158    return fN;
02159 }
02160 //____________________________________________________________________________
02161 #ifdef StreamElelement
02162 #define __StreamElelement__ StreamElelement
02163 #undef StreamElelement
02164 #endif
02165 
02166 #define StreamElementIn(type)  case TTableDescriptor::_NAME2_(k,type):        \
02167  if (evolutionOn) {                                  \
02168      if (nextCol->fDimensions)  {                    \
02169        if (nextCol->fOffset != UInt_t(-1)) {         \
02170           R__b.ReadFastArray((_NAME2_(type,_t) *)(row+nextCol->fOffset),nextCol->fSize/sizeof(_NAME2_(type,_t)));   \
02171        } else {                                        \
02172            _NAME2_(type,_t) *readPtrV = new _NAME2_(type,_t)[nextCol->fSize/sizeof(_NAME2_(type,_t))];              \
02173            R__b.ReadFastArray((_NAME2_(type,_t) *)(row+nextCol->fOffset),nextCol->fSize/sizeof(_NAME2_(type,_t)));  \
02174            delete [] readPtrV;                       \
02175            readPtrV = 0;                             \
02176        }                                             \
02177      }                                               \
02178      else  {                                         \
02179        _NAME2_(type,_t) skipBuffer;                  \
02180        _NAME2_(type,_t) *readPtr =  (_NAME2_(type,_t) *)(row+nextCol->fOffset); \
02181        if (nextCol->fOffset == UInt_t(-1)) readPtr = &skipBuffer;               \
02182        R__b >> *readPtr;                             \
02183      }                                               \
02184  } else {                                            \
02185    if (nextCol->fDimensions)  {                      \
02186      R__b.ReadFastArray  ((_NAME2_(type,_t) *)(row+nextCol->fOffset),nextCol->fSize/sizeof(_NAME2_(type,_t)));  \
02187    } else                                                       \
02188      R__b >> *(_NAME2_(type,_t) *)(row+nextCol->fOffset);       \
02189  }                                                              \
02190  break
02191 
02192 #define StreamElementOut(type) case TTableDescriptor::_NAME2_(k,type):    \
02193  if (nextCol->fDimensions)                                    \
02194     R__b.WriteFastArray((_NAME2_(type,_t) *)(row+nextCol->fOffset), nextCol->fSize/sizeof(_NAME2_(type,_t))); \
02195  else                                                         \
02196     R__b << *(_NAME2_(type,_t) *)(row+nextCol->fOffset);      \
02197  break
02198 
02199 //______________________________________________________________________________
02200 TTableDescriptor  *TTable::GetRowDescriptors() const
02201 {
02202    //to be documented
02203    TTableDescriptor *dsc = 0;
02204    if (IsA()) dsc = GetDescriptorPointer();
02205    if (!dsc) {
02206       Error("GetRowDescriptors()","%s has no dictionary !",GetName());
02207       dsc = GetTableDescriptors();
02208       ((TTableDescriptor *)this)->SetDescriptorPointer(dsc);
02209    }
02210    return dsc;
02211 }
02212 //______________________________________________________________________________
02213 TTableDescriptor *TTable::GetDescriptorPointer() const
02214 {
02215    //to be documented
02216    assert(0);
02217    return 0;
02218 }
02219 
02220 //______________________________________________________________________________
02221 void TTable::SetDescriptorPointer(TTableDescriptor *)
02222 {
02223    //to be documented
02224    assert(0);
02225 }
02226 
02227 //______________________________________________________________________________
02228 void TTable::Streamer(TBuffer &R__b)
02229 {
02230    // Stream an array of the "plain" C-structures
02231    TTableDescriptor *ioDescriptor = GetRowDescriptors();
02232    TTableDescriptor *currentDescriptor = ioDescriptor;
02233    Version_t R__v = 0;
02234    if (R__b.IsReading()) {
02235       // Check whether the file is the "obsolete" one
02236       R__v = R__b.ReadVersion();
02237       Bool_t evolutionOn = kFALSE;
02238       if (R__v>=2) {
02239          if (IsA() != TTableDescriptor::Class()) {
02240             if (R__v>3) {
02241                R__b >> ioDescriptor;
02242             } else {  // backward compatibility
02243                ioDescriptor =  new TTableDescriptor();
02244                ioDescriptor->Streamer(R__b);
02245             }
02246             if (!currentDescriptor) {
02247                currentDescriptor = ioDescriptor;
02248                SetDescriptorPointer(currentDescriptor);
02249             }
02250             if (currentDescriptor->fSecondDescriptor != ioDescriptor) {
02251                // Protection against of memory leak.
02252                delete currentDescriptor->fSecondDescriptor;
02253                currentDescriptor->fSecondDescriptor = ioDescriptor;
02254             }
02255 
02256             // compare two descriptors
02257             evolutionOn = (Bool_t)ioDescriptor->UpdateOffsets(currentDescriptor);
02258          }
02259       }
02260       TTable::StreamerTable(R__b,R__v);
02261       if (fMaxIndex <= 0) return;
02262       char *row= fTable;
02263       Int_t maxColumns = ioDescriptor->NumberOfColumns();
02264       Int_t rowSize = GetRowSize();
02265       if (evolutionOn) Reset(0); // Clean table
02266       for (Int_t indx=0;indx<fMaxIndex;indx++,row += rowSize) {
02267          tableDescriptor_st *nextCol = ioDescriptor->GetTable();
02268          for (Int_t colCounter=0; colCounter < maxColumns; colCounter++,nextCol++) {
02269             // Stream one table row supplied
02270             switch(nextCol->fType) {
02271                StreamElementIn(Float);
02272                StreamElementIn(Int);
02273                StreamElementIn(Long);
02274                StreamElementIn(Short);
02275                StreamElementIn(Double);
02276                StreamElementIn(UInt);
02277                StreamElementIn(ULong);
02278                StreamElementIn(UChar);
02279                StreamElementIn(Char);
02280                StreamElementIn(Bool);
02281                case TTableDescriptor::kPtr: {
02282                   Ptr_t readPtr;
02283                   R__b >> readPtr;
02284                   if (evolutionOn) {
02285                      // TTableMap skipBuffer;
02286                      //  R__b >> readPtr;
02287                      if (nextCol->fOffset == UInt_t(-1)) delete readPtr; // skip this member
02288                      else *(Ptr_t *)(row+nextCol->fOffset) = readPtr;
02289                   } else {
02290                      *(Ptr_t *)(row+nextCol->fOffset) = readPtr;
02291                }
02292                break;
02293                                        }
02294                default:
02295                break;
02296             };
02297          }
02298       }
02299    } else {
02300       TSeqCollection *save = fList;
02301       R__b.WriteVersion(TTable::IsA());
02302 
02303       //      if (Class_Version()==2)
02304       if (IsA() != TTableDescriptor::Class()) {
02305          if ( Class_Version()>3 ) {
02306             R__b << ioDescriptor;
02307          } else {  // backward compatibility
02308             ioDescriptor->Streamer(R__b);
02309          }
02310       } else {
02311          if ( Class_Version()<=3 ) fList = 0;
02312       }
02313 
02314       TTable::StreamerTable(R__b);
02315       if (fMaxIndex <= 0) return;
02316       char *row= fTable;
02317       Int_t maxColumns = ioDescriptor->NumberOfColumns();
02318       Int_t rowSize = GetRowSize();
02319       for (Int_t indx=0;indx<fMaxIndex;indx++,row += rowSize) {
02320          tableDescriptor_st *nextCol = ioDescriptor->GetTable();
02321          for (Int_t colCounter=0; colCounter < maxColumns; colCounter++,nextCol++) {
02322             // Stream one table row supplied
02323             switch(nextCol->fType) {
02324                StreamElementOut(Float);
02325                StreamElementOut(Int);
02326                StreamElementOut(Long);
02327                StreamElementOut(Short);
02328                StreamElementOut(Double);
02329                StreamElementOut(UInt);
02330                StreamElementOut(ULong);
02331                StreamElementOut(UChar);
02332                StreamElementOut(Char);
02333                StreamElementOut(Bool);
02334                case TTableDescriptor::kPtr:
02335                   R__b << *(Ptr_t *)(row+nextCol->fOffset);
02336                   break;
02337                default:
02338                   break;
02339             };
02340          }
02341       }
02342       fList = save;
02343    }
02344 }
02345 #ifdef __StreamElelement__
02346 #define StreamElelement __StreamElelement__
02347 #undef __StreamElelement__
02348 #endif
02349 
02350 //_______________________________________________________________________
02351 void TTable::Update()
02352 {
02353    //to be documented
02354 }
02355 
02356 //_______________________________________________________________________
02357 void TTable::Update(TDataSet *set, UInt_t opt)
02358 {
02359  // Kill the table current data
02360  // and adopt those from set
02361    if (set->HasData()) {
02362       // Check whether the new table has the same type
02363       if (strcmp(GetTitle(),set->GetTitle()) == 0 ) {
02364          TTable *table =  (TTable *)set;
02365          Adopt(table->GetSize(),table->GetArray());
02366          // Adopt can not distniguish the "allocated" and "used"
02367          // rows,
02368          // correct the corrupted number of the "used" rows
02369          SetUsedRows(table->GetNRows());
02370          // mark that object lost the STAF table and can not delete it anymore
02371          table->SetBit(kIsNotOwn);
02372          // mark we took over of this STAF table
02373          ResetBit(kIsNotOwn);
02374       } else
02375          Error("Update",
02376              "This table is <%s> but the updating one has a wrong type <%s>",GetTitle(),set->GetTitle());
02377    }
02378    TDataSet::Update(set,opt);
02379 }
02380 //_______________________________________________________________________
02381 const char *TTable::TableDictionary(const char *className,const char *structName,TTableDescriptor *&ColDescriptors)
02382 {
02383    // Query the TClass instance for the C-stucture dicitonary
02384    // This method is to be used  with TableImp CPP macro (see $ROOTSYS/table/inc/Ttypes.h
02385    if (className){/*NotUsed*/};
02386    TClass *r = TClass::GetClass(structName,1);
02387    ColDescriptors = new TTableDescriptor(r);
02388    return structName;
02389 }
02390 
02391 
02392  //  ----   Table descriptor service   ------
02393 Int_t        TTable::GetColumnIndex(const Char_t *columnName) const {return GetRowDescriptors()->ColumnByName(columnName);}
02394 const Char_t *TTable::GetColumnName(Int_t columnIndex) const {return GetRowDescriptors()->ColumnName(columnIndex); }
02395 const UInt_t *TTable::GetIndexArray(Int_t columnIndex) const {return GetRowDescriptors()->IndexArray(columnIndex); }
02396 UInt_t       TTable::GetNumberOfColumns()              const {return GetRowDescriptors()->NumberOfColumns();       }
02397 
02398 UInt_t       TTable::GetOffset(Int_t columnIndex)      const {return GetRowDescriptors()->Offset(columnIndex); }
02399 Int_t        TTable::GetOffset(const Char_t *columnName) const {return GetRowDescriptors()->Offset(columnName); }
02400 
02401 UInt_t       TTable::GetColumnSize(Int_t columnIndex)  const {return GetRowDescriptors()->ColumnSize(columnIndex); }
02402 Int_t        TTable::GetColumnSize(const Char_t *columnName) const {return GetRowDescriptors()->ColumnSize(columnName); }
02403 
02404 UInt_t       TTable::GetTypeSize(Int_t columnIndex)    const {return GetRowDescriptors()->TypeSize(columnIndex); }
02405 Int_t        TTable::GetTypeSize(const Char_t *columnName) const {return GetRowDescriptors()->TypeSize(columnName); }
02406 
02407 UInt_t       TTable::GetDimensions(Int_t columnIndex)  const {return GetRowDescriptors()->Dimensions(columnIndex); }
02408 Int_t        TTable::GetDimensions(const Char_t *columnName) const {return GetRowDescriptors()->Dimensions(columnName); }
02409 
02410 TTable::EColumnType  TTable::GetColumnType(Int_t columnIndex)  const {return GetRowDescriptors()->ColumnType(columnIndex); }
02411 TTable::EColumnType  TTable::GetColumnType(const Char_t *columnName) const {return GetRowDescriptors()->ColumnType(columnName); }
02412 
02413 //  pointer iterator
02414 //________________________________________________________________________________________________________________
02415 TTable::piterator::piterator(const TTable *t,EColumnType type): fCurrentRowIndex(0),fCurrentColIndex(0),fRowSize(0),fCurrentRowPtr(0),fCurrentColPtr(0)
02416 {
02417    //to be documented
02418    Int_t sz = 0;
02419    if (t) sz = t->GetNRows();
02420    if (sz) {
02421       fRowSize       = t->GetRowSize();
02422       fCurrentRowPtr = (const Char_t *)t->GetArray();
02423 
02424       TTableDescriptor    *tabsDsc       = t->GetRowDescriptors();
02425       TTableDescriptor::iterator ptr     = tabsDsc->begin();
02426       TTableDescriptor::iterator lastPtr = tabsDsc->end();
02427       UInt_t i =0;
02428       for( i = 0; ptr != lastPtr; ptr++,i++)
02429          if ( tabsDsc->ColumnType(i) == type ) fPtrs.push_back(tabsDsc->Offset(i));
02430       if (fPtrs.size()==0) {
02431          MakeEnd(t->GetNRows());
02432       } else {
02433          column();
02434       }
02435    } else {
02436       MakeEnd(0);
02437    }
02438 } // piterator(TTable *)
02439 

Generated on Tue Jul 5 14:45:09 2011 for ROOT_528-00b_version by  doxygen 1.5.1