Logo ROOT   6.10/09
Reference Guide
TTable.cxx
Go to the documentation of this file.
1 // @(#)root/table:$Id$
2 // Author: Valery Fine(fine@bnl.gov) 03/07/98
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 
13 ////////////////////////////////////////////////////////////////////////////
14 // //
15 // TTable //
16 // //
17 // Wraps the array of the plain C-structures (one C-structure per element)//
18 // //
19 // class TTable provides the automatic schema evolution for //
20 // the derived "table" classes saved with ROOT format. //
21 // //
22 // "Automatic Schema evolution" provides: //
23 // - skipping data-member if it is not present for the current //
24 // implementation of the "table" but was present at the time the //
25 // table was written; //
26 // - assign a default value ZERO for the brand-new data-members, //
27 // those were not in the structure when the object was written but //
28 // present now; //
29 // - trace propely any change in the order of the data-members //
30 // //
31 // To enjoy this class one has to derive one's own custom class: //
32 // //
33 // St_dst_track_Table.h: //
34 // --------------------- //
35 // #ifndef STAF_St_dst_track_Table //
36 // #define STAF_St_dst_track_Table //
37 // //
38 // #include "TTable.h" //
39 // //
40 // // #include "dst_track.h" the C-structure defintion may be kept //
41 // separately //
42 // typedef struct dst_track_st { //
43 // float r0; /* radius at start (cm). See also comments*/
44 // float phi0; /* azimuthal angle at start (deg) */
45 // float z0; /* z-coord. at start (cm) */
46 // float psi; /* azimuthal angle of pT vector (deg) */
47 // float tanl; /* tan(dip) =pz/pt at start */
48 // float invpt; /* 1/pt at start (GeV/c)^(-1) */
49 // float curvature; /* Track curvature (1/cm) */
50 // float covar[15]; /* full covariance matrix */
51 // float chisq[2]; /* Chi-square per degree of freedom */
52 // float x_first[3]; /* coord. of first measured point (cm) */
53 // float x_last[3]; /* coord. of last measured point (cm) */
54 // float length; /* from first to last point (cm) */
55 // float impact; /* primary vertex (cm) */
56 // unsigned long map[2]; /* extrap. info. (see preceding comments)*/
57 // int id; /* Primary key (see comments) */
58 // int iflag; /* bitmask quality info. (see comments) */
59 // int det_id; /* Detector id information */
60 // int method; /* Track finding/fitting method, packed */
61 // int pid; /* Geant particle ID for assumed mass */
62 // int n_point; /* SVT, TPC, FTPC component #s are packed */
63 // int n_max_point; /* SVT, TPC, FTPC component #s are packed */
64 // int n_fit_point; /* SVT, TPC, FTPC component #s are packed */
65 // int icharge; /* Particle charge in units of |e| */
66 // int id_start_vertex; /* final fit and primary track candidates */
67 // } DST_TRACK_ST; //
68 // //
69 // class St_dst_track : public TTable //
70 // { //
71 // public: //
72 // ClassDefTable(St_dst_track,dst_track_st) //
73 // ClassDef(St_dst_track,2) //C++ wrapper for <dst_track> StAF table //
74 // }; //
75 // #endif //
76 // --------------------- //
77 // //
78 // where the CPP macro defines several convinient methods for the //
79 // "table" class (see: $ROOTSYS/table/inc/Ttypes.h for details: //
80 // //
81 // #define ClassDefTable(className,structName)
82 // protected:
83 // static TTableDescriptor *fgColDescriptors;
84 // virtual TTableDescriptor *GetDescriptorPointer() const { return fgColDescriptors;}
85 // virtual void SetDescriptorPointer(TTableDescriptor *list) const { fgColDescriptors = list;}
86 // public:
87 // typedef structName* iterator;
88 // className() : TTable(_QUOTE_(className),sizeof(structName)) {SetType(_QUOTE_(structName));}
89 // className(const char *name) : TTable(name,sizeof(structName)) {SetType(_QUOTE_(structName));}
90 // className(Int_t n) : TTable(_QUOTE_(className),n,sizeof(structName)) {SetType(_QUOTE_(structName));}
91 // className(const char *name,Int_t n) : TTable(name,n,sizeof(structName)) {SetType(_QUOTE_(structName));}
92 // structName *GetTable(Int_t i=0) const { return ((structName *)GetArray())+i;}
93 // structName &operator[](Int_t i){ assert(i>=0 && i < GetNRows()); return *GetTable(i); }
94 // const structName &operator[](Int_t i) const { assert(i>=0 && i < GetNRows()); return *((const structName *)(GetTable(i))); }
95 // structName *begin() const { return GetNRows()? GetTable(0):0;}
96 // structName *end() const {Int_t i = GetNRows(); return i? GetTable(i):0;}
97 // //
98 // The class implementation file may 2 lines and look as follows: //
99 // (for the example above): //
100 // //
101 // St_dst_track_Table.cxx: //
102 // ----------------------- //
103 // #include "St_dst_track_Table.h" //
104 // TableClassImpl(St_dst_track, dst_track_st) //
105 // ----------------------- //
106 // LinkDef.h //
107 // ----------------------- //
108 // To provide ROOT I/O for this class TWO CINT dictonary entries //
109 // should be defined with your custom LinkDef.h file //
110 // 1. First entry (as usually) for the class derived from TTable //
111 // for example: //
112 // #pragma C++ class St_dst_track //
113 // 2. Second entry for the C-structure wrapped into the class. //
114 // Since C-structuire is not derived from TObject it must be //
115 // properly defined as "foreign" ROOT class //
116 // #pragma C++ class dst_track_st+; //
117 // ----------------------- //
118 // meta-variables i$ and n$ introduced //
119 // where "i$" stands for the current row index //
120 // "n$" stands for the total number of rows //
121 // meta-variable can be used along the normal //
122 // table column names in the expressions (see for example //
123 // method TTable::Draw //
124 // //
125 ////////////////////////////////////////////////////////////////////////////
126 
127 #include <assert.h>
128 
129 #ifdef WIN32
130 # include <float.h>
131 #endif
132 
133 #include "RConfigure.h"
134 
135 //#if ROOT_VERSION_CODE >= ROOT_VERSION(3,03,5)
136 #include "Riostream.h"
137 //#include <iomanip.h>
138 
139 // #endif
140 
141 #include "TROOT.h"
142 #include "TBaseClass.h"
143 #include "TSystem.h"
144 #include "TBuffer.h"
145 #include "TMath.h"
146 #include "TClass.h"
147 #include "TBrowser.h"
148 #include "TString.h"
149 #include "TInterpreter.h"
150 #include "TDataSetIter.h"
151 #include "TTable.h"
152 #include "TTableDescriptor.h"
153 #include "TColumnView.h"
154 
155 #include "TGaxis.h"
156 #include "TH1.h"
157 #include "TH2.h"
158 #include "TProfile.h"
159 #include "TVirtualPad.h"
160 #include "TEventList.h"
161 #include "TPolyMarker.h"
162 #include "TView.h"
163 #include "TGaxis.h"
164 #include "TPolyMarker3D.h"
165 
166 #include "THLimitsFinder.h"
167 
168 #include "TTableMap.h"
169 
170 static TH1 *gCurrentTableHist = 0;
171 
172 static const char *gDtorName = "dtor";
173 static Int_t gNbins[4] = {100,100,100,100}; //Number of bins per dimension
174 static Float_t gVmin[4] = {0,0,0,0}; //Minima of varexp columns
175 static Float_t gVmax[4] = {20,20,20,20}; //Maxima of varexp columns
176 
177 const char *TTable::fgTypeName[] = {
178  "NAN", "float", "int", "long", "short", "double",
179  "unsigned int", "unsigned long","unsigned short",
180  "unsigned char", "char", "Ptr_t"
181 };
182 
183 ////////////////////////////////////////////////////////////////////////////////
184 ///
185 /// ArrayLayout - calculates the array layout recursively
186 ///
187 /// Input:
188 /// -----
189 /// dim - dimension of the targeted array
190 /// size - the max index for each dimension
191 ///
192 /// Output:
193 /// ------
194 /// layout - the "start index" for each dimension of an array
195 ///
196 
197 static void ArrayLayout(UInt_t *layout,const UInt_t *size, Int_t dim)
198 {
199  if (dim && layout && size) {
200  if (++layout[dim-1] >= size[dim-1]) {
201  layout[dim-1] = 0;
202  dim--;
203  ArrayLayout(layout,size, dim);
204  }
205  }
206 }
207 
209 
210 ////////////////////////////////////////////////////////////////////////////////
211 /// protected: create a new TTableDescriptor descriptor for this table
212 
213 TTableDescriptor *TTable::GetTableDescriptors() const {
214  assert(0);
215  return new TTableDescriptor(this);
216 }
217 
218 ////////////////////////////////////////////////////////////////////////////////
219 ///
220 /// AsString represents the value provided via "void *b" with type defined
221 /// by "name"
222 ///
223 /// void *buf - the pointer to the value to be printed out.
224 /// type - the basic data type for the value above
225 /// width - the number of psotion to be used to print the value out
226 ///
227 
228 void TTable::AsString(void *buf, EColumnType type, Int_t width,std::ostream &out) const
229 {
230  int prevPrec = out.precision();
231  const std::ios_base::fmtflags prevFmt = out.flags();
232 
233  switch (type) {
234  case kFloat:
235  out << std::dec << std::setw(width) << std::setprecision(width-3) << *(float *)buf;
236  break;
237  case kInt:
238  out << std::dec << std::setw(width) << *(int *)buf;
239  break;
240  case kLong:
241  out << std::dec << std::setw(width) << *(Long_t *)buf;
242  break;
243  case kShort:
244  out << std::dec << std::setw(width) << *(short *)buf;
245  break;
246  case kDouble:
247  out << std::dec << std::setw(width) << std::setprecision(width-3) << *(double *)buf;
248  break;
249  case kUInt:
250  out << std::dec << std::setw(width) << *(unsigned int *)buf;
251  break;
252  case kULong:
253  out << std::dec << std::setw(width) << *(ULong_t *)buf;
254  break;
255  case kUShort:
256  out << std::setw(width) << "0x" << std::hex << *(unsigned short *)buf;
257  break;
258  case kUChar:
259  out << std::setw(width) << "0x" << std::hex << int(*(unsigned char *)buf);
260  break;
261  case kChar:
262  out << std::setw(width) << *(char *)buf;
263  break;
264  case kBool:
265  out << std::setw(width) << *(Bool_t *)buf;
266  break;
267  case kPtr:
268  out << "->" << std::setw(width) << *(void **)buf;
269  break;
270  default:
271  out << "\"NaN\"";
272  break;
273  };
274  out.precision(prevPrec);
275  out.setf(prevFmt);
276 }
277 
278 ////////////////////////////////////////////////////////////////////////////////
279 ///return table type name
280 
282 {
283  return fgTypeName[type];
284 }
285 
286 ////////////////////////////////////////////////////////////////////////////////
287 /// return the Id of the C basic type by given name
288 /// return kNAN if the name provided fits no knwn basic name.
289 ///
290 
292 {
293  Int_t allTypes = sizeof(fgTypeName)/sizeof(const char *);
294  for (int i = 0; i < allTypes; i++)
295  if (!strcmp(fgTypeName[i],typeName)) return EColumnType(i);
296  return kNAN;
297 }
298 
299 ////////////////////////////////////////////////////////////////////////////////
300 /// Returns a pointer to the i-th row of the table
301 
302 const void *TTable::At(Int_t i) const
303 {
304  if (!BoundsOk("TTable::At", i)) {
305  Warning("TTable::At","%s.%s",GetName(),GetType());
306  i = 0;
307  }
308  return (const void *)(fTable+i*fSize);
309 }
310 
311 ////////////////////////////////////////////////////////////////////////////////
312 /// CopyRows copies nRows from starting from the srcRow of srcTable
313 /// to the dstRow in this table upto nRows or by the end of this table.
314 ///
315 /// This table if automaticaly increased if expand = kTRUE.
316 /// The old values of this table rows are to be destroyed and
317 /// replaced with the new ones.
318 ///
319 /// PARAMETERS:
320 /// srcTable - a pointer to the table "donor"
321 /// srcRow - the index of the first row of the table donor to copy from
322 /// dstRow - the index of the first row of this table to copy to
323 /// nRows - the total number of rows to be copied. This table will be expanded
324 /// as needed if expand = kTRUE (it is kFALSE "by default")
325 /// = 0 to copy ALL remain rows from the srcTable.
326 /// expand - flag whether this table should reallocated if needed.
327 ///
328 /// RETURN:
329 /// the number of the rows been copied
330 
331 Int_t TTable::CopyRows(const TTable *srcTable, Long_t srcRow, Long_t dstRow, Long_t nRows, Bool_t expand)
332 {
333  assert(!TestBit(kIsNotOwn));
334  if (!(srcTable && srcTable->GetNRows()) || srcRow > srcTable->GetNRows()-1 ) return 0;
335  if (strcmp(GetType(),srcTable->GetType())) {
336  // check this table current capacity
337  if (!nRows) nRows = srcTable->GetNRows();
338  Long_t tSize = GetTableSize();
339  Long_t extraRows = (tSize - dstRow) - nRows;
340  if (extraRows < 0) {
341  if (expand) {
342  ReAllocate(tSize - extraRows);
343  extraRows = 0;
344  }
345  nRows += extraRows;
346  }
347  if (dstRow+nRows > GetNRows()) SetNRows(dstRow+nRows);
348  ::memmove((*this)[dstRow],(*srcTable)[srcRow],(size_t)GetRowSize()*nRows);
349  return nRows;
350  } else
351  Error("CopyRows",
352  "This table is <%s> but the src table has a wrong type <%s>",GetType()
353  ,srcTable->GetType());
354  return 0;
355 }
356 ////////////////////////////////////////////////////////////////////////////////
357 /// Delete one or several rows from the table
358 ///
359 /// Int_t indx - index of the first row to be deleted
360 /// Int_t nRows - the total number of rows to be deleted
361 /// = 1 "by default
362 
364 {
365  if (CopyRows(this, indx+nRows, indx, GetNRows()-indx-nRows))
366  SetUsedRows(GetNRows() - nRows);
367 }
368 ////////////////////////////////////////////////////////////////////////////////
369 ///*-*-*-*-*-*-*-*-*-*-*Draw expression varexp for specified entries-*-*-*-*-*
370 ///*-* ===========================================
371 ///
372 /// This function accepts TCut objects as arguments.
373 /// Useful to use the string operator +
374 /// example:
375 /// table.Draw("x",cut1+cut2+cut3);
376 ///
377 /// TCutG object with "CUTG" name can be created via the graphics editor.
378 ///
379 
380 TH1 *TTable::Draw(TCut varexp, TCut selection, Option_t *option, Int_t nentries, Int_t firstentry)
381 {
382  return TTable::Draw(varexp.GetTitle(), selection.GetTitle(), option, nentries, firstentry);
383 }
384 
385 ////////////////////////////////////////////////////////////////////////////////
386 ///*-*-*-*-*-*-*-*-*-*-*Draw expression varexp for specified entries-*-*-*-*-*
387 ///*-* ===========================================
388 ///
389 /// varexp is an expression of the general form e1:e2:e3
390 /// where e1,etc is a C++ expression referencing a combination of the TTable columns
391 /// One can use two extra meta variable "i$" and "n$" along with the table
392 /// column names.
393 /// i$ is to involve the current row number
394 /// n$ refers the total num,ber of rows of this table provided by TTable::GetNRows()
395 ///
396 /// Example:
397 /// varexp = x simplest case: draw a 1-Dim distribution of column named x
398 /// = sqrt(x) : draw distribution of sqrt(x)
399 /// = x*y/z
400 /// = y:sqrt(x) 2-Dim dsitribution of y versus sqrt(x)
401 /// = i$:sqrt(x) 2-Dim dsitribution of i versus sqrt(x[i])
402 /// = phep[0]:sqrt(phep[3]) 2-Dim dsitribution of phep[0] versus sqrt(phep[3])
403 ///
404 /// Note that the variables e1, e2 or e3 may contain a boolean expression as well.
405 /// example, if e1= x*(y<0), the value histogrammed will be x if y<0
406 /// and will be 0 otherwise.
407 ///
408 /// selection is a C++ expression with a combination of the columns.
409 /// The value corresponding to the selection expression is used as a weight
410 /// to fill the histogram.
411 /// If the expression includes only boolean operations, the result
412 /// is 0 or 1. If the result is 0, the histogram is not filled.
413 /// In general, the expression may be of the form:
414 ///
415 /// value*(boolean expression)
416 ///
417 /// if boolean expression is true, the histogram is filled with
418 /// a weight = value.
419 /// Examples:
420 /// selection1 = "x<y && sqrt(z)>3.2 && 6 < i$ && i$ < n$"
421 /// selection2 = "(x+y)*(sqrt(z)>3.2"
422 /// selection3 = "signal*(log(signal)>1.2)"
423 /// selection1 returns a weigth = 0 or 1
424 /// selection2 returns a weight = x+y if sqrt(z)>3.2
425 /// returns a weight = 0 otherwise.
426 /// selection3 returns a weight = signal if log(signal)>1.2
427 ///
428 /// option is the drawing option
429 /// see TH1::Draw for the list of all drawing options.
430 /// If option contains the string "goff", no graphics is generated.
431 ///
432 /// nentries is the number of entries to process (default is all)
433 /// first is the first entry to process (default is 0)
434 ///
435 /// Saving the result of Draw to an histogram
436 /// =========================================
437 /// By default the temporary histogram created is called htemp.
438 /// If varexp0 contains >>hnew (following the variable(s) name(s),
439 /// the new histogram created is called hnew and it is kept in the current
440 /// directory.
441 /// Example:
442 /// tree.Draw("sqrt(x)>>hsqrt","y>0")
443 /// will draw sqrt(x) and save the histogram as "hsqrt" in the current
444 /// directory.
445 ///
446 /// By default, the specified histogram is reset.
447 /// To continue to append data to an existing histogram, use "+" in front
448 /// of the histogram name;
449 /// table.Draw("sqrt(x)>>+hsqrt","y>0")
450 /// will not reset hsqrt, but will continue filling.
451 ///
452 /// Making a Profile histogram
453 /// ==========================
454 /// In case of a 2-Dim expression, one can generate a TProfile histogram
455 /// instead of a TH2F histogram by specyfying option=prof or option=profs.
456 /// The option=prof is automatically selected in case of y:x>>pf
457 /// where pf is an existing TProfile histogram.
458 ///
459 /// Saving the result of Draw to a TEventList
460 /// =========================================
461 /// TTable::Draw can be used to fill a TEventList object (list of entry numbers)
462 /// instead of histogramming one variable.
463 /// If varexp0 has the form >>elist , a TEventList object named "elist"
464 /// is created in the current directory. elist will contain the list
465 /// of entry numbers satisfying the current selection.
466 /// Example:
467 /// tree.Draw(">>yplus","y>0")
468 /// will create a TEventList object named "yplus" in the current directory.
469 /// In an interactive session, one can type (after TTable::Draw)
470 /// yplus.Print("all")
471 /// to print the list of entry numbers in the list.
472 ///
473 /// By default, the specified entry list is reset.
474 /// To continue to append data to an existing list, use "+" in front
475 /// of the list name;
476 /// table.Draw(">>+yplus","y>0")
477 /// will not reset yplus, but will enter the selected entries at the end
478 /// of the existing list.
479 ///
480 
481 TH1 *TTable::Draw(const char *varexp00, const char *selection, Option_t *option,Int_t nentries, Int_t firstentry)
482 {
483  if (GetNRows() == 0 || varexp00 == 0 || varexp00[0]==0) return 0;
484  TString opt;
485 // char *hdefault = (char *)"htemp";
486  const char *hdefault = "htemp";
487  Int_t i,j,action;
488  Int_t hkeep = 0;
489  opt = option;
490  opt.ToLower();
491  char *varexp0 = StrDup(varexp00);
492  char *hname = strstr(varexp0,">>");
493  TH1 *oldh1 = 0;
494  TEventList *elist = 0;
495  Bool_t profile = kFALSE;
496 
497  gCurrentTableHist = 0;
498  if (hname) {
499  *hname = 0;
500  hname += 2;
501  hkeep = 1;
502  i = strcspn(varexp0,">>");
503  Bool_t hnameplus = kFALSE;
504  while (*hname == ' ') hname++;
505  if (*hname == '+') {
506  hnameplus = kTRUE;
507  hname++;
508  while (*hname == ' ') hname++;
509  j = strlen(hname)-1;
510  while (j) {
511  if (hname[j] != ' ') break;
512  hname[j] = 0;
513  j--;
514  }
515  }
516  if (i) {
517  oldh1 = (TH1*)gDirectory->Get(hname);
518  if (oldh1 && !hnameplus) oldh1->Reset();
519  } else {
520  elist = (TEventList*)gDirectory->Get(hname);
521  if (!elist) {
522  elist = new TEventList(hname,selection,1000,0);
523  }
524  if (elist && !hnameplus) elist->Reset();
525  }
526  }
527  if (!hname || *hname==0) {
528  hkeep = 0;
529  if (gDirectory) {
530  oldh1 = (TH1*)gDirectory->Get(hdefault);
531  if (oldh1 ) { oldh1->Delete(); oldh1 = 0;}
532  }
533  }
534 
535  // Look for colons
536  const Char_t *expressions[] ={varexp0,0,0,0,selection};
537  Int_t maxExpressions = sizeof(expressions)/sizeof(Char_t *);
538  Char_t *nextColon = varexp0;
539  Int_t colIndex = 1;
540  while ((nextColon = strchr(nextColon,':')) && ( colIndex < maxExpressions - 1 ) ) {
541  *nextColon = 0;
542  nextColon++;
543  expressions[colIndex] = nextColon;
544  colIndex++;
545  }
546 
547  expressions[colIndex] = selection;
548 
549 
550 ////////////////////////////////////////////////////////////////////////////////
551 
552  Printf(" Draw %s for <%s>\n", varexp00, selection);
553  Char_t *exprFileName = MakeExpression(expressions,colIndex+1);
554  if (!exprFileName) {
555  delete [] varexp0;
556  return 0;
557  }
558 
559 //--------------------------------------------------
560 // if (!fVar1 && !elist) return 0;
561 
562 //*-*- In case oldh1 exists, check dimensionality
563  Int_t dimension = colIndex;
564 
565  TString title = expressions[0];
566  for (i=1;i<colIndex;i++) {
567  title += ":";
568  title += expressions[i];
569  }
570  Int_t nsel = strlen(selection);
571  if (nsel > 1) {
572  if (nsel < 80-title.Length()) {
573  title += "{";
574  title += selection;
575  title += "}";
576  } else
577  title += "{...}";
578  }
579 
580  const Char_t *htitle = title.Data();
581 
582  if (oldh1) {
583  Int_t mustdelete = 0;
584  if (oldh1->InheritsFrom(TProfile::Class())) profile = kTRUE;
585  if (opt.Contains("prof")) {
586  if (!profile) mustdelete = 1;
587  } else {
588  if (oldh1->GetDimension() != dimension) mustdelete = 1;
589  }
590  if (mustdelete) {
591  Warning("Draw","Deleting old histogram with different dimensions");
592  delete oldh1; oldh1 = 0;
593  }
594  }
595 //*-*- Create a default canvas if none exists
596  if (!gPad && !opt.Contains("goff") && dimension > 0) {
597  gROOT->MakeDefCanvas();
598  }
599 //*-*- 1-D distribution
600  if (dimension == 1) {
601  action = 1;
602  if (!oldh1) {
603  gNbins[0] = 100;
604  if (gPad && opt.Contains("same")) {
605  TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
606  if (oldhtemp) {
607  gNbins[0] = oldhtemp->GetXaxis()->GetNbins();
608  gVmin[0] = oldhtemp->GetXaxis()->GetXmin();
609  gVmax[0] = oldhtemp->GetXaxis()->GetXmax();
610  } else {
611  gVmin[0] = gPad->GetUxmin();
612  gVmax[0] = gPad->GetUxmax();
613  }
614  } else {
615  action = -1;
616  }
617  }
618  TH1F *h1;
619  if (oldh1) {
620  h1 = (TH1F*)oldh1;
621  gNbins[0] = h1->GetXaxis()->GetNbins(); // for proofserv
622  } else {
623  h1 = new TH1F(hname,htitle,gNbins[0],gVmin[0],gVmax[0]);
624  if (!hkeep) {
625  h1->SetBit(kCanDelete);
626  h1->SetDirectory(0);
627  }
628  if (opt.Length() && opt[0] == 'e') h1->Sumw2();
629  }
630 
631  EntryLoop(exprFileName,action, h1, nentries, firstentry, option);
632 
633 // if (!fDraw && !opt.Contains("goff")) h1->Draw(option);
634  if (!opt.Contains("goff")) h1->Draw(option);
635 
636 //*-*- 2-D distribution
637  } else if (dimension == 2) {
638  action = 2;
639  if (!opt.Contains("same") && gPad) gPad->Clear();
640  if (!oldh1 || !opt.Contains("same")) {
641  gNbins[0] = 40;
642  gNbins[1] = 40;
643  if (opt.Contains("prof")) gNbins[1] = 100;
644  if (opt.Contains("same")) {
645  TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
646  if (oldhtemp) {
647  gNbins[1] = oldhtemp->GetXaxis()->GetNbins();
648  gVmin[1] = oldhtemp->GetXaxis()->GetXmin();
649  gVmax[1] = oldhtemp->GetXaxis()->GetXmax();
650  gNbins[0] = oldhtemp->GetYaxis()->GetNbins();
651  gVmin[0] = oldhtemp->GetYaxis()->GetXmin();
652  gVmax[0] = oldhtemp->GetYaxis()->GetXmax();
653  } else {
654  gNbins[1] = 40;
655  gVmin[1] = gPad->GetUxmin();
656  gVmax[1] = gPad->GetUxmax();
657  gNbins[0] = 40;
658  gVmin[0] = gPad->GetUymin();
659  gVmax[0] = gPad->GetUymax();
660  }
661  } else {
662  action = -2;
663  }
664  }
665  if (profile || opt.Contains("prof")) {
666  TProfile *hp;
667  if (oldh1) {
668  action = 4;
669  hp = (TProfile*)oldh1;
670  } else {
671  if (action < 0) action = -4;
672  if (opt.Contains("profs"))
673  hp = new TProfile(hname,htitle,gNbins[1],gVmin[1], gVmax[1],"s");
674  else
675  hp = new TProfile(hname,htitle,gNbins[1],gVmin[1], gVmax[1],"");
676  if (!hkeep) {
677  hp->SetBit(kCanDelete);
678  hp->SetDirectory(0);
679  }
680  }
681 
682  EntryLoop(exprFileName,action,hp,nentries, firstentry, option);
683 
684  if (!opt.Contains("goff")) hp->Draw(option);
685  } else {
686  TH2F *h2;
687  if (oldh1) {
688  h2 = (TH2F*)oldh1;
689  } else {
690  h2 = new TH2F(hname,htitle,gNbins[1],gVmin[1], gVmax[1], gNbins[0], gVmin[0], gVmax[0]);
691  if (!hkeep) {
692  const Int_t kNoStats = BIT(9);
693  h2->SetBit(kCanDelete);
694  h2->SetBit(kNoStats);
695  h2->SetDirectory(0);
696  }
697  }
698  Int_t noscat = strlen(option);
699  if (opt.Contains("same")) noscat -= 4;
700  if (noscat) {
701  EntryLoop(exprFileName,action,h2,nentries, firstentry, option);
702  // if (!fDraw && !opt.Contains("goff")) h2->Draw(option);
703  if (!opt.Contains("goff")) h2->Draw(option);
704  } else {
705  action = 12;
706  if (!oldh1 && !opt.Contains("same")) action = -12;
707  EntryLoop(exprFileName,action,h2,nentries, firstentry, option);
708 // if (oldh1 && !fDraw && !opt.Contains("goff")) h2->Draw(option);
709  if (oldh1 && !opt.Contains("goff")) h2->Draw(option);
710  }
711  }
712 
713 //*-*- 3-D distribution
714  } else if (dimension == 3) {
715  action = 13;
716  if (!opt.Contains("same")) action = -13;
717  EntryLoop(exprFileName,action,0,nentries, firstentry, option);
718 
719 //*-* an Event List
720  //} else if (elist) {
721  // action = 5;
722 // Int_t oldEstimate = fEstimate;
723 // SetEstimate(1);
724  // EntryLoop(exprFileName,action,elist,nentries, firstentry, option);
725 // SetEstimate(oldEstimate);
726  }
727  delete [] exprFileName;
728  delete [] varexp0;
729  return gCurrentTableHist;
730 }
731 
732 ////////////////////////////////////////////////////////////////////////////////
733 ///*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
734 ///*-* ==========================
735 ///*-* This mathod is a straight copy of void TTree::FindGoodLimits method
736 ///*-*
737 
739 {
740  Double_t binlow=0,binhigh=0,binwidth=0;
741  Int_t n;
742  Double_t dx = 0.1*(xmax-xmin);
743  Double_t umin = xmin - dx;
744  Double_t umax = xmax + dx;
745  if (umin < 0 && xmin >= 0) umin = 0;
746  if (umax > 0 && xmax <= 0) umax = 0;
747 
748 #if ROOT_VERSION_CODE >= ROOT_VERSION(3,03,5)
749  THLimitsFinder::Optimize(umin,umax,nbins,binlow,binhigh,n,binwidth,"");
750 #else
751  TGaxis::Optimize(umin,umax,nbins,binlow,binhigh,n,binwidth,"");
752 #endif
753 
754  if (binwidth <= 0 || binwidth > 1.e+39) {
755  xmin = -1;
756  xmax = 1;
757  } else {
758  xmin = binlow;
759  xmax = binhigh;
760  }
761 
762  newbins = nbins;
763 }
764 
765 ////////////////////////////////////////////////////////////////////////////////
766 ///
767 /// EntryLoop creates a CINT bytecode to evaluate the given expressions for
768 /// all table rows in loop and fill the appropriated histograms.
769 ///
770 /// Solution for Byte code
771 /// From: Masaharu Goto <MXJ02154@nifty.ne.jp>
772 /// To: <fine@bnl.gov>
773 /// Cc: <rootdev@hpsalo.cern.ch>
774 /// Sent: 13-th august 1999 year 23:01
775 ///
776 /// action = 1 Fill 1-D histogram obj
777 /// = 2 Fill 2-D histogram obj
778 /// = 3 Fill 3-D histogram obj
779 /// = 4 Fill Profile histogram obj
780 /// = 5 Fill a TEventlist
781 /// = 11 Estimate Limits
782 /// = 12 Fill 2-D PolyMarker obj
783 /// = 13 Fill 3-D PolyMarker obj
784 /// action < 0 Evaluate Limits for case abs(action)
785 ///
786 /// Load file
787 
788 Bool_t TTable::EntryLoop(const Char_t *exprFileName,Int_t &action, TObject *obj
789  ,Int_t nentries, Int_t firstentry, Option_t *option)
790 {
791  Double_t rmin[3],rmax[3];
792  if (gInterpreter->LoadFile((Char_t *)exprFileName) < 0) {
793  Error("EntryLoop","Error: loading file %s",exprFileName);
794  return kFALSE; // can not load file
795  }
796 
797  // Float_t Selection(Float_t *results[], void *address[], int& i$, int n$)
798  // where i$ - meta variable to set current row index
799  // n$ - meta variable to set the total table size
800  const Char_t *funcName = "SelectionQWERTY";
801 #ifdef BYTECODE
802  const Char_t *argtypes = "Float_t *,float **, int&, int& ";
803  Long_t offset;
804  G__ClassInfo globals;
805  G__MethodInfo func = globals.GetMethod(funcName,argtypes,&offset);
806 
807  // Compile bytecode
808  struct G__bytecodefunc *pbc = func.GetBytecode();
809  if(!pbc) {
810  Error("EntryLoop","Bytecode compilation %s",funcName);
811  G__unloadfile((Char_t *)exprFileName);
812  return kFALSE; // can not get bytecode
813  }
814 #endif
815  // Prepare callfunc object
816  int i;
817  int nRows = GetNRows();
818  TTableDescriptor *tabsDsc = GetRowDescriptors();
819  tableDescriptor_st *descTable = tabsDsc->GetTable();
820  Float_t results[] = {1,1,1,1,1};
821  Char_t **addressArray = (Char_t **)new ULong_t[tabsDsc->GetNRows()];
822  Char_t *thisTable = (Char_t *)GetArray();
823  const Char_t *argtypes = "Float_t *,float **, int&, int& ";
824  Long_t offset;
825  ClassInfo_t *globals = gInterpreter->ClassInfo_Factory();
826  CallFunc_t *callfunc = gInterpreter->CallFunc_Factory();
827  gInterpreter->CallFunc_SetFunc(callfunc,globals,funcName,argtypes,&offset);
828 
829  gInterpreter->CallFunc_SetArg(callfunc,(Long_t)(&results[0])); // give 'Float_t *results[5]' as 1st argument
830  gInterpreter->CallFunc_SetArg(callfunc,(Long_t)(addressArray)); // give 'void *addressArray[]' as 2nd argument
831  gInterpreter->CallFunc_SetArg(callfunc,(Long_t)(&i)); // give 'int& i$' as 3nd argument
832  gInterpreter->CallFunc_SetArg(callfunc,(Long_t)(&nRows)); // give 'int& n$= nRows as 4th argument
833 
834  // Call bytecode in loop
835 
836 #define CALLMETHOD gInterpreter->CallFunc_Exec(callfunc,0);
837 
838 #define TAKEACTION_BEGIN \
839  descTable = tabsDsc->GetTable(); \
840  for (i=0; i < tabsDsc->GetNRows(); i++,descTable++ ) \
841  addressArray[i] = addressEntry + descTable->fOffset; \
842  for(i=firstentry;i<lastEntry;i++) { \
843  CALLMETHOD
844 
845 #define TAKEACTION_END for (int j=0; j < tabsDsc->GetNRows(); j++ ) addressArray[j] += rSize;}
846 
847 
848  if (firstentry < nRows ) {
849  Long_t rSize = GetRowSize();
850  Char_t *addressEntry = thisTable + rSize*firstentry;
851  Int_t lastEntry = TMath::Min(UInt_t(firstentry+nentries),UInt_t(nRows));
852  if (action < 0) {
853  gVmin[0] = gVmin[1] = gVmin[2] = 1e30;
854  gVmax[0] = gVmax[1] = gVmax[2] = -gVmin[0];
855  }
856  Int_t nchans = 0;
857  switch ( action ) {
858  case -1: {
860  if (results[1]) {
861  if (gVmin[0] > results[0]) gVmin[0] = results[0];
862  if (gVmax[0] < results[0]) gVmax[0] = results[0];
863  }
865 
866  nchans = gNbins[0];
867  if (gVmin[0] >= gVmax[0]) { gVmin[0] -= 1; gVmax[0] += 1;}
868  FindGoodLimits(nchans,gNbins[0],gVmin[0],gVmax[0]);
869  ((TH1 *)obj)->SetBins(gNbins[0],gVmin[0],gVmax[0]);
870  }
871  // Intentional fall though
872  case 1:
874  if (results[1]) ((TH1 *)obj)->Fill(Axis_t(results[0]),Stat_t(results[1]));
876  gCurrentTableHist = ((TH1 *)obj);
877  break;
878  case -2:
880  if (results[2]) {
881  if (gVmin[0] > results[1]) gVmin[0] = results[1];
882  if (gVmax[0] < results[1]) gVmax[0] = results[1];
883  if (gVmin[1] > results[0]) gVmin[1] = results[0];
884  if (gVmax[1] < results[0]) gVmax[1] = results[0];
885  }
887  nchans = gNbins[0];
888  if (gVmin[0] >= gVmax[0]) { gVmin[0] -= 1; gVmax[0] += 1;}
889  FindGoodLimits(nchans,gNbins[0],gVmin[0],gVmax[0]);
890  if (gVmin[1] >= gVmax[1]) { gVmin[1] -= 1; gVmax[1] += 1;}
891  FindGoodLimits(nchans,gNbins[1],gVmin[1],gVmax[1]);
892  ((TH1*)obj)->SetBins(gNbins[1],gVmin[1],gVmax[1],gNbins[0],gVmin[0],gVmax[0]);
893  // Intentional fall though
894  case 2:
895  if (obj->IsA() == TH2F::Class()) {
897  if (results[2]) ((TH2F*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
899  } else if (obj->IsA() == TH2S::Class()) {
901  if (results[2]) ((TH2S*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
903  } else if (obj->IsA() == TH2C::Class()) {
905  if (results[2]) ((TH2C*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
907  } else if (obj->IsA() == TH2D::Class()) {
909  if (results[2]) ((TH2D*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
911  }
912  gCurrentTableHist = ((TH1 *)obj);
913  break;
914  case -4:
916  if (results[2]) {
917  if (gVmin[0] > results[1]) gVmin[0] = results[1];
918  if (gVmax[0] < results[1]) gVmax[0] = results[1];
919  if (gVmin[1] > results[0]) gVmin[1] = results[0];
920  if (gVmax[1] < results[0]) gVmax[1] = results[0];
921  }
923  nchans = gNbins[1];
924  if (gVmin[1] >= gVmax[1]) { gVmin[1] -= 1; gVmax[1] += 1;}
925  FindGoodLimits(nchans,gNbins[1],gVmin[1],gVmax[1]);
926  ((TProfile*)obj)->SetBins(gNbins[1],gVmin[1],gVmax[1]);
927  // Intentional fall though
928  case 4:
930  if (results[2]) ((TProfile*)obj)->Fill(Axis_t(results[0]),Axis_t(results[1]),Stat_t(results[2]));
932  break;
933  case -12:
935  if (results[2]) {
936  if (gVmin[0] > results[1]) gVmin[0] = results[1];
937  if (gVmax[0] < results[1]) gVmax[0] = results[1];
938  if (gVmin[1] > results[0]) gVmin[1] = results[0];
939  if (gVmax[1] < results[0]) gVmax[1] = results[0];
940  }
942  nchans = gNbins[0];
943  if (gVmin[0] >= gVmax[0]) { gVmin[0] -= 1; gVmax[0] += 1;}
944  FindGoodLimits(nchans,gNbins[0],gVmin[0],gVmax[0]);
945  if (gVmin[1] >= gVmax[1]) { gVmin[1] -= 1; gVmax[1] += 1;}
946  FindGoodLimits(nchans,gNbins[1],gVmin[1],gVmax[1]);
947  ((TH2F*)obj)->SetBins(gNbins[1],gVmin[1],gVmax[1],gNbins[0],gVmin[0],gVmax[0]);
948  // Intentional fall though
949  case 12: {
950  if (!strstr(option,"same") && !strstr(option,"goff")) {
951  ((TH2F*)obj)->DrawCopy(option);
952  gPad->Update();
953  }
954 // pm->SetMarkerStyle(GetMarkerStyle());
955 // pm->SetMarkerColor(GetMarkerColor());
956 // pm->SetMarkerSize(GetMarkerSize());
957  Float_t *x = new Float_t[lastEntry-firstentry]; // pm->GetX();
958  Float_t *y = new Float_t[lastEntry-firstentry]; // pm->GetY();
959  Float_t u, v;
960  Float_t umin = gPad->GetUxmin();
961  Float_t umax = gPad->GetUxmax();
962  Float_t vmin = gPad->GetUymin();
963  Float_t vmax = gPad->GetUymax();
964  Int_t pointIndex = 0;
966  if (results[2]) {
967  u = gPad->XtoPad(results[0]);
968  v = gPad->YtoPad(results[1]);
969  if (u < umin) u = umin;
970  if (u > umax) u = umax;
971  if (v < vmin) v = vmin;
972  if (v > vmax) v = vmax;
973  x[pointIndex] = u;
974  y[pointIndex] = v;
975  pointIndex++;
976  }
978  if (pointIndex && !strstr(option,"goff")) {
979  TPolyMarker *pm = new TPolyMarker(pointIndex,x,y);
980  pm->Draw();
981  pm->SetBit(kCanDelete);
982  }
983  if (!((TH2F*)obj)->TestBit(kCanDelete))
984  if (pointIndex)
985  for(i=0;i<pointIndex;i++) ((TH2F*)obj)->Fill(x[i], y[i]);
986  delete [] x; delete [] y;
987  gCurrentTableHist = ((TH1*)obj);
988  }
989  break;
990  case -13:
992  if (results[3]) {
993  if (gVmin[0] > results[2]) gVmin[0] = results[2];
994  if (gVmax[0] < results[2]) gVmax[0] = results[2];
995  if (gVmin[1] > results[1]) gVmin[1] = results[1];
996  if (gVmax[1] < results[1]) gVmax[1] = results[1];
997  if (gVmin[2] > results[0]) gVmin[2] = results[0];
998  if (gVmax[2] < results[0]) gVmax[2] = results[0];
999  }
1001  rmin[0] = gVmin[2]; rmin[1] = gVmin[1]; rmin[2] = gVmin[0];
1002  rmax[0] = gVmax[2]; rmax[1] = gVmax[1]; rmax[2] = gVmax[0];
1003  gPad->Clear();
1004  gPad->Range(-1,-1,1,1);
1005  TView::CreateView(1,rmin,rmax);
1006  // Intentional fall though
1007  case 13: {
1008  TPolyMarker3D *pm3d = new TPolyMarker3D(lastEntry-firstentry);
1009  pm3d->SetBit(kCanDelete);
1010 // pm3d->SetMarkerStyle(GetMarkerStyle());
1011 // pm3d->SetMarkerColor(GetMarkerColor());
1012 // pm3d->SetMarkerSize(GetMarkerSize());
1014  if (results[3]) pm3d->SetNextPoint(results[0],results[1],results[2]);
1016  pm3d->Draw();
1017  }
1018  break;
1019  default:
1020  Error("EntryLoop","unknown action \"%d\" for table <%s>", action, GetName());
1021  break;
1022  };
1023  }
1024  gInterpreter->CallFunc_Delete(callfunc);
1025  gInterpreter->ClassInfo_Delete(globals);
1026  gInterpreter->UnloadFile((Char_t *)exprFileName);
1027  delete [] addressArray;
1028  return kTRUE;
1029 }
1030 
1031 ////////////////////////////////////////////////////////////////////////////////
1032 /// Default TTable ctor.
1033 
1034 TTable::TTable(const char *name, Int_t size) : TDataSet(name),
1035  fSize(size),fN(0), fTable(0),fMaxIndex(0)
1036 {
1037  if (size == 0) Warning("TTable(0)","Wrong table format");
1038 }
1039 
1040 ////////////////////////////////////////////////////////////////////////////////
1041 /// Create TTable object and set array size to n longs.
1042 
1043 TTable::TTable(const char *name, Int_t n,Int_t size) : TDataSet(name),
1044  fSize(size),fN(0),fTable(0),fMaxIndex(0)
1045 {
1046  if (n > 0) Set(n);
1047 }
1048 
1049 ////////////////////////////////////////////////////////////////////////////////
1050 /// Create TTable object and initialize it with values of array.
1051 
1052 TTable::TTable(const char *name, Int_t n, Char_t *table,Int_t size) : TDataSet(name),
1053  fSize(size),fN(0),fTable(0),fMaxIndex(0)
1054 {
1055  Set(n, table);
1056 }
1057 
1058 ////////////////////////////////////////////////////////////////////////////////
1059 /// Create TTable object and initialize it with values of array.
1060 
1061 TTable::TTable(const char *name, const char *type, Int_t n, Char_t *array, Int_t size)
1062  : TDataSet(name),fSize(size),fTable(0),fMaxIndex(0)
1063 {
1064  fTable = array;
1065  SetType(type);
1066  SetfN(n);
1067 }
1068 
1069 ////////////////////////////////////////////////////////////////////////////////
1070 /// Copy constructor.
1071 
1072 TTable::TTable(const TTable &table):TDataSet(table)
1073 {
1074  fTable = 0;
1075  SetUsedRows(table.GetNRows());
1076  fSize = table.GetRowSize();
1077  SetfN(table.fN);
1078  Set(table.fN, table.fTable);
1079 }
1080 
1081 ////////////////////////////////////////////////////////////////////////////////
1082 /// TTable assignment operator.
1083 /// This operator REALLOCATEs this table to fit the number of
1084 /// the USED rows of the source table if any
1085 
1087 {
1088  if (strcmp(GetType(),rhs.GetType()) == 0) {
1089  if (this != &rhs && rhs.GetNRows() >0 ){
1090  Set(rhs.GetNRows(), rhs.fTable);
1091  SetUsedRows(rhs.GetNRows());
1092  }
1093  } else
1094  Error("operator=","Can not copy <%s> table into <%s> table", rhs.GetType(),GetType());
1095  return *this;
1096 }
1097 
1098 ////////////////////////////////////////////////////////////////////////////////
1099 /// Delete TTable object.
1100 
1102 {
1103  Delete();
1104 }
1105 
1106 ////////////////////////////////////////////////////////////////////////////////
1107 /// Adopt array arr into TTable, i.e. don't copy arr but use it directly
1108 /// in TTable. User may not delete arr, TTable dtor will do it.
1109 
1110 void TTable::Adopt(Int_t n, void *arr)
1111 {
1112  Clear();
1113 
1114  SetfN(n); SetUsedRows(n);
1115  fTable = (char *)arr;
1116 }
1117 
1118 ////////////////////////////////////////////////////////////////////////////////
1119 /// Add the "row" at the GetNRows() position, and
1120 /// reallocate the table if neccesary, and
1121 /// return the row index the "row" has occupied.
1122 ///
1123 /// row == 0 see method TTable::AddAt(const void *row, Int_t i)
1124 
1125 Int_t TTable::AddAt(const void *row)
1126 {
1127  Int_t gap = GetTableSize() - GetNRows();
1128  // do we need to add an extra space?
1129  if (gap < 1) ReAllocate(GetTableSize() + TMath::Max(1,Int_t(0.3*GetTableSize())));
1130  Int_t indx = GetNRows();
1131  AddAt(row,indx);
1132  return indx;
1133 }
1134 ////////////////////////////////////////////////////////////////////////////////
1135 /// Add one element ("row") of structure at position "i".
1136 /// Check for out of bounds.
1137 ///
1138 /// If the row == 0 the "i" cell is still occupied and
1139 /// filled with the pattern "ff"
1140 
1141 void TTable::AddAt(const void *row, Int_t i)
1142 {
1143  if (!BoundsOk("TTable::AddAt", i))
1144  i = 0;
1145  if (row) memcpy(fTable+i*fSize,row,fSize);
1146  else memset(fTable+i*fSize,127,fSize);
1148 }
1149 
1150 ////////////////////////////////////////////////////////////////////////////////
1151 /// Copy the C-structure src into the new location
1152 /// the length of the strucutre is defined by this class descriptor
1153 
1155 {
1156  ::memcpy(dest,src,fSize*fN);
1157 }
1158 ////////////////////////////////////////////////////////////////////////////////
1159 ///to be documented
1160 
1162 {
1163  array.Set(fN);
1164  CopyStruct(array.fTable,fTable);
1165 }
1166 ////////////////////////////////////////////////////////////////////////////////
1167 /// Get a comment from the table descriptor
1168 
1169 const Char_t *TTable::GetColumnComment(Int_t columnIndex) const {
1170  TDataSetIter nextComment(GetRowDescriptors()->MakeCommentField(kFALSE));
1171  TDataSet *nxc = 0;
1172  for (int i=0; i<= columnIndex; i++) nxc = nextComment();
1173  return nxc ? nxc->GetTitle() : 0;
1174 }
1175 ////////////////////////////////////////////////////////////////////////////////
1176 /// Append nRows row of the array "row" to the table
1177 /// return
1178 /// - the new table size (# of table rows)
1179 /// - 0 if the object doesn't own the internal array and can not expand it
1180 
1181 Long_t TTable::AppendRows(const void *row, UInt_t nRows)
1182 {
1183  if (!TestBit(kIsNotOwn) && row && nRows ) {
1184  Int_t indx = GetNRows();
1185  ReAllocate(nRows);
1186  // Copy (insert) the extra staff in
1187  ::memmove(fTable+indx*fSize,row,fSize*nRows);
1188  }
1189  return TestBit(kIsNotOwn) ? 0 : GetSize();
1190 }
1191 ////////////////////////////////////////////////////////////////////////////////
1192 /// void InsertRows(cons void *row, Long_t indx, UInt_t nRows)
1193 ///
1194 /// Insert one or several rows into the table at "indx" position
1195 /// The rest table stuff is shifted down
1196 ///
1197 /// cons void - a pointer to the array of rows to be inserted
1198 /// Long_t indx = The position these rows will be inserted to
1199 /// Int_t nRows - the total number of rows to be inserted
1200 /// = 1 "by default
1201 /// return:
1202 /// The number of the rows has been shifted to accomodate
1203 /// the new rows.
1204 ///
1205 
1206 Long_t TTable::InsertRows(const void *row, Long_t indx, UInt_t nRows)
1207 {
1208  Long_t nShifted = 0;
1209  if (nRows > 0) {
1210  // Shift the table down
1211  nShifted = CopyRows(this, indx, indx+nRows, GetNRows()+nRows);
1212  // Copy (insert) the extra staff in
1213  ::memmove(fTable+indx*fSize,row,fSize*nRows);
1214  }
1215  return nShifted;
1216 
1217 }
1218 ////////////////////////////////////////////////////////////////////////////////
1219 /// Reallocate this table leaving only (used rows)+1 allocated
1220 /// GetTableSize() = GetNRows() + 1
1221 /// returns a pointer to the first row of the reallocated table
1222 /// Note:
1223 /// The table is reallocated if it is an owner of the internal array
1224 
1226 {
1227  ReAlloc(GetNRows()+1);
1228  return (void *)fTable;
1229 }
1230 ////////////////////////////////////////////////////////////////////////////////
1231 /// Reallocate this table leaving only <newsize> allocated
1232 /// GetTableSize() = newsize;
1233 /// returns a pointer to the first row of the reallocated table
1234 /// Note:
1235 /// The table is reallocated if it is an owner of the internal array
1236 
1238 {
1239  if (newsize > fN) ReAlloc(newsize);
1240  return (void *)fTable;
1241 }
1242 
1243 ////////////////////////////////////////////////////////////////////////////////
1244 /// The table is reallocated if it is an owner of the internal array
1245 
1246 void TTable::ReAlloc(Int_t newsize)
1247 {
1248  if (!TestBit(kIsNotOwn) && newsize > 0) {
1249  void *arr = 0;
1250  Int_t sleepCounter = 0;
1251  while (!(arr = realloc(fTable,fSize*newsize))) {
1252  sleepCounter++;
1253  Warning("ReAlloc",
1254  "Not enough memory to Reallocate %d bytes for table <%s::%s>. Please cancel some jobs",
1255  newsize, GetType(),GetName());
1256  gSystem->Sleep(1000*600);
1257  if (sleepCounter > 30) {
1258  Error("ReAlloc","I can not wait anymore. Good bye");
1259  assert(0);
1260  }
1261  }
1262  SetfN(newsize);
1263  fTable = (char *)arr;
1264  }
1265 }
1266 
1267 ////////////////////////////////////////////////////////////////////////////////
1268 /// Allocate a space for the new table, if any
1269 /// Sleep for a while if space is not available and try again
1270 
1272 {
1273  if (!fTable) {
1274  void *ptr = 0;
1275  Int_t sleepCounter = 0;
1276  while (!(ptr = malloc(fSize*fN))) {
1277  sleepCounter++;
1278  Warning("Create",
1279  "Not enough memory to allocate %d rows for table <%s::%s>. Please cancel some jobs",
1280  fN, GetType(),GetName());
1281  gSystem->Sleep(1000*600);
1282  if (sleepCounter > 30){
1283  Error("Create","I can not wait anymore. Good bye");
1284  assert(0);
1285  }
1286  }
1287  fTable = (Char_t *)ptr;
1288  // make sure all link-columns are zero
1289  memset(fTable,0,fSize*fN);
1290  }
1291  return fTable;
1292 }
1293 
1294 ////////////////////////////////////////////////////////////////////////////////
1295 /// Wrap each table coulumn with TColumnView object to browse.
1296 
1298  if (!b) return;
1299  TDataSet::Browse(b);
1300  Int_t nrows = TMath::Min(Int_t(GetNRows()),6);
1301  if (nrows == 0) nrows = 1;
1302  Print(0,nrows);
1303  // Add the table columns to the browser
1304  UInt_t nCol = GetNumberOfColumns();
1305  for (UInt_t i = 0;i<nCol;i++){
1306  TColumnView *view = 0;
1307  UInt_t nDim = GetDimensions(i);
1308  const Char_t *colName = GetColumnName(i);
1309  if (!nDim) { // scalar
1310  // This will cause a small memory leak
1311  // unless TBrowser recognizes kCanDelete bit
1312  if( GetColumnType(i)== kPtr) {
1313  UInt_t offset = GetOffset(i);
1314  TTableMap *m = *(TTableMap **)(((char *)GetArray())+offset);
1315  if (m) {
1316  TString nameMap = "*";
1317  nameMap += m->Table()->GetName();
1318  b->Add(m,nameMap.Data());
1319  }
1320  } else {
1321  view = new TColumnView(GetColumnName(i),this);
1322  view->SetBit(kCanDelete);
1323  b->Add(view,view->GetName());
1324  }
1325  } else { // array
1326  const UInt_t *indx = GetIndexArray(i);
1327  UInt_t totalSize = 1;
1328  UInt_t k;
1329  for (k=0;k<nDim; k++) totalSize *= indx[k];
1330  for (k=0;k<totalSize;k++) {
1331  TString buffer;
1332  buffer.Form("%s[%d]",colName,k);
1333  view = new TColumnView(buffer,this);
1334  view->SetBit(kCanDelete);
1335  b->Add(view,view->GetName());
1336  }
1337  }
1338  }
1339 }
1340 
1341 ////////////////////////////////////////////////////////////////////////////////
1342 /// Deletes the internal array of this class
1343 /// if this object does own its internal table
1344 
1346 {
1347  if (!fTable) return;
1348  Bool_t dtor = kFALSE;
1349  dtor = opt && (strcmp(opt,gDtorName)==0);
1350  if (!opt || !opt[0] || dtor ) {
1351  if (! TestBit(kIsNotOwn)) {
1352  if (!dtor) ResetMap();
1353  free(fTable);
1354  }
1355  fTable = 0;
1356  fMaxIndex = 0;
1357  SetfN(0);
1358  return;
1359  }
1360 }
1361 
1362 ////////////////////////////////////////////////////////////////////////////////
1363 ///
1364 /// Delete the internal array and free the memory it occupied
1365 /// if this object did own this array
1366 ///
1367 /// Then perform TDataSet::Delete(opt)
1368 
1370 {
1371  Clear(gDtorName);
1372  TDataSet::Delete(opt);
1373 }
1374 
1375 ////////////////////////////////////////////////////////////////////////////////
1376 ///to be documented
1377 
1379 {
1380  TClass *cl = 0;
1382  if (dsc) cl = dsc->RowClass();
1383  else Error("GetRowClass()","Table descriptor of <%s::%s> table lost",
1384  GetName(),GetType());
1385  return cl;
1386 }
1387 
1388 ////////////////////////////////////////////////////////////////////////////////
1389 /// Returns the number of the used rows for the wrapped table
1390 
1392  return fMaxIndex;
1393 }
1394 
1395 ////////////////////////////////////////////////////////////////////////////////
1396 /// Returns the size (in bytes) of one table row
1397 
1399  return fSize;
1400 }
1401 
1402 ////////////////////////////////////////////////////////////////////////////////
1403 /// Returns the number of the allocated rows
1404 
1406  return fN;
1407 }
1408 
1409 ////////////////////////////////////////////////////////////////////////////////
1410 ///*-*-*-*-*-*-*-*-*Fit a projected item(s) from a TTable*-*-*-*-*-*-*-*-*-*
1411 ///*-* =======================================
1412 ///
1413 /// formula is a TF1 expression.
1414 ///
1415 /// See TTable::Draw for explanations of the other parameters.
1416 ///
1417 /// By default the temporary histogram created is called htemp.
1418 /// If varexp contains >>hnew , the new histogram created is called hnew
1419 /// and it is kept in the current directory.
1420 /// Example:
1421 /// table.Fit(pol4,"sqrt(x)>>hsqrt","y>0")
1422 /// will fit sqrt(x) and save the histogram as "hsqrt" in the current
1423 /// directory.
1424 ///
1425 
1426 void TTable::Fit(const char *formula ,const char *varexp, const char *selection,Option_t *option ,Option_t *goption,Int_t nentries, Int_t firstentry)
1427 {
1428  TString opt(option);
1429  opt += "goff";
1430 
1431  Draw(varexp,selection,opt,nentries,firstentry);
1432 
1433  TH1 *hfit = gCurrentTableHist;
1434  if (hfit) {
1435  Printf("hname=%s, formula=%s, option=%s, goption=%s\n",hfit->GetName(),formula,option,goption);
1436  // remove bit temporary
1437  Bool_t canDeleteBit = hfit->TestBit(kCanDelete);
1438  if (canDeleteBit) hfit->ResetBit(kCanDelete);
1439  hfit->Fit(formula,option,goption);
1440  if (TestBit(canDeleteBit)) hfit->SetBit(kCanDelete);
1441  }
1442  else Printf("ERROR hfit=0\n");
1443 }
1444 
1445 ////////////////////////////////////////////////////////////////////////////////
1446 ///Returns the type of the wrapped C-structure kept as the TNamed title
1447 
1448 const Char_t *TTable::GetType() const
1449 {
1450  return GetTitle();
1451 }
1452 
1453 ////////////////////////////////////////////////////////////////////////////////
1454 /// return Folder flag to be used by TBrowse object
1455 /// The table is a folder if
1456 /// - it has sub-dataset
1457 /// or
1458 /// - GetNRows > 0
1459 
1461  return kTRUE; // to provide the "fake" folder bit to workaround TKey::Browse()
1462 
1463 #if 0
1464  // this became useless due TKey::Browse new implementation
1465  return
1466  (fList && fList->Last() ? kTRUE : kFALSE)
1467  ||
1468  (GetNRows() > 0);
1469 #endif
1470 }
1471 
1472 ////////////////////////////////////////////////////////////////////////////////
1473 ///
1474 /// return the total number of the NaN for float/double cells of this table
1475 /// Thanks Victor Perevoztchikov
1476 ///
1477 
1479 {
1480  EColumnType code;
1481  char const *cell,*colname,*table;
1482  double word;
1483  int icol,irow,colsize,wordsize,nwords,iword,nerr,offset;
1484 
1485  TTableDescriptor *rowDes = GetRowDescriptors();
1486  assert(rowDes!=0);
1487  table = (const char*)GetArray();
1488 
1489  int ncols = rowDes->GetNumberOfColumns();
1490 
1491  int lrow = GetRowSize();
1492  int nrows = GetNRows ();
1493  nerr =0;
1494  for (icol=0; icol < ncols; icol++) {// loop over cols
1495  code = rowDes->GetColumnType(icol);
1496  if (code!=kFloat && code!=kDouble) continue;
1497 
1498  offset = rowDes->GetOffset (icol);
1499  colsize = rowDes->GetColumnSize(icol);
1500  wordsize = rowDes->GetTypeSize (icol);
1501  nwords = colsize/wordsize;
1502  for (irow=0; irow < nrows; irow++) { //loop over rows
1503  cell = table + offset + irow*lrow;
1504  for (iword=0;iword<nwords; iword++,cell+=wordsize) { //words in col
1505  word = (code==kDouble) ? *(double*)cell : *(float*)cell;
1506  if (TMath::Finite(word)) continue;
1507 // ERROR FOUND
1508  nerr++; colname = rowDes->GetColumnName(icol);
1509  Warning("NaN"," Table %s.%s.%d\n",GetName(),colname,irow);
1510  }
1511  }
1512  }
1513  return nerr;
1514 }
1515 
1516 ////////////////////////////////////////////////////////////////////////////////
1517 /// This static method creates a new TTable object if provided
1518 
1519 TTable *TTable::New(const Char_t *name, const Char_t *type, void *array, UInt_t size)
1520 {
1521  TTable *table = 0;
1522  if (type && name) {
1523  TString tableType(type);
1524  TString t = tableType.Strip();
1525 
1526  TString classname("St_");
1527  classname += t;
1528  TClass *cl = TClass::GetClass(classname);
1529  if (cl) {
1530  table = (TTable *)cl->New();
1531  if (table) {
1532  table->SetTablePointer(array);
1533  table->SetName(name);
1534  table->SetfN(size);
1535  table->SetUsedRows(size);
1536  }
1537  }
1538  }
1539  return table;
1540 }
1541 ////////////////////////////////////////////////////////////////////////////////
1542 /// Generate an out-of-bounds error. Always returns false.
1543 
1544 Bool_t TTable::OutOfBoundsError(const char *where, Int_t i) const
1545 {
1546  Error(where, "index %d out of bounds (size: %d, this: 0x%lx)", i, fN, (Long_t)this);
1547  return kFALSE;
1548 }
1549 ////////////////////////////////////////////////////////////////////////////////
1550 /// Create IDL table defintion (to be used for XDF I/O)
1551 
1552 Char_t *TTable::Print(Char_t *strbuf,Int_t lenbuf) const
1553 {
1554  Int_t iOut = 0;
1555 
1557  if (!dscT ) {
1558  Error("Print"," No dictionary entry for <%s> structure", GetTitle());
1559  if (lenbuf>0) iOut += snprintf(strbuf,lenbuf," *** Errror ***");
1560  return strbuf;
1561  }
1562 
1564  if (lenbuf>0) {
1565  // cut of the "_st" suffix
1566  Char_t *typenam = new Char_t [strlen(dscT->GetName())+1];
1567  strlcpy(typenam,dscT->GetName(),strlen(dscT->GetName())+1);
1568  // look for the last "_"
1569  Char_t *last = strrchr(typenam,'_');
1570  // Check whether it is "_st"
1571  Char_t *eon = 0;
1572  if (last) eon = strstr(last,"_st");
1573  // Cut it off if any
1574  if (eon) *eon = '\0';
1575  iOut += snprintf(strbuf+iOut,lenbuf-iOut,"struct %s {",typenam);
1576  delete [] typenam;
1577  } else {
1578  std::cout << "struct " << dscT->GetName() << " {" << std::endl;
1579  }
1580 
1581  TTableDescriptor::iterator dsc = dscT->begin();
1582  TTableDescriptor::iterator dscE = dscT->end();
1583  TDataSetIter nextComment(dscT->MakeCommentField(kFALSE));
1584  for (;dsc != dscE; dsc++) {
1586  TString name = GetTypeName(EColumnType((*dsc).fType));
1587  if (lenbuf>0) {
1588  // convert C type names to CORBA type names
1589  name.ReplaceAll("unsigned char","octet");
1590  name.ReplaceAll("int","long");
1591  iOut += snprintf(strbuf+iOut,lenbuf-iOut," %s %s",name.Data(),(*dsc).fColumnName);
1592  } else
1593  std::cout << '\t'<< name.Data() << '\t'<< (*dsc).fColumnName;
1594 
1595  Int_t indx;
1596  Int_t dim = (*dsc).fDimensions;
1597  for (indx = 0; indx < dim; indx++) {
1598  if (lenbuf>0)
1599  iOut += snprintf(strbuf+iOut,lenbuf-iOut,"[%d]",(*dsc).fIndexArray[indx]);
1600  else
1601  std::cout << "[" << std::dec << (*dsc).fIndexArray[indx]<<"]";
1602  }
1603  // print comment if any
1604  TDataSet *nxc = nextComment();
1605  if (lenbuf>0)
1606  iOut += snprintf(strbuf+iOut,lenbuf-iOut, ";");
1607  else {
1608  const char *title = nxc ? nxc->GetTitle() : " ";
1609  std::cout << ";\t//" << title << std::endl;
1610  }
1611  } /* dsc */
1612 
1614  if (lenbuf>0)
1615  iOut += snprintf(strbuf+iOut,lenbuf-iOut, "}");
1616  else
1617  std::cout << "}" << std::endl;
1618  return strbuf;
1619 }
1620 
1621 ////////////////////////////////////////////////////////////////////////////////
1622 /// Print general table inforamtion
1623 
1625 {
1626  std::cout << std::endl << " ---------------------------------------------------------------------------------------" << std::endl
1627  << " " << Path()
1628  <<" Allocated rows: "<<fN
1629  <<"\t Used rows: "<<fMaxIndex
1630  <<"\t Row size: " << fSize << " bytes"
1631  <<std::endl;
1632  return 0;
1633 }
1634 
1635 ////////////////////////////////////////////////////////////////////////////////
1636 ///const Char_t *TTable::Print(Int_t row, Int_t rownumber, const Char_t *colfirst, const Char_t *collast) const
1637 ///
1638 /// Print the contents of internal table per COLUMN.
1639 ///
1640 /// row - the index of the first row to print (counting from ZERO)
1641 /// rownumber - the total number of rows to print out (=10 by default)
1642 ///
1643 /// (No use !) Char_t *colfirst, *collast - the names of the first/last
1644 /// to print out (not implemented yet)
1645 ///
1646 ///--------------------------------------------------------------
1647 /// Check bounds and adjust it
1648 
1649 const Char_t *TTable::Print(Int_t row, Int_t rownumber, const Char_t *, const Char_t *) const
1650 {
1651  Int_t const width = 8;
1652  Int_t rowStep = 10; // The maximun values to print per line
1653  Int_t rowNumber = rownumber;
1654  if (row > Int_t(GetSize()) || GetSize() == UInt_t(0)) {
1655  PrintHeader();
1656  std::cout << " ======================================================================================" << std::endl
1657  << " There are " << GetSize() << " allocated rows for this table only" << std::endl
1658  << " ======================================================================================" << std::endl;
1659  return 0;
1660  }
1661  if (rowNumber > Int_t(GetSize()-row)) rowNumber = GetSize()-row;
1662  if (!rowNumber) return 0;
1663  rowStep = TMath::Min(rowStep,rowNumber);
1664 
1665  Int_t cdate = 0;
1666  Int_t ctime = 0;
1667  UInt_t *cdatime = 0;
1668  Bool_t isdate = kFALSE;
1669 
1671  if (!dscT ) return 0;
1672 
1673  // 3. Loop by "rowStep x lines"
1674 
1675  const Char_t *startRow = (const Char_t *)GetArray() + row*GetRowSize();
1676  Int_t rowCount = rowNumber;
1677  Int_t thisLoopLenth = 0;
1678  const Char_t *nextRow = 0;
1679  while (rowCount) {
1680  PrintHeader();
1681  if (GetNRows() == 0) {// to Print empty table header
1682  std::cout << " ======================================================================================" << std::endl
1683  << " There is NO filled row in this table" << std::endl
1684  << " ======================================================================================" << std::endl;
1685  return 0;
1686  }
1687  std::cout << " Table: " << dscT->GetName()<< "\t";
1688  for (Int_t j = row+rowNumber-rowCount; j<row+rowNumber-rowCount+rowStep && j < row+rowNumber ;j++) {
1689  Int_t hW = width-2;
1690  if (j>=10) hW -= (int)TMath::Log10(float(j))-1;
1691  std::cout << std::setw(hW) << "["<<j<<"]";
1692  std::cout << " :" ;
1693  }
1694  std::cout << std::endl
1695  << " ======================================================================================" << std::endl;
1696  TTableDescriptor::iterator member = dscT->begin();
1697  TTableDescriptor::iterator dscE = dscT->end();
1698  TDataSetIter nextComment(dscT->MakeCommentField(kFALSE));
1699 
1700  for (; member != dscE; member++){
1701  TString membertype = GetTypeName(EColumnType((*member).fType));
1702  isdate = kFALSE;
1703  if (strcmp((*member).fColumnName,"fDatime") == 0 && EColumnType((*member).fType) == kUInt)
1704  isdate = kTRUE;
1705  std::cout << membertype.Data();
1706 
1707  // Add the dimensions to "array" members
1708  Int_t dim = (*member).fDimensions;
1709  Int_t indx = 0;
1710  UInt_t *arrayLayout = 0;
1711  if (dim) {
1712  arrayLayout = new UInt_t[dim];
1713  memset(arrayLayout,0,dim*sizeof(Int_t));
1714  }
1715  Int_t arrayLength = 1;
1716  while (indx < dim ){ // Take in account the room this index will occupy
1717  arrayLength *= (*member).fIndexArray[indx];
1718  indx++;
1719  }
1720  // Encode data value or pointer value
1721  Int_t offset = (*member).fOffset;
1722  Int_t thisStepRows;
1723  thisLoopLenth = TMath::Min(rowCount,rowStep);
1724  Int_t indexOffset;
1725  Bool_t breakLoop = kFALSE;
1726 
1727  for (indexOffset=0; indexOffset < arrayLength && !breakLoop; indexOffset++) {
1728  nextRow = startRow;
1729 
1730  if (!indexOffset) std::cout << "\t" << (*member).fColumnName;
1731  else std::cout << "\t" << std::setw(strlen((*member).fColumnName)) << " ";
1732 
1733  if (dim) {
1734  for (Int_t i=0;i<dim;i++) std::cout << "["<<std::dec<<arrayLayout[i]<<"]";
1735  ArrayLayout(arrayLayout,(*member).fIndexArray,dim);
1736  }
1737  std::cout << "\t";
1738  if ( strlen((*member).fColumnName)+3*dim < 8) std::cout << "\t";
1739 
1740  for (thisStepRows = 0;thisStepRows < thisLoopLenth; thisStepRows++,nextRow += GetRowSize()) {
1741  const char *pointer = nextRow + offset + indexOffset*(*member).fTypeSize;
1742  if (isdate) {
1743  cdatime = (UInt_t*)pointer;
1744  TDatime::GetDateTime(cdatime[0],cdate,ctime);
1745  std::cout << cdate << "/" << ctime;
1746  } else if ((*member).fType == kChar && dim == 1) {
1747  char charbuffer[11];
1748  strlcpy(charbuffer,pointer,TMath::Min(10,arrayLength)+1);
1749  charbuffer[10] = 0;
1750  std::cout << "\"" << charbuffer;
1751  if (arrayLength > 10)
1752  std::cout << " . . . ";
1753  std::cout << "\"";
1754  breakLoop = kTRUE;
1755  } else {
1756  AsString((void *)pointer,EColumnType((*member).fType),width,std::cout);
1757  std::cout << " :";
1758  }
1759  }
1760  // Encode the column's comment
1761  if (indexOffset==0) {
1762  TDataSet *nxc = nextComment();
1763  std::cout << " " << (const char *)(nxc ? nxc->GetTitle() : "no comment");
1764  }
1765  std::cout << std::endl;
1766  }
1767  if (arrayLayout) delete [] arrayLayout;
1768  }
1769  rowCount -= thisLoopLenth;
1770  startRow = nextRow;
1771  }
1772  std::cout << "---------------------------------------------------------------------------------------" << std::endl;
1773  return 0;
1774 }
1775 ////////////////////////////////////////////////////////////////////////////////
1776 ///to be documented
1777 
1779 {
1782  Printf("\tclass %s: public TTable\t --> Allocated rows: %d\t Used rows: %d\t Row size: %d bytes\n",
1783  IsA()->GetName(),int(fN),int(fMaxIndex),int(fSize));
1784 
1785 }
1786 
1787 ////////////////////////////////////////////////////////////////////////////////
1788 ///*-*-*-*-*-*-*-*-*Make a projection of a TTable using selections*-*-*-*-*-*-*
1789 ///*-* =============================================
1790 ///
1791 /// Depending on the value of varexp (described in Draw) a 1-D,2-D,etc
1792 /// projection of the TTable will be filled in histogram hname.
1793 /// Note that the dimension of hname must match with the dimension of varexp.
1794 ///
1795 
1796 void TTable::Project(const char *hname, const char *varexp, const char *selection, Option_t *option,Int_t nentries, Int_t firstentry)
1797 {
1798  TString var;
1799  var.Form("%s>>%s",varexp,hname);
1800 
1801  TString opt(option);
1802  opt += "goff";
1803 
1804  Draw(var,selection,opt,nentries,firstentry);
1805 }
1806 
1807 ////////////////////////////////////////////////////////////////////////////////
1808 /// Shrink the table to free the unused but still allocated rows
1809 
1811 {
1812  ReAllocate();
1813  return TDataSet::Purge(opt);
1814 }
1815 
1816 ////////////////////////////////////////////////////////////////////////////////
1817 /// Save a primitive as a C++ statement(s) on output stream "out".
1818 
1819 void TTable::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
1820 {
1821  UInt_t arrayLayout[10],arraySize[10];
1822  const unsigned char *pointer=0,*startRow=0;
1823  int i,rowCount;unsigned char ic;
1824 
1825  out << "TDataSet *CreateTable() { " << std::endl;
1826 
1827  Int_t rowNumber = GetNRows();
1829 
1830 // Is anything Wrong??
1831  if (!rowNumber || !dscT ) {//
1832  out << "// The output table was bad-defined!" << std::endl
1833  << " fprintf(stderr, \"Bad table found. Please remove me\\n\");" << std::endl
1834  << " return 0; } " << std::endl;
1835  return;
1836  }
1837 
1838  startRow = (const UChar_t *)GetArray();
1839  assert(startRow!=0);
1840 
1841  const Char_t *rowId = "row";
1842  const Char_t *tableId = "tableSet";
1843 
1844 // Generate the header
1845 
1846  const char *className = IsA()->GetName();
1847 
1848  out << "// -----------------------------------------------------------------" << std::endl;
1849  out << "// " << Path()
1850  << " Allocated rows: "<< rowNumber
1851  <<" Used rows: "<< rowNumber
1852  <<" Row size: " << fSize << " bytes" << std::endl;
1853  out << "// " << " Table: " << dscT->GetName()<<"[0]--> "
1854  << dscT->GetName()<<"["<<rowNumber-1 <<"]" << std::endl;
1855  out << "// ====================================================================" << std::endl;
1856  out << "// ------ Test whether this table share library was loaded ------" << std::endl;
1857  out << " if (!TClass::GetClass(\"" << className << "\")) return 0;" << std::endl;
1858  out << dscT->GetName() << " " << rowId << ";" << std::endl
1859  << className << " *" << tableId << " = new "
1860  << className
1861  << "(\""<<GetName()<<"\"," << GetNRows() << ");" << std::endl
1862  << "//" <<std::endl ;
1863 
1864 // Row loop
1865  TDataSetIter nextComment(dscT->MakeCommentField(kFALSE));
1866  for (rowCount=0;rowCount<rowNumber; rowCount++,startRow += fSize, nextComment.Reset()) { //row loop
1867  out << "memset(" << "&" << rowId << ",0," << tableId << "->GetRowSize()" << ");" << std::endl ;
1868 
1869 // Member loop
1870  TTableDescriptor::iterator member = dscT->begin();
1871  TTableDescriptor::iterator dscE = dscT->end();
1872  for (; member != dscE; member++) { //LOOP over members
1873  TString memberType = GetTypeName(EColumnType((*member).fType));
1874  TString memberName((*member).fColumnName);
1875 
1876  // Encode the column's comment
1877  TDataSet *nxc = nextComment();
1878  TString memberTitle(nxc ? nxc->GetTitle() : "no comment");
1879 
1880  Int_t offset = (*member).fOffset;
1881  int mayBeName = 0;
1882  if (memberName.Index("name",0,TString::kIgnoreCase)>=0) mayBeName=1999;
1883  if (memberName.Index("file",0,TString::kIgnoreCase)>=0) mayBeName=1999;
1884  int typeSize = (*member).fTypeSize;
1885 
1886 // Add the dimensions to "array" members
1887  Int_t dim = (*member).fDimensions;
1888  if (dim) memset(arrayLayout,0,dim*sizeof(Int_t));
1889  Int_t arrayLength = 1;
1890  for (int indx=0;indx < dim ;indx++){
1891  arraySize[indx] = (*member).fIndexArray[indx];;
1892  arrayLength *= arraySize[indx];
1893  }
1894 
1895 // Special case, character array
1896  int charLen = (memberType.CompareTo("char")==0);
1897  if (charLen) { //Char case
1898  charLen=arrayLength;
1899  pointer = startRow + offset;
1900 // Actual size of char array
1901  if (mayBeName) {
1902  charLen = strlen((const char*)pointer)+1;
1903  if (charLen>arrayLength) charLen = arrayLength;
1904  } else {
1905  for(;charLen && !pointer[charLen-1];charLen--){;}
1906  if (!charLen) charLen=1;
1907  }
1908 
1909  out << " memcpy(&" << rowId << "." << (const char*)memberName;
1910  out << ",\"";
1911  for (int ii=0; ii<charLen;ii++) {
1912  ic = pointer[ii];
1913  if (ic && (isalnum(ic)
1914  || strchr("!#$%&()*+-,./:;<>=?@{}[]_|~",ic))) {//printable
1915  out << ic;
1916  } else { //nonprintable
1917  out << "\\x" << std::setw(2) << std::setfill('0') << std::hex << (unsigned)ic ;
1918  out << std::setw(1) << std::setfill(' ') << std::dec;
1919  }
1920  }
1921  out << "\"," << std::dec << charLen << ");";
1922  out << "// " << (const char*)memberTitle << std::endl;
1923  continue;
1924  } //EndIf of char case
1925 
1926 // Normal member
1927  Int_t indexOffset;
1928  for (indexOffset=0; indexOffset < arrayLength ; indexOffset++) {//array loop
1929  out << std::setw(3) << " " ;
1930  out << " " << rowId << "." << (const char*)memberName;
1931 
1932  if (dim) {
1933  for (i=0;i<dim;i++) {out << "["<<std::dec<<arrayLayout[i]<<"]";}
1934  ArrayLayout(arrayLayout,arraySize,dim);}
1935 
1936 // Generate "="
1937  out << "\t = ";
1938 
1939  pointer = startRow + offset + indexOffset*typeSize;
1940 
1941  AsString((void *)pointer,EColumnType((*member).fType),10,out);
1942 
1943 // Encode data member title
1944  if (indexOffset==0) out << "; // " << (const char*)memberTitle;
1945  out << ";" << std::endl;
1946  }//end array loop
1947  }//end of member loop
1948 
1949  out << tableId << "->AddAt(&" << rowId <<");" << std::endl;
1950 
1951  }//end of row loop
1952  out << "// ----------------- end of code ---------------" << std::endl
1953  << " return (TDataSet *)tableSet;" << std::endl
1954  << "}" << std::endl;
1955  return;
1956 }
1957 
1958 ////////////////////////////////////////////////////////////////////////////////
1959 /// Set array size of TTable object to n longs. If n<0 leave array unchanged.
1960 
1962 {
1963  if (n < 0) return;
1964  if (fN != n) Clear();
1965  SetfN(n);
1966  if (fN == 0) return;
1967  Create();
1968  if (TTable::GetNRows()) Reset();
1969 }
1970 ////////////////////////////////////////////////////////////////////////////////
1971 ///to be documented
1972 
1973 void TTable::SetTablePointer(void *table)
1974 {
1975  if (fTable) free(fTable);
1976  fTable = (Char_t *)table;
1977 }
1978 
1979 ////////////////////////////////////////////////////////////////////////////////
1980 ///to be documented
1981 
1982 void TTable::SetType(const char *const type)
1983 {
1984  SetTitle(type);
1985 }
1986 
1987 ////////////////////////////////////////////////////////////////////////////////
1988 /// Create a name of the file in the temporary directory if any
1989 
1991 {
1992  const Char_t *tempDirs = gSystem->Getenv("TEMP");
1993  if (!tempDirs) tempDirs = gSystem->Getenv("TMP");
1994  if (!tempDirs) tempDirs = "/tmp";
1995  if (gSystem->AccessPathName(tempDirs)) tempDirs = ".";
1996  if (gSystem->AccessPathName(tempDirs)) return 0;
1997  TString fileName;
1998  fileName.Form("Selection.C.%d.tmp",gSystem->GetPid());
1999  return gSystem->ConcatFileName(tempDirs,fileName.Data());
2000 }
2001 
2002 ////////////////////////////////////////////////////////////////////////////////
2003 /// Create CINT macro to evaluate the user-provided expresssion
2004 /// Expression may contains:
2005 /// - the table columen names
2006 /// - 2 meta names: i$ - the current column index,
2007 /// n$ - the total table size provided by TTable::GetNRows() method
2008 ///
2009 /// return the name of temporary file with the current expressions
2010 ///
2011 
2012 Char_t *TTable::MakeExpression(const Char_t *expressions[],Int_t nExpressions)
2013 {
2014  const Char_t *typeNames[] = {"NAN","float", "int", "long", "short", "double"
2015  ,"unsigned int","unsigned long", "unsigned short","unsigned char"
2016  ,"char", "TTableMap &"};
2017  const char *resID = "results";
2018  const char *addressID = "address";
2019  Char_t *fileName = GetExpressionFileName();
2020  if (!fileName) {
2021  Error("MakeExpression","Can not create a temporary file");
2022  return 0;
2023  }
2024 
2025  std::ofstream str;
2026  str.open(fileName);
2027  if (str.bad() ) {
2028  Error("MakeExpression","Can not open the temporary file <%s>",fileName);
2029  delete [] fileName;
2030  return 0;
2031  }
2032 
2034  const tableDescriptor_st *descTable = dsc->GetTable();
2035  // Create function
2036  str << "void SelectionQWERTY(float *"<<resID<<", float **"<<addressID<< ", int& i$, int& n$ )" << std::endl;
2037  str << "{" << std::endl;
2038  int i = 0;
2039  for (i=0; i < dsc->GetNRows(); i++,descTable++ ) {
2040  // Take the column name
2041  const Char_t *columnName = descTable->fColumnName;
2042  const Char_t *type = 0;
2043  // First check whether we do need this column
2044  for (Int_t exCount = 0; exCount < nExpressions; exCount++) {
2045  if (expressions[exCount] && expressions[exCount][0] && strstr(expressions[exCount],columnName)) goto LETSTRY;
2046  }
2047  continue;
2048 LETSTRY:
2049  Bool_t isScalar = !(descTable->fDimensions);
2050  Bool_t isFloat = descTable->fType == kFloat;
2051  type = typeNames[descTable->fType];
2052  str << type << " ";
2053  if (!isScalar) str << "*";
2054 
2055  str << columnName << " = " ;
2056  if (isScalar) str << "*(";
2057  if (!isFloat) str << "(" << type << "*)";
2058  str << addressID << "[" << i << "]";
2059  if (isScalar) str << ")" ;
2060  str << ";" << std::endl;
2061  }
2062  // Create expressions
2063  for (i=0; i < nExpressions; i++ ) {
2064  if (expressions[i] && expressions[i][0])
2065  str << " "<<resID<<"["<<i<<"]=(float)(" << expressions[i] << ");" << std::endl;
2066 // if (i == nExpressions-1 && i !=0 )
2067 // str << " if ("<<resID<<"["<<i<<"] == 0){ return; }" << std::endl;
2068  };
2069  str << "}" << std::endl;
2070  str.close();
2071  // Create byte code and check syntax
2072  if (str.good()) return fileName;
2073  delete [] fileName;
2074  return 0;
2075 }
2076 
2077 ////////////////////////////////////////////////////////////////////////////////
2078 /// Fill the entire table with byte "c" ;
2079 //// c=0 "be default"
2080 
2082 {
2083  if (fTable) {
2084  ResetMap(kTRUE);
2085  ::memset(fTable,c,fSize*fN);
2086  if (c) ResetMap(kFALSE);
2087  }
2088 }
2089 
2090 ////////////////////////////////////////////////////////////////////////////////
2091 /// Clean all filled columns with the pointers to TTableMap
2092 /// if any
2093 /// wipe = kTRUE - delete all object the Map's point to
2094 /// kFALSE - zero pointer, do not call "delete" though
2095 
2097 {
2098  piterator links = pbegin();
2099  piterator lastLinks = pend();
2100  for (;links != lastLinks;links++) {
2101  TTableMap **mp = (TTableMap **)(*links);
2102  if (wipe) delete *mp;
2103  *mp = 0;
2104  }
2105 }
2106 ////////////////////////////////////////////////////////////////////////////////
2107 /// Set array size of TTable object to n longs and copy array.
2108 /// If n<0 leave array unchanged.
2109 
2110 void TTable::Set(Int_t n, Char_t *array)
2111 {
2112  if (n < 0) return;
2113  if (fN < n) Clear();
2114 
2115  SetfN(n);
2116 
2117  if (fN == 0) return;
2118  Create();
2119  CopyStruct(fTable,array);
2120  fMaxIndex = n;
2121 }
2122 
2123 ////////////////////////////////////////////////////////////////////////////////
2124 /// Stream an object of class TTable.
2125 
2127 {
2128  if (b.IsReading()) {
2129  TDataSet::Streamer(b);
2130  b >> fN;
2131  StreamerHeader(b,version);
2132  // Create a table to fit nok rows
2133  Set(fMaxIndex);
2134  } else {
2135  TDataSet::Streamer(b);
2136  b << fN;
2137  StreamerHeader(b,version);
2138  }
2139 }
2140 
2141 ////////////////////////////////////////////////////////////////////////////////
2142 /// Read "table parameters first"
2143 
2145 {
2146  if (b.IsReading()) {
2147  Long_t rbytes;
2148  if (version) { } // version to remove compiler warning
2149 #ifdef __STAR__
2150  if (version < 3) {
2151  // skip obsolete STAR fields (for the sake of the backward compatibility)
2152  // char name[20]; /* table name */
2153  // char type[20]; /* table type */
2154  // long maxlen; /* # rows allocated */
2155  long len = b.Length() + (20+4) + (20+4) + 4;
2156  b.SetBufferOffset(len);
2157  }
2158 #endif
2159  b >> fMaxIndex; // fTableHeader->nok; /* # rows filled */
2160  b >> rbytes; /* number of bytes per row */
2161  if (GetRowSize() == -1) fSize = rbytes;
2162  if (rbytes - GetRowSize()) {
2163  Warning("StreamerHeader","Schema evolution warning: row size mismatch: expected %ld, read %ld bytes\n",GetRowSize(),rbytes);
2164  }
2165 
2166 #ifdef __STAR__
2167  if (version < 3) {
2168  // skip obsolete STAR fields (for the sake of the backward compatibility)
2169  // long dsl_pointer; /* swizzled (DS_DATASET_T*) */
2170  // long data_pointer; /* swizzled (char*) */
2171  long len = b.Length() + (4) + (4);
2172  b.SetBufferOffset(len);
2173  }
2174 #endif
2175  } else {
2176  b << fMaxIndex; //fTableHeader->nok; /* # rows filled */
2177  b << fSize; // fTableHeader->rbytes; /* number of bytes per row */
2178  }
2179 }
2180 ////////////////////////////////////////////////////////////////////////////////
2181 ///to be documented
2182 
2184 {
2185  fN = len;
2186  return fN;
2187 }
2188 ////////////////////////////////////////////////////////////////////////////////
2189 
2190 #ifdef StreamElelement
2191 #define __StreamElelement__ StreamElelement
2192 #undef StreamElelement
2193 #endif
2194 
2195 #define StreamElementIn(type) case TTableDescriptor::_NAME2_(k,type): \
2196 if (evolutionOn) { \
2197  if (nextCol->fDimensions) { \
2198  if (nextCol->fOffset != UInt_t(-1)) { \
2199  R__b.ReadFastArray((_NAME2_(type,_t) *)(row+nextCol->fOffset),nextCol->fSize/sizeof(_NAME2_(type,_t))); \
2200  } else { \
2201  _NAME2_(type,_t) *readPtrV = new _NAME2_(type,_t)[nextCol->fSize/sizeof(_NAME2_(type,_t))]; \
2202  R__b.ReadFastArray((_NAME2_(type,_t) *)(row+nextCol->fOffset),nextCol->fSize/sizeof(_NAME2_(type,_t))); \
2203  delete [] readPtrV; \
2204  readPtrV = 0; \
2205  } \
2206  } \
2207  else { \
2208  _NAME2_(type,_t) skipBuffer; \
2209  _NAME2_(type,_t) *readPtr = (_NAME2_(type,_t) *)(row+nextCol->fOffset); \
2210  if (nextCol->fOffset == UInt_t(-1)) readPtr = &skipBuffer; \
2211  R__b >> *readPtr; \
2212  } \
2213 } else { \
2214  if (nextCol->fDimensions) { \
2215  R__b.ReadFastArray ((_NAME2_(type,_t) *)(row+nextCol->fOffset),nextCol->fSize/sizeof(_NAME2_(type,_t))); \
2216  } else \
2217  R__b >> *(_NAME2_(type,_t) *)(row+nextCol->fOffset); \
2218 } \
2219 break
2220 
2221 #define StreamElementOut(type) case TTableDescriptor::_NAME2_(k,type): \
2222 if (nextCol->fDimensions) \
2223  R__b.WriteFastArray((_NAME2_(type,_t) *)(row+nextCol->fOffset), nextCol->fSize/sizeof(_NAME2_(type,_t))); \
2224 else \
2225  R__b << *(_NAME2_(type,_t) *)(row+nextCol->fOffset); \
2226 break
2227 
2228 ////////////////////////////////////////////////////////////////////////////////
2229 ///to be documented
2230 
2232 {
2233  TTableDescriptor *dsc = 0;
2234  if (IsA()) dsc = GetDescriptorPointer();
2235  if (!dsc) {
2236  Error("GetRowDescriptors()","%s has no dictionary !",GetName());
2237  dsc = GetTableDescriptors();
2238  ((TTableDescriptor *)this)->SetDescriptorPointer(dsc);
2239  }
2240  return dsc;
2241 }
2242 ////////////////////////////////////////////////////////////////////////////////
2243 ///to be documented
2244 
2246 {
2247  assert(0);
2248  return 0;
2249 }
2250 
2251 ////////////////////////////////////////////////////////////////////////////////
2252 ///to be documented
2253 
2255 {
2256  assert(0);
2257 }
2258 
2259 ////////////////////////////////////////////////////////////////////////////////
2260 /// Stream an array of the "plain" C-structures
2261 
2262 void TTable::Streamer(TBuffer &R__b)
2263 {
2264  TTableDescriptor *ioDescriptor = GetRowDescriptors();
2265  TTableDescriptor *currentDescriptor = ioDescriptor;
2266  Version_t R__v = 0;
2267  if (R__b.IsReading()) {
2268  // Check whether the file is the "obsolete" one
2269  R__v = R__b.ReadVersion();
2270  Bool_t evolutionOn = kFALSE;
2271  if (R__v>=2) {
2272  if (IsA() != TTableDescriptor::Class()) {
2273  if (R__v>3) {
2274  R__b >> ioDescriptor;
2275  } else { // backward compatibility
2276  ioDescriptor = new TTableDescriptor();
2277  ioDescriptor->Streamer(R__b);
2278  }
2279  if (!currentDescriptor) {
2280  currentDescriptor = ioDescriptor;
2281  SetDescriptorPointer(currentDescriptor);
2282  }
2283  if (currentDescriptor->fSecondDescriptor != ioDescriptor) {
2284  // Protection against of memory leak.
2285  delete currentDescriptor->fSecondDescriptor;
2286  currentDescriptor->fSecondDescriptor = ioDescriptor;
2287  }
2288 
2289  // compare two descriptors
2290  evolutionOn = (Bool_t)ioDescriptor->UpdateOffsets(currentDescriptor);
2291  }
2292  }
2293  TTable::StreamerTable(R__b,R__v);
2294  if (fMaxIndex <= 0) return;
2295  char *row= fTable;
2296  Int_t maxColumns = ioDescriptor->NumberOfColumns();
2297  Int_t rowSize = GetRowSize();
2298  if (evolutionOn) Reset(0); // Clean table
2299  for (Int_t indx=0;indx<fMaxIndex;indx++,row += rowSize) {
2300  tableDescriptor_st *nextCol = ioDescriptor->GetTable();
2301  for (Int_t colCounter=0; colCounter < maxColumns; colCounter++,nextCol++) {
2302  // Stream one table row supplied
2303  switch(nextCol->fType) {
2314  case TTableDescriptor::kPtr: {
2315  Ptr_t readPtr;
2316  R__b >> readPtr;
2317  if (evolutionOn) {
2318  // TTableMap skipBuffer;
2319  // R__b >> readPtr;
2320  if (nextCol->fOffset == UInt_t(-1)) delete readPtr; // skip this member
2321  else *(Ptr_t *)(row+nextCol->fOffset) = readPtr;
2322  } else {
2323  *(Ptr_t *)(row+nextCol->fOffset) = readPtr;
2324  }
2325  break;
2326  }
2327  default:
2328  break;
2329  };
2330  }
2331  }
2332  } else {
2333  TSeqCollection *save = fList;
2334  R__b.WriteVersion(TTable::IsA());
2335 
2336  // if (Class_Version()==2)
2337  if (IsA() != TTableDescriptor::Class()) {
2338  if ( Class_Version()>3 ) {
2339  R__b << ioDescriptor;
2340  } else { // backward compatibility
2341  ioDescriptor->Streamer(R__b);
2342  }
2343  } else {
2344  if ( Class_Version()<=3 ) fList = 0;
2345  }
2346 
2347  TTable::StreamerTable(R__b);
2348  if (fMaxIndex <= 0) return;
2349  char *row= fTable;
2350  Int_t maxColumns = ioDescriptor->NumberOfColumns();
2351  Int_t rowSize = GetRowSize();
2352  for (Int_t indx=0;indx<fMaxIndex;indx++,row += rowSize) {
2353  tableDescriptor_st *nextCol = ioDescriptor->GetTable();
2354  for (Int_t colCounter=0; colCounter < maxColumns; colCounter++,nextCol++) {
2355  // Stream one table row supplied
2356  switch(nextCol->fType) {
2368  R__b << *(Ptr_t *)(row+nextCol->fOffset);
2369  break;
2370  default:
2371  break;
2372  };
2373  }
2374  }
2375  fList = save;
2376  }
2377 }
2378 #ifdef __StreamElelement__
2379 #define StreamElelement __StreamElelement__
2380 #undef __StreamElelement__
2381 #endif
2382 
2383 ////////////////////////////////////////////////////////////////////////////////
2384 ///to be documented
2385 
2387 {
2388 }
2389 
2390 ////////////////////////////////////////////////////////////////////////////////
2391 /// Kill the table current data
2392 /// and adopt those from set
2393 
2395 {
2396  if (set->HasData()) {
2397  // Check whether the new table has the same type
2398  if (strcmp(GetTitle(),set->GetTitle()) == 0 ) {
2399  TTable *table = (TTable *)set;
2400  Adopt(table->GetSize(),table->GetArray());
2401  // Adopt can not distniguish the "allocated" and "used"
2402  // rows,
2403  // correct the corrupted number of the "used" rows
2404  SetUsedRows(table->GetNRows());
2405  // mark that object lost the STAF table and can not delete it anymore
2406  table->SetBit(kIsNotOwn);
2407  // mark we took over of this STAF table
2409  } else
2410  Error("Update",
2411  "This table is <%s> but the updating one has a wrong type <%s>",GetTitle(),set->GetTitle());
2412  }
2413  TDataSet::Update(set,opt);
2414 }
2415 ////////////////////////////////////////////////////////////////////////////////
2416 /// Query the TClass instance for the C-stucture dicitonary
2417 /// This method is to be used with TableImp CPP macro (see $ROOTSYS/table/inc/Ttypes.h
2418 
2419 const char *TTable::TableDictionary(const char *className,const char *structName,TTableDescriptor *&ColDescriptors)
2420 {
2421  if (className){/*NotUsed*/};
2422  TClass *r = TClass::GetClass(structName,1);
2423  ColDescriptors = new TTableDescriptor(r);
2424  return structName;
2425 }
2426 
2427 
2428  // ---- Table descriptor service ------
2429 Int_t TTable::GetColumnIndex(const Char_t *columnName) const {return GetRowDescriptors()->ColumnByName(columnName);}
2430 const Char_t *TTable::GetColumnName(Int_t columnIndex) const {return GetRowDescriptors()->ColumnName(columnIndex); }
2431 const UInt_t *TTable::GetIndexArray(Int_t columnIndex) const {return GetRowDescriptors()->IndexArray(columnIndex); }
2433 
2434 UInt_t TTable::GetOffset(Int_t columnIndex) const {return GetRowDescriptors()->Offset(columnIndex); }
2435 Int_t TTable::GetOffset(const Char_t *columnName) const {return GetRowDescriptors()->Offset(columnName); }
2436 
2437 UInt_t TTable::GetColumnSize(Int_t columnIndex) const {return GetRowDescriptors()->ColumnSize(columnIndex); }
2438 Int_t TTable::GetColumnSize(const Char_t *columnName) const {return GetRowDescriptors()->ColumnSize(columnName); }
2439 
2440 UInt_t TTable::GetTypeSize(Int_t columnIndex) const {return GetRowDescriptors()->TypeSize(columnIndex); }
2441 Int_t TTable::GetTypeSize(const Char_t *columnName) const {return GetRowDescriptors()->TypeSize(columnName); }
2442 
2443 UInt_t TTable::GetDimensions(Int_t columnIndex) const {return GetRowDescriptors()->Dimensions(columnIndex); }
2444 Int_t TTable::GetDimensions(const Char_t *columnName) const {return GetRowDescriptors()->Dimensions(columnName); }
2445 
2446 TTable::EColumnType TTable::GetColumnType(Int_t columnIndex) const {return GetRowDescriptors()->ColumnType(columnIndex); }
2447 TTable::EColumnType TTable::GetColumnType(const Char_t *columnName) const {return GetRowDescriptors()->ColumnType(columnName); }
2448 
2449 // pointer iterator
2450 ////////////////////////////////////////////////////////////////////////////////
2451 ///to be documented
2452 
2453 TTable::piterator::piterator(const TTable *t,EColumnType type): fCurrentRowIndex(0),fCurrentColIndex(0),fRowSize(0),fCurrentRowPtr(0),fCurrentColPtr(0)
2454 {
2455  Int_t sz = 0;
2456  if (t) sz = t->GetNRows();
2457  if (sz) {
2458  fRowSize = t->GetRowSize();
2459  fCurrentRowPtr = (const Char_t *)t->GetArray();
2460 
2461  TTableDescriptor *tabsDsc = t->GetRowDescriptors();
2462  TTableDescriptor::iterator ptr = tabsDsc->begin();
2463  TTableDescriptor::iterator lastPtr = tabsDsc->end();
2464  UInt_t i =0;
2465  for( i = 0; ptr != lastPtr; ptr++,i++)
2466  if ( tabsDsc->ColumnType(i) == type ) fPtrs.push_back(tabsDsc->Offset(i));
2467  if (fPtrs.size()==0) {
2468  MakeEnd(t->GetNRows());
2469  } else {
2470  column();
2471  }
2472  } else {
2473  MakeEnd(0);
2474  }
2475 } // piterator(TTable *)
2476 
void Add(TObject *obj, const char *name=0, Int_t check=-1)
Add object with name to browser.
Definition: TBrowser.cxx:261
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1272
void SetBufferOffset(Int_t offset=0)
Definition: TBuffer.h:88
virtual void AsString(void *buf, EColumnType type, Int_t width, std::ostream &out) const
AsString represents the value provided via "void *b" with type defined by "name". ...
Definition: TTable.cxx:228
virtual void Draw(Option_t *option="")
Draw.
Bool_t IsReading() const
Definition: TBuffer.h:81
virtual Long_t InsertRows(const void *rows, Long_t indx, UInt_t nRows=1)
void InsertRows(cons void *row, Long_t indx, UInt_t nRows)
Definition: TTable.cxx:1206
RooCmdArg Optimize(Int_t flag=2)
tableDescriptor_st * end() const
virtual int GetPid()
Get process id.
Definition: TSystem.cxx:714
float xmin
Definition: THbookFile.cxx:93
Int_t SetfN(Long_t len)
to be documented
Definition: TTable.cxx:2183
Long_t fSize
Definition: TTable.h:52
void * ReAllocate()
Reallocate this table leaving only (used rows)+1 allocated GetTableSize() = GetNRows() + 1 returns a ...
Definition: TTable.cxx:1225
short Version_t
Definition: RtypesCore.h:61
piterator pbegin()
Definition: TTable.h:260
void SetUsedRows(Int_t n)
Definition: TTable.h:286
float Float_t
Definition: RtypesCore.h:53
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8053
const char * Int
const char Option_t
Definition: RtypesCore.h:62
double Axis_t
Definition: RtypesCore.h:72
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:640
tableDescriptor_st * begin() const
unsigned int fOffset
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:131
#define TAKEACTION_BEGIN
#define BIT(n)
Definition: Rtypes.h:75
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:311
UInt_t fRowSize
Definition: TTable.h:231
virtual Bool_t IsFolder() const
return Folder flag to be used by TBrowse object The table is a folder if
Definition: TTable.cxx:1460
static void Optimize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BWID, Option_t *option="")
static function to compute reasonable axis limits
virtual TClass * GetRowClass() const
to be documented
Definition: TTable.cxx:1378
Int_t GetSize() const
Definition: TTable.h:116
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
virtual const Char_t * PrintHeader() const
Print general table inforamtion.
Definition: TTable.cxx:1624
#define gROOT
Definition: TROOT.h:375
void StreamerHeader(TBuffer &b, Version_t version=3)
Read "table parameters first".
Definition: TTable.cxx:2144
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:587
const char * Long
Basic string class.
Definition: TString.h:129
tomato 1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:551
static const char * TableDictionary()
Definition: TTable.h:263
tomato 2-D histogram with a byte per channel (see TH1 documentation)
Definition: TH2.h:132
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:168
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual void Delete(Option_t *opt="")
Delete - deletes the list of the TDataSet objects and all "Structural Members" as well This method do...
Definition: TDataSet.cxx:320
const char * Char
#define gInterpreter
Definition: TInterpreter.h:499
static const char * gDtorName
Definition: TTable.cxx:172
static Float_t gVmax[4]
Definition: TTable.cxx:175
int nbins[3]
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)=0
Profile Histogram.
Definition: TProfile.h:32
virtual Int_t UpdateOffsets(const TTableDescriptor *newDesciptor)
"Schema evolution" Method updates the offsets with a new ones from another descriptor ...
const char * UInt
#define malloc
Definition: civetweb.c:818
#define StreamElementIn(type)
Definition: TTable.cxx:2195
TTable & operator=(const TTable &rhs)
TTable assignment operator.
Definition: TTable.cxx:1086
virtual const Char_t * GetType() const
Returns the type of the wrapped C-structure kept as the TNamed title.
Definition: TTable.cxx:1448
virtual void Browse(TBrowser *b)
Browse this dataset (called by TBrowser).
Definition: TDataSet.cxx:297
virtual void Browse(TBrowser *b)
Wrap each table coulumn with TColumnView object to browse.
Definition: TTable.cxx:1297
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:687
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
Definition: TH1.cxx:6419
virtual TTableDescriptor * GetDescriptorPointer() const
to be documented
Definition: TTable.cxx:2245
if object in a list can be deleted
Definition: TObject.h:58
Int_t Length() const
Definition: TBuffer.h:94
virtual Char_t * Print(Char_t *buf, Int_t n) const
Create IDL table defintion (to be used for XDF I/O)
Definition: TTable.cxx:1552
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TTable.cxx:1819
virtual Int_t GetDimension() const
Definition: TH1.h:263
const UInt_t * IndexArray(Int_t columnIndex) const
Sequenceable collection abstract base class.
Double_t GetXmin() const
Definition: TAxis.h:133
const Char_t * ColumnName(Int_t columnIndex) const
static Int_t gNbins[4]
Definition: TTable.cxx:173
Double_t x[n]
Definition: legend1.C:17
TTable(const char *name=0, Int_t size=0)
Default TTable ctor.
Definition: TTable.cxx:1034
void Class()
Definition: Class.C:29
virtual void Sleep(UInt_t milliSec)
Sleep milliSec milli seconds.
Definition: TSystem.cxx:444
virtual const UInt_t * GetIndexArray(Int_t columnIndex) const
Definition: TTable.cxx:2431
static void ArrayLayout(UInt_t *layout, const UInt_t *size, Int_t dim)
ArrayLayout - calculates the array layout recursively.
Definition: TTable.cxx:197
virtual void SetNRows(Int_t n)
Definition: TTable.h:288
Double_t Log10(Double_t x)
Definition: TMath.h:652
UInt_t Offset(Int_t columnIndex) const
static void GetDateTime(UInt_t datetime, Int_t &date, Int_t &time)
Static function that returns the date and time.
Definition: TDatime.cxx:434
virtual const char * Getenv(const char *env)
Get environment variable.
Definition: TSystem.cxx:1634
EColumnType
Definition: TTable.h:82
Int_t Finite(Double_t x)
Definition: TMath.h:655
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTable.h:308
TH1F * h1
Definition: legend1.C:5
#define realloc
Definition: civetweb.c:820
TSeqCollection * fList
Definition: TDataSet.h:59
UInt_t ColumnSize(Int_t columnIndex) const
static Float_t gVmin[4]
Definition: TTable.cxx:174
TClass * RowClass() const
virtual void Set(Int_t n)
Set array size of TTable object to n longs. If n<0 leave array unchanged.
Definition: TTable.cxx:1961
Bool_t OutOfBoundsError(const char *where, Int_t i) const
Generate an out-of-bounds error. Always returns false.
Definition: TTable.cxx:1544
Int_t CopyRows(const TTable *srcTable, Long_t srcRow=0, Long_t dstRow=0, Long_t nRows=0, Bool_t expand=kFALSE)
CopyRows copies nRows from starting from the srcRow of srcTable to the dstRow in this table upto nRow...
Definition: TTable.cxx:331
const char * UChar
const char * Float
A specialized string object used for TTree selections.
Definition: TCut.h:25
const void * At(Int_t i) const
Returns a pointer to the i-th row of the table.
Definition: TTable.cxx:302
virtual Bool_t EntryLoop(const Char_t *exprFileName, Int_t &action, TObject *obj, Int_t nentries=1000000000, Int_t firstentry=0, Option_t *option="")
EntryLoop creates a CINT bytecode to evaluate the given expressions for all table rows in loop and fi...
Definition: TTable.cxx:788
virtual void Delete(Option_t *option="")
Delete this object.
Definition: TObject.cxx:176
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
TTableDescriptor * fSecondDescriptor
virtual TTableDescriptor * GetRowDescriptors() const
to be documented
Definition: TTable.cxx:2231
UInt_t TypeSize(Int_t columnIndex) const
virtual void SetType(const char *const type)
to be documented
Definition: TTable.cxx:1982
piterator pend()
Definition: TTable.h:261
TRandom2 r(17)
virtual const Char_t * GetColumnComment(Int_t columnIndex) const
Get a comment from the table descriptor.
Definition: TTable.cxx:1169
R__EXTERN TSystem * gSystem
Definition: TSystem.h:539
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2851
SVector< double, 2 > v
Definition: Dict.h:5
#define TAKEACTION_END
tomato 2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:249
virtual Int_t AddAt(const void *c)
Add the "row" at the GetNRows() position, and reallocate the table if neccesary, and return the row i...
Definition: TTable.cxx:1125
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:436
virtual Long_t AppendRows(const void *row, UInt_t nRows)
Append nRows row of the array "row" to the table return.
Definition: TTable.cxx:1181
virtual UInt_t GetDimensions(Int_t columnIndex) const
Definition: TTable.cxx:2443
virtual void Fit(const char *formula, const char *varexp, const char *selection="", Option_t *option="", Option_t *goption="", Int_t nentries=1000000000, Int_t firstentry=0)
-*-*-*-*-*-*-*-*Fit a projected item(s) from a TTable-*-*-*-*-*-*-*-*-* *-* =========================...
Definition: TTable.cxx:1426
tableDescriptor_st * GetTable(Int_t i=0) const
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2332
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
Char_t * fTable
Definition: TTable.h:57
void ** column()
Definition: TTable.h:236
Ssiz_t Length() const
Definition: TString.h:388
const char * ULong
virtual const Char_t * GetColumnName(Int_t columnIndex) const
Definition: TTable.cxx:2430
virtual void Update()
to be documented
Definition: TTable.cxx:2386
Int_t NaN()
return the total number of the NaN for float/double cells of this table Thanks Victor Perevoztchikov ...
Definition: TTable.cxx:1478
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition: TString.cxx:1080
virtual void ResetMap(Bool_t wipe=kTRUE)
Clean all filled columns with the pointers to TTableMap if any wipe = kTRUE - delete all object the M...
Definition: TTable.cxx:2096
virtual void DeleteRows(Long_t indx, UInt_t nRows=1)
Delete one or several rows from the table.
Definition: TTable.cxx:363
A TEventList object is a list of selected events (entries) in a TTree.
Definition: TEventList.h:31
std::vector< ULong_t > fPtrs
Definition: TTable.h:228
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2632
const char * Double
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:71
TAxis * GetYaxis()
Definition: TH1.h:301
float xmax
Definition: THbookFile.cxx:93
#define StreamElementOut(type)
Definition: TTable.cxx:2221
virtual void Draw(Option_t *option="")
Draws 3-D polymarker with its current attributes.
static TH1 * gCurrentTableHist
Definition: TTable.cxx:170
void * GetArray() const
Definition: TTable.h:280
#define Printf
Definition: TGeoToOCC.h:18
char * StrDup(const char *str)
Duplicate the string str.
Definition: TString.cxx:2524
virtual TString Path() const
return the full path of this data set
Definition: TDataSet.cxx:626
const Bool_t kFALSE
Definition: RtypesCore.h:92
virtual ~TTable()
Delete TTable object.
Definition: TTable.cxx:1101
long Long_t
Definition: RtypesCore.h:50
virtual Int_t SetNextPoint(Double_t x, Double_t y, Double_t z)
Set point following LastPoint to x, y, z.
virtual UInt_t GetColumnSize(Int_t columnIndex) const
Definition: TTable.cxx:2437
A PolyMarker is defined by an array on N points in a 2-D space.
Definition: TPolyMarker.h:31
#define ClassImp(name)
Definition: Rtypes.h:336
void StreamerTable(TBuffer &b, Version_t version=3)
Stream an object of class TTable.
Definition: TTable.cxx:2126
virtual TTableDescriptor * GetTableDescriptors() const
protected: create a new TTableDescriptor descriptor for this table
Definition: TTable.cxx:213
Bool_t BoundsOk(const char *where, Int_t at) const
Definition: TTable.h:272
double Double_t
Definition: RtypesCore.h:55
virtual Int_t Purge(Option_t *opt="")
Purge - deletes all "dummy" "Structural Members" those are not ended up with some dataset with data i...
Definition: TDataSet.cxx:758
int type
Definition: TGX11.cxx:120
#define free
Definition: civetweb.c:821
unsigned long ULong_t
Definition: RtypesCore.h:51
int nentries
Definition: THbookFile.cxx:89
void SetTablePointer(void *table)
to be documented
Definition: TTable.cxx:1973
TTable::EColumnType ColumnType(Int_t columnIndex) const
Double_t y[n]
Definition: legend1.C:17
double func(double *x, double *p)
Definition: stressTF1.cxx:213
void CopyStruct(Char_t *dest, const Char_t *src)
Copy the C-structure src into the new location the length of the strucutre is defined by this class d...
Definition: TTable.cxx:1154
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:572
The TH1 histogram class.
Definition: TH1.h:56
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
static EColumnType GetTypeId(const char *typeName)
return the Id of the C basic type by given name return kNAN if the name provided fits no knwn basic n...
Definition: TTable.cxx:291
static TTable * New(const Char_t *name, const Char_t *type, void *array, UInt_t size)
This static method creates a new TTable object if provided.
Definition: TTable.cxx:1519
tomato 2-D histogram with a short per channel (see TH1 documentation)
Definition: TH2.h:171
unsigned int fDimensions
static Char_t * GetExpressionFileName()
Create a name of the file in the temporary directory if any.
Definition: TTable.cxx:1990
int CompareTo(const char *cs, ECaseCompare cmp=kExact) const
Compare a string to char *cs2.
Definition: TString.cxx:396
virtual Int_t GetColumnIndex(const Char_t *columnName) const
Definition: TTable.cxx:2429
UInt_t NumberOfColumns() const
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition: TClass.cxx:2885
double Stat_t
Definition: RtypesCore.h:73
Definition: TTable.h:48
static void FindGoodLimits(Int_t nbins, Int_t &newbins, Float_t &xmin, Float_t &xmax)
-*-*-*-*-*-*-*-*Find reasonable bin values-*-*-*-*-*-*-*-*-*-*-*-*-*-* *-* ==========================...
Definition: TTable.cxx:738
virtual Long_t GetTableSize() const
Returns the number of the allocated rows.
Definition: TTable.cxx:1405
static TView * CreateView(Int_t system=1, const Double_t *rmin=0, const Double_t *rmax=0)
Create a concrete default 3-d view via the plug-in manager.
Definition: TView.cxx:39
const Char_t * fCurrentRowPtr
Definition: TTable.h:232
const char * Bool
Mother of all ROOT objects.
Definition: TObject.h:37
char Char_t
Definition: RtypesCore.h:29
virtual Char_t * MakeExpression(const Char_t *expressions[], Int_t nExpressions)
Create CINT macro to evaluate the user-provided expresssion Expression may contains: ...
Definition: TTable.cxx:2012
virtual void PrintContents(Option_t *opt="") const
to be documented
Definition: TTable.cxx:1778
A 3D polymarker.
Definition: TPolyMarker3D.h:32
static const char * GetTypeName(EColumnType type)
return table type name
Definition: TTable.cxx:281
virtual void Update()
Update()
Definition: TDataSet.cxx:864
virtual void Project(const char *hname, const char *varexp, const char *selection="", Option_t *option="", Int_t nentries=1000000000, Int_t firstentry=0)
-*-*-*-*-*-*-*-*Make a projection of a TTable using selections-*-*-*-*-*-* *-* ======================...
Definition: TTable.cxx:1796
const TTable * Table()
Definition: TTableMap.h:51
#define dest(otri, vertexptr)
Definition: triangle.c:1040
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
virtual Long_t GetNRows() const
Returns the number of the used rows for the wrapped table.
Definition: TTable.cxx:1391
Char_t * Create()
Allocate a space for the new table, if any Sleep for a while if space is not available and try again...
Definition: TTable.cxx:1271
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition: TH1.cxx:8132
void MakeEnd(UInt_t lastRowIndex)
Definition: TTable.h:374
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
#define snprintf
Definition: civetweb.c:822
virtual UInt_t GetTypeSize(Int_t columnIndex) const
Definition: TTable.cxx:2440
#define gPad
Definition: TVirtualPad.h:284
virtual Long_t GetRowSize() const
Returns the size (in bytes) of one table row.
Definition: TTable.cxx:1398
Int_t fN
Definition: TTable.h:56
virtual void SetDescriptorPointer(TTableDescriptor *list)
to be documented
Definition: TTable.cxx:2254
piterator(const TTable *t=0, EColumnType type=kPtr)
to be documented
Definition: TTable.cxx:2453
virtual void Delete(Option_t *opt="")
Delete the internal array and free the memory it occupied if this object did own this array...
Definition: TTable.cxx:1369
#define gDirectory
Definition: TDirectory.h:211
void ResetBit(UInt_t f)
Definition: TObject.h:158
unsigned char UChar_t
Definition: RtypesCore.h:34
virtual EColumnType GetColumnType(Int_t columnIndex) const
Definition: TTable.cxx:2446
Int_t GetNbins() const
Definition: TAxis.h:121
Int_t ColumnByName(const Char_t *columnName=0) const
Find the column index but the column name.
virtual void Clear(Option_t *opt="")
Deletes the internal array of this class if this object does own its internal table.
Definition: TTable.cxx:1345
static const char * fgTypeName[kEndColumnType]
Definition: TTable.h:89
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
THist< 2, float, THistStatContent, THistStatUncertainty > TH2F
Definition: THist.hxx:317
virtual void PrintContents(Option_t *opt="") const
Callback method to complete ls() method recursive loop This is to allow to sepoarate navigation and t...
Definition: TDataSet.cxx:618
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3564
const Bool_t kTRUE
Definition: RtypesCore.h:91
UInt_t Dimensions(Int_t columnIndex) const
virtual char * ConcatFileName(const char *dir, const char *name)
Concatenate a directory and a file name. User must delete returned string.
Definition: TSystem.cxx:1051
virtual void Adopt(Int_t n, void *array)
Adopt array arr into TTable, i.e.
Definition: TTable.cxx:1110
Double_t GetXmax() const
Definition: TAxis.h:134
TDataSet * MakeCommentField(Bool_t createFlag=kTRUE)
Instantiate a comment dataset if any.
const Int_t n
Definition: legend1.C:16
virtual UInt_t GetNumberOfColumns() const
Definition: TTable.cxx:2432
Long_t fMaxIndex
Definition: TTable.h:58
virtual void CopySet(TTable &array)
to be documented
Definition: TTable.cxx:1161
virtual UInt_t GetOffset(Int_t columnIndex) const
Definition: TTable.cxx:2434
void ReAlloc(Int_t newsize)
The table is reallocated if it is an owner of the internal array.
Definition: TTable.cxx:1246
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:859
TAxis * GetXaxis()
Definition: TH1.h:300
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual void Reset(Int_t c=0)
Fill the entire table with byte "c" ; / c=0 "be default".
Definition: TTable.cxx:2081
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual TObject * Last() const =0
virtual void Reset(Option_t *option="")
Reset number of entries in event list.
Definition: TEventList.cxx:328
virtual Int_t Purge(Option_t *opt="")
Shrink the table to free the unused but still allocated rows.
Definition: TTable.cxx:1810
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition: TClass.cxx:4706
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:290
const char * Data() const
Definition: TString.h:347
const char * Short