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