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