Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TSelectorDraw.cxx
Go to the documentation of this file.
1// @(#)root/treeplayer:$Id$
2// Author: Rene Brun 08/01/2003
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/** \class TSelectorDraw
13A specialized TSelector for TTree::Draw.
14*/
15
16#include "TSelectorDraw.h"
17#include "TROOT.h"
18#include "TH2.h"
19#include "TH3.h"
20#include "TView.h"
21#include "TGraph.h"
22#include "TPolyMarker3D.h"
23#include "TDirectory.h"
24#include "TVirtualPad.h"
25#include "TProfile.h"
26#include "TProfile2D.h"
27#include "TTreeFormulaManager.h"
28#include "TEnv.h"
29#include "TTree.h"
30#include "TCut.h"
31#include "TEntryList.h"
32#include "TEventList.h"
33#include "TEntryListArray.h"
34#include "THLimitsFinder.h"
35#include "TStyle.h"
36#include "TClass.h"
37#include "TColor.h"
38#include "strlcpy.h"
39
41
43
44////////////////////////////////////////////////////////////////////////////////
45/// Default selector constructor.
46
48{
49 fTree = nullptr;
50 fW = nullptr;
51 fValSize = 4;
52 fVal = new Double_t*[fValSize];
53 fVmin = new Double_t[fValSize];
54 fVmax = new Double_t[fValSize];
55 fNbins = new Int_t[fValSize];
56 fVarMultiple = new bool[fValSize];
58 for (Int_t i = 0; i < fValSize; ++i) {
59 fVal[i] = nullptr;
60 fVar[i] = nullptr;
61 }
62 fManager = nullptr;
63 fMultiplicity = 0;
64 fSelect = nullptr;
65 fSelectedRows = 0;
66 fDraw = 0;
67 fObject = nullptr;
68 fOldHistogram = nullptr;
69 fObjEval = false;
70 fSelectMultiple = false;
71 fCleanElist = false;
72 fTreeElist = nullptr;
73 fAction = 0;
74 fNfill = 0;
75 fDimension = 0;
76 fOldEstimate = 0;
77 fForceRead = 0;
78 fWeight = 1;
80 fTreeElistArray = nullptr;
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Selector destructor.
85
87{
89 delete [] fVar;
90 if (fVal) {
91 for (Int_t i = 0; i < fValSize; ++i)
92 delete [] fVal[i];
93 delete [] fVal;
94 }
95 if (fVmin) delete [] fVmin;
96 if (fVmax) delete [] fVmax;
97 if (fNbins) delete [] fNbins;
98 if (fVarMultiple) delete [] fVarMultiple;
99 if (fW) delete [] fW;
100}
101
102////////////////////////////////////////////////////////////////////////////////
103/// Called every time a loop on the tree(s) starts.
104
106{
107 SetStatus(0);
108 ResetAbort();
110 fSelectedRows = 0;
111 fTree = tree;
112 fDimension = 0;
113 fAction = 0;
114
115 TObject *obj = fInput->FindObject("varexp");
116 const char *varexp0 = obj ? obj->GetTitle() : "";
117 obj = fInput->FindObject("selection");
118 const char *selection = obj ? obj->GetTitle() : "";
119 const char *option = GetOption();
120
121 TString opt, abrt;
122 char *hdefault = (char *)"htemp";
123 char *varexp = nullptr;
124 Int_t i, j, hkeep;
125 opt = option;
126 opt.ToLower();
127 fOldHistogram = nullptr;
128 TEntryList *enlist = nullptr;
129 TEventList *evlist = nullptr;
130 TString htitle;
131 bool profile = false;
132 bool optSame = false;
133 bool optEnlist = false;
134 bool optEnlistArray = false;
135 bool optpara = false;
136 bool optcandle = false;
137 bool opt5d = false;
138 if (opt.Contains("same")) {
139 optSame = true;
140 opt.ReplaceAll("same", "");
141 }
142 if (opt.Contains("entrylist")) {
143 optEnlist = true;
144 if (opt.Contains("entrylistarray")) {
145 optEnlistArray = true;
146 opt.ReplaceAll("entrylistarray", "");
147 } else {
148 opt.ReplaceAll("entrylist", "");
149 }
150 }
151 if (opt.Contains("para")) {
152 optpara = true;
153 opt.ReplaceAll("para", "");
154 }
155 if (opt.Contains("candle")) {
156 optcandle = true;
157 opt.ReplaceAll("candle", "");
158 }
159 if (opt.Contains("gl5d")) {
160 opt5d = true;
161 opt.ReplaceAll("gl5d", "");
162 }
163 TCut realSelection(selection);
164 //input list - only TEntryList
165 TEntryList *inElist = fTree->GetEntryList();
166 evlist = fTree->GetEventList();
167 if (evlist && inElist) {
168 //this is needed because the input entry list was created
169 //by the fTree from the input TEventList and is owned by the fTree.
170 //Calling GetEntryList function changes ownership and here
171 //we want fTree to still delete this entry list
172
173 inElist->SetBit(kCanDelete, true);
174 }
175 fCleanElist = false;
176 fTreeElist = inElist;
177
178 fTreeElistArray = inElist ? dynamic_cast<TEntryListArray*>(fTreeElist) : nullptr;
179
180
181 if (inElist && inElist->GetReapplyCut()) {
182 realSelection *= inElist->GetTitle();
183 }
184
185 // what each variable should contain:
186 // varexp0 - original expression eg "a:b>>htest"
187 // hname - name of new or old histogram
188 // hkeep - flag if to keep new histogram
189 // hnameplus - flag if to add to current histo
190 // i - length of variable expression stipped of everything after ">>"
191 // varexp - variable expression stipped of everything after ">>"
192 // fOldHistogram - pointer to hist hname
193 // elist - pointer to selection list of hname
194
195 bool canExtend = true;
196 if (optSame) canExtend = false;
197
198 Int_t nbinsx = 0, nbinsy = 0, nbinsz = 0;
199 Double_t xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
200
201 fObject = nullptr;
202 char *hname = nullptr;
203 TString hnamealloc;
204 i = 0;
205 if (varexp0 && strlen(varexp0)) {
206 for (UInt_t k = strlen(varexp0) - 1; k > 0; k--) {
207 if (varexp0[k] == '>' && varexp0[k-1] == '>') {
208 i = (int)(&(varexp0[k-1]) - varexp0); // length of varexp0 before ">>"
209 hnamealloc = &(varexp0[k+1]); // use TString to copy value of the
210 hname = (char *) hnamealloc.Data();
211 break;
212 }
213 }
214 }
215 // char *hname = (char*)strstr(varexp0,">>");
216 if (hname) {
217 hkeep = 1;
218 varexp = new char[i+1];
219 varexp[0] = 0; //necessary if i=0
220 bool hnameplus = false;
221 while (*hname == ' ') hname++;
222 if (*hname == '+') {
223 hnameplus = true;
224 hname++;
225 while (*hname == ' ') hname++; //skip ' '
226 }
227 j = strlen(hname) - 1; // skip ' ' at the end
228 while (j) {
229 if (hname[j] != ' ') break;
230 hname[j] = 0;
231 j--;
232 }
233
234 if (i) {
235 strlcpy(varexp,varexp0,i+1);
236
237 Int_t mustdelete = 0;
239
240 // parse things that follow the name of the histo between '(' and ')'.
241 // At this point hname contains the name of the specified histogram.
242 // Now the syntax is extended to handle an hname of the following format
243 // hname(nBIN [[,[xlow]][,xhigh]],...)
244 // so enclosed in brackets is the binning information, xlow, xhigh, and
245 // the same for the other dimensions
246
247 char *pstart; // pointer to '('
248 char *pend; // pointer to ')'
249 char *cdummy; // dummy pointer
250 int ncomma; // number of commas between '(' and ')', later number of arguments
251 int ncols; // number of columns in varexpr
252 Double_t value; // parsed value (by sscanf)
253
254 const Int_t maxvalues = 9;
255
256 pstart = strchr(hname, '(');
257 pend = strchr(hname, ')');
258 if (pstart != nullptr) { // found the bracket
259
260 mustdelete = 1;
261
262 // check that there is only one open and close bracket
263 if (pstart == strrchr(hname, '(') && pend == strrchr(hname, ')')) {
264
265 // count number of ',' between '(' and ')'
266 ncomma = 0;
267 cdummy = pstart;
268 cdummy = strchr(&cdummy[1], ',');
269 while (cdummy != nullptr) {
270 cdummy = strchr(&cdummy[1], ',');
271 ncomma++;
272 }
273
274 if (ncomma + 1 > maxvalues) {
275 Error("DrawSelect", "ncomma+1>maxvalues, ncomma=%d, maxvalues=%d", ncomma, maxvalues);
276 ncomma = maxvalues - 1;
277 }
278
279 ncomma++; // number of arguments
280 cdummy = pstart;
281
282 // number of columns
283 ncols = 1;
284 for (j = 0; j < i; j++) {
285 if (varexp[j] == ':'
286 && !((j > 0 && varexp[j-1] == ':') || varexp[j+1] == ':')
287 ) {
288 ncols++;
289 }
290 }
291 if (ncols > 3) { // max 3 columns
292 Error("DrawSelect", "ncols > 3, ncols=%d", ncols);
293 ncols = 0;
294 }
295
296 // check dimensions before and after ">>"
297 if (ncols * 3 < ncomma) {
298 Error("DrawSelect", "ncols*3 < ncomma ncols=%d, ncomma=%d", ncols, ncomma);
299 ncomma = ncols * 3;
300 }
301
302 // scan the values one after the other
303 for (j = 0; j < ncomma; j++) {
304 cdummy++; // skip '(' or ','
305 if (sscanf(cdummy, " %lf ", &value) == 1) {
306 cdummy = strchr(&cdummy[1], ',');
307
308 switch (j) { // do certain settings depending on position of argument
309 case 0: // binning x-axis
310 nbinsx = (Int_t)value;
311 if (ncols < 2) {
312 gEnv->SetValue("Hist.Binning.1D.x", nbinsx);
313 } else if (ncols < 3) {
314 gEnv->SetValue("Hist.Binning.2D.x", nbinsx);
315 gEnv->SetValue("Hist.Binning.2D.Prof", nbinsx);
316 } else {
317 gEnv->SetValue("Hist.Binning.3D.x", nbinsx);
318 gEnv->SetValue("Hist.Binning.3D.Profx", nbinsx);
319 }
320
321 break;
322 case 1: // lower limit x-axis
323 xmin = value;
324 break;
325 case 2: // upper limit x-axis
326 xmax = value;
327 break;
328 case 3: // binning y-axis
329 nbinsy = (Int_t)value;
330 if (ncols < 3) gEnv->SetValue("Hist.Binning.2D.y", nbinsy);
331 else {
332 gEnv->SetValue("Hist.Binning.3D.y", nbinsy);
333 gEnv->SetValue("Hist.Binning.3D.Profy", nbinsy);
334 }
335 break;
336 case 4: // lower limit y-axis
337 ymin = value;
338 break;
339 case 5: // upper limit y-axis
340 ymax = value;
341 break;
342 case 6: // binning z-axis
343 nbinsz = (Int_t)value;
344 gEnv->SetValue("Hist.Binning.3D.z", nbinsz);
345 break;
346 case 7: // lower limit z-axis
347 zmin = value;
348 break;
349 case 8: // upper limit z-axis
350 zmax = value;
351 break;
352 default:
353 Error("DrawSelect", "j>8");
354 break;
355 }
356 } // if sscanf == 1
357 } // for j=0;j<ncomma;j++
358 } else {
359 Error("Begin", "Two open or close brackets found, hname=%s", hname);
360 }
361
362 // fix up hname
363 pstart[0] = '\0'; // removes things after (and including) '('
364 } // if '(' is found
365
366 j = strlen(hname) - 1; // skip ' ' at the end
367 while (j > 0) {
368 if (hname[j] != ' ') break; // skip ' ' at the end
369 hname[j] = 0;
370 j--;
371 }
372
373 TObject *oldObject = gDirectory->Get(hname); // if hname contains '(...)' the return values is NULL, which is what we want
374 fOldHistogram = oldObject ? dynamic_cast<TH1*>(oldObject) : nullptr;
375
376 if (!fOldHistogram && oldObject && !oldObject->InheritsFrom(TH1::Class())) {
377 abrt.Form("An object of type '%s' has the same name as the requested histo (%s)", oldObject->IsA()->GetName(), hname);
378 Abort(abrt);
379 delete[] varexp;
380 return;
381 }
382 if (fOldHistogram && !hnameplus) fOldHistogram->Reset(); // reset unless adding is wanted
383
384 if (mustdelete) {
385 if (gDebug) {
386 Warning("Begin", "Deleting old histogram, since (possibly new) limits and binnings have been given");
387 }
388 delete fOldHistogram; fOldHistogram=nullptr;
389 }
390
391 } else {
392 // make selection list (i.e. varexp0 starts with ">>")
393 TObject *oldObject = gDirectory->Get(hname);
394 if (optEnlist) {
395 //write into a TEntryList
396 enlist = oldObject ? dynamic_cast<TEntryList*>(oldObject) : nullptr;
397
398 if (!enlist && oldObject) {
399 abrt.Form("An object of type '%s' has the same name as the requested event list (%s)",
400 oldObject->IsA()->GetName(), hname);
401 Abort(abrt);
402 delete[] varexp;
403 return;
404 }
405 if (!enlist) {
406 if (optEnlistArray) {
407 enlist = new TEntryListArray(hname, realSelection.GetTitle());
408 } else {
409 enlist = new TEntryList(hname, realSelection.GetTitle());
410 }
411 }
412 if (enlist) {
413 if (!hnameplus) {
414 if (enlist == inElist) {
415 // We have been asked to reset the input list!!
416 // Let's set it aside for now ...
417 if (optEnlistArray) {
418 inElist = new TEntryListArray(*enlist);
419 } else {
420 inElist = new TEntryList(*enlist);
421 }
422 fCleanElist = true;
423 fTree->SetEntryList(inElist);
424 }
425 enlist->Reset();
426 enlist->SetTitle(realSelection.GetTitle());
427 } else {
428 TCut old = enlist->GetTitle();
429 TCut upd = old || realSelection.GetTitle();
430 enlist->SetTitle(upd.GetTitle());
431 }
432 }
433 } else {
434 //write into a TEventList
435 evlist = oldObject ? dynamic_cast<TEventList*>(oldObject) : nullptr;
436
437 if (!evlist && oldObject) {
438 abrt.Form("An object of type '%s' has the same name as the requested event list (%s)",
439 oldObject->IsA()->GetName(), hname);
440 Abort(abrt);
441 delete[] varexp;
442 return;
443 }
444 if (!evlist) {
445 evlist = new TEventList(hname, realSelection.GetTitle(), 1000, 0);
446 }
447 if (evlist) {
448 if (!hnameplus) {
449 if (evlist == fTree->GetEventList()) {
450 // We have been asked to reset the input list!!
451 // Let's set it aside for now ...
452 Abort("Input and output lists are the same!");
453 delete[] varexp;
454 return;
455 }
456 evlist->Reset();
457 evlist->SetTitle(realSelection.GetTitle());
458 } else {
459 TCut old = evlist->GetTitle();
460 TCut upd = old || realSelection.GetTitle();
461 evlist->SetTitle(upd.GetTitle());
462 }
463 }
464 }
465
466 } // if (i)
467 } else { // if (hname)
468 hname = hdefault;
469 hkeep = 0;
470 const size_t varexpLen = strlen(varexp0) + 1;
471 varexp = new char[varexpLen];
472 strlcpy(varexp, varexp0, varexpLen);
473 if (gDirectory) {
474 fOldHistogram = (TH1*)gDirectory->Get(hname);
475 if (fOldHistogram) { fOldHistogram->Delete(); fOldHistogram = nullptr;}
476 }
477 }
478
479 // Decode varexp and selection
480 if (!CompileVariables(varexp, realSelection.GetTitle())) {
481 abrt.Form("Variable compilation failed: {%s,%s}", varexp, realSelection.GetTitle());
482 Abort(abrt);
483 delete[] varexp;
484 return;
485 }
486 if (fDimension > 4 && !(optpara || optcandle || opt5d || opt.Contains("goff"))) {
487 Abort("Too many variables. Use the option \"para\", \"gl5d\" or \"candle\" to display more than 4 variables.");
488 delete[] varexp;
489 return;
490 }
491 if (fDimension < 2 && (optpara || optcandle)) {
492 Abort("The options \"para\" and \"candle\" require at least 2 variables.");
493 delete[] varexp;
494 return;
495 }
496
497 // In case fOldHistogram exists, check dimensionality
498 Int_t nsel = strlen(selection);
499 if (nsel > 1) {
500 htitle.Form("%s {%s}", varexp, selection);
501 } else {
502 htitle = varexp;
503 }
504 if (fOldHistogram) {
505 Int_t olddim = fOldHistogram->GetDimension();
506 Int_t mustdelete = 0;
508 profile = true;
509 olddim = 2;
510 }
512 profile = true;
513 olddim = 3;
514 }
515 if (opt.Contains("prof") && fDimension > 1) {
516 // ignore "prof" for 1D.
517 if (!profile || olddim != fDimension) mustdelete = 1;
518 } else if (opt.Contains("col") && fDimension>2) {
519 if (olddim+1 != fDimension) mustdelete = 1;
520 } else {
521 if (olddim != fDimension) mustdelete = 1;
522 }
523 if (mustdelete) {
524 Warning("Begin", "Deleting old histogram with different dimensions");
525 delete fOldHistogram;
526 fOldHistogram = nullptr;
527 }
528 }
529
530 // Create a default canvas if none exists
531 fDraw = 0;
532 if (!gPad && !opt.Contains("goff") && fDimension > 0) {
533 gROOT->MakeDefCanvas();
534 if (!gPad) {
535 Abort("Creation of default canvas failed");
536 delete[] varexp;
537 return;
538 }
539 }
540
541 // 1-D distribution
542 TH1 *hist;
543 if (fDimension == 1) {
544 fAction = 1;
545 if (!fOldHistogram) {
546 fNbins[0] = gEnv->GetValue("Hist.Binning.1D.x", 100);
547 if (gPad && optSame) {
548 TListIter np(gPad->GetListOfPrimitives());
549 TObject *op;
550 TH1 *oldhtemp = nullptr;
551 while ((op = np()) && !oldhtemp) {
552 if (op->InheritsFrom(TH1::Class())) oldhtemp = (TH1 *)op;
553 }
554 if (oldhtemp) {
555 fNbins[0] = oldhtemp->GetXaxis()->GetNbins();
556 fVmin[0] = oldhtemp->GetXaxis()->GetXmin();
557 fVmax[0] = oldhtemp->GetXaxis()->GetXmax();
558 } else {
559 fVmin[0] = gPad->GetUxmin();
560 fVmax[0] = gPad->GetUxmax();
561 }
562 } else {
563 fAction = -1;
564 fVmin[0] = xmin;
565 fVmax[0] = xmax;
566 if (xmin < xmax) canExtend = false;
567 }
568 }
569 if (fOldHistogram) {
570 hist = fOldHistogram;
571 fNbins[0] = hist->GetXaxis()->GetNbins();
572 } else {
573 TString precision = gEnv->GetValue("Hist.Precision.1D", "float");
574 if (precision.Contains("float")) {
575 hist = new TH1F(hname, htitle.Data(), fNbins[0], fVmin[0], fVmax[0]);
576 } else {
577 hist = new TH1D(hname, htitle.Data(), fNbins[0], fVmin[0], fVmax[0]);
578 }
587 if (canExtend) hist->SetCanExtend(TH1::kAllAxes);
588 if (!hkeep) {
589 hist->GetXaxis()->SetTitle(fVar[0]->GetTitle());
590 hist->SetBit(kCanDelete);
591 if (!opt.Contains("goff")) hist->SetDirectory(nullptr);
592 }
593 if (opt.Length() && opt.Contains("e")) hist->Sumw2();
594 }
595 fVar[0]->SetAxis(hist->GetXaxis());
596 fObject = hist;
597
598 // 2-D distribution
599 } else if (fDimension == 2 && !(optpara || optcandle)) {
600 fAction = 2;
601 if (!fOldHistogram || !optSame) {
602 fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y", 40);
603 fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x", 40);
604 if (opt.Contains("prof")) fNbins[1] = gEnv->GetValue("Hist.Binning.2D.Prof", 100);
605 if (optSame) {
606 TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
607 if (oldhtemp) {
608 fNbins[1] = oldhtemp->GetXaxis()->GetNbins();
609 fVmin[1] = oldhtemp->GetXaxis()->GetXmin();
610 fVmax[1] = oldhtemp->GetXaxis()->GetXmax();
611 fNbins[0] = oldhtemp->GetYaxis()->GetNbins();
612 fVmin[0] = oldhtemp->GetYaxis()->GetXmin();
613 fVmax[0] = oldhtemp->GetYaxis()->GetXmax();
614 } else {
615 fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x", 40);
616 fVmin[1] = gPad->GetUxmin();
617 fVmax[1] = gPad->GetUxmax();
618 fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y", 40);
619 fVmin[0] = gPad->GetUymin();
620 fVmax[0] = gPad->GetUymax();
621 }
622 } else {
623 if (!fOldHistogram) fAction = -2;
624 fVmin[1] = xmin;
625 fVmax[1] = xmax;
626 fVmin[0] = ymin;
627 fVmax[0] = ymax;
628 if (xmin < xmax && ymin < ymax) canExtend = false;
629 }
630 }
631 if (profile || opt.Contains("prof")) {
632 TProfile *hp;
633 if (fOldHistogram) {
634 fAction = 4;
635 hp = (TProfile*)fOldHistogram;
636 } else {
637 if (fAction < 0) {
638 fAction = -4;
639 fVmin[1] = xmin;
640 fVmax[1] = xmax;
641 if (xmin < xmax) canExtend = false;
642 }
643 if (fAction == 2) {
644 //we come here when option = "same prof"
645 fAction = -4;
646 TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
647 if (oldhtemp) {
648 fNbins[1] = oldhtemp->GetXaxis()->GetNbins();
649 fVmin[1] = oldhtemp->GetXaxis()->GetXmin();
650 fVmax[1] = oldhtemp->GetXaxis()->GetXmax();
651 }
652 }
653 if (opt.Contains("profs")) {
654 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "s");
655 } else if (opt.Contains("profi")) {
656 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "i");
657 } else if (opt.Contains("profg")) {
658 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "g");
659 } else {
660 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "");
661 }
662 if (!hkeep) {
663 hp->SetBit(kCanDelete);
664 if (!opt.Contains("goff")) hp->SetDirectory(nullptr);
665 }
674 if (canExtend) hp->SetCanExtend(TH1::kAllAxes);
675 }
676 fVar[1]->SetAxis(hp->GetXaxis());
677 fObject = hp;
678
679 } else {
680 TH2 *h2;
681 if (fOldHistogram) {
682 h2 = (TH2F*)fOldHistogram;
683 } else {
684 TString precision = gEnv->GetValue("Hist.Precision.2D", "float");
685 if (precision.Contains("float")) {
686 h2 = new TH2F(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
687 } else {
688 h2 = new TH2D(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
689 }
690 h2->SetLineColor(fTree->GetLineColor());
691 h2->SetLineWidth(fTree->GetLineWidth());
692 h2->SetLineStyle(fTree->GetLineStyle());
693 h2->SetFillColor(fTree->GetFillColor());
694 h2->SetFillStyle(fTree->GetFillStyle());
695 h2->SetMarkerStyle(fTree->GetMarkerStyle());
696 h2->SetMarkerColor(fTree->GetMarkerColor());
697 h2->SetMarkerSize(fTree->GetMarkerSize());
698 if (canExtend) h2->SetCanExtend(TH1::kAllAxes);
699 if (!hkeep) {
700 h2->GetXaxis()->SetTitle(fVar[1]->GetTitle());
701 h2->GetYaxis()->SetTitle(fVar[0]->GetTitle());
702 h2->SetBit(TH1::kNoStats);
703 h2->SetBit(kCanDelete);
704 if (!opt.Contains("goff")) h2->SetDirectory(nullptr);
705 }
706 }
707 fVar[0]->SetAxis(h2->GetYaxis());
708 fVar[1]->SetAxis(h2->GetXaxis());
709 bool graph = false;
710 Int_t l = opt.Length();
711 if (l == 0 || optSame) graph = true;
712 if (opt.Contains("p") || opt.Contains("*") || opt.Contains("l")) graph = true;
713 if (opt.Contains("surf") || opt.Contains("lego") || opt.Contains("cont")) graph = false;
714 if (opt.Contains("col") || opt.Contains("hist") || opt.Contains("scat")) graph = false;
715 if (opt.Contains("box")) graph = false;
716 fObject = h2;
717 if (graph) {
718 fAction = 12;
719 if (!fOldHistogram && !optSame) fAction = -12;
720 }
721 }
722
723 // 3-D distribution
724 } else if ((fDimension == 3 || fDimension == 4) && !(optpara || optcandle)) {
725 fAction = 3;
726 if (fDimension == 4) fAction = 40;
727 if (!fOldHistogram || !optSame) {
728 fNbins[0] = gEnv->GetValue("Hist.Binning.3D.z", 20);
729 fNbins[1] = gEnv->GetValue("Hist.Binning.3D.y", 20);
730 fNbins[2] = gEnv->GetValue("Hist.Binning.3D.x", 20);
731 if (fDimension == 3 && opt.Contains("prof")) {
732 fNbins[1] = gEnv->GetValue("Hist.Binning.3D.Profy", 20);
733 fNbins[2] = gEnv->GetValue("Hist.Binning.3D.Profx", 20);
734 } else if (fDimension == 3 && opt.Contains("col")) {
735 fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y", 40);
736 fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x", 40);
737 }
738 if (optSame) {
739 TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
740 if (oldhtemp) {
741 fNbins[2] = oldhtemp->GetXaxis()->GetNbins();
742 fVmin[2] = oldhtemp->GetXaxis()->GetXmin();
743 fVmax[2] = oldhtemp->GetXaxis()->GetXmax();
744 fNbins[1] = oldhtemp->GetYaxis()->GetNbins();
745 fVmin[1] = oldhtemp->GetYaxis()->GetXmin();
746 fVmax[1] = oldhtemp->GetYaxis()->GetXmax();
747 fNbins[0] = oldhtemp->GetZaxis()->GetNbins();
748 fVmin[0] = oldhtemp->GetZaxis()->GetXmin();
749 fVmax[0] = oldhtemp->GetZaxis()->GetXmax();
750 } else {
751 TView *view = gPad->GetView();
752 if (!view) {
753 Error("Begin", "You cannot use option same when no 3D view exists");
754 fVmin[0] = fVmin[1] = fVmin[2] = -1;
755 fVmax[0] = fVmax[1] = fVmax[2] = 1;
756 view = TView::CreateView(1, fVmin, fVmax);
757 }
758 Double_t *rmin = view->GetRmin();
759 Double_t *rmax = view->GetRmax();
760 fNbins[2] = gEnv->GetValue("Hist.Binning.3D.z", 20);
761 fVmin[2] = rmin[0];
762 fVmax[2] = rmax[0];
763 fNbins[1] = gEnv->GetValue("Hist.Binning.3D.y", 20);
764 fVmin[1] = rmin[1];
765 fVmax[1] = rmax[1];
766 fNbins[0] = gEnv->GetValue("Hist.Binning.3D.x", 20);
767 fVmin[0] = rmin[2];
768 fVmax[0] = rmax[2];
769 }
770 } else {
771 if (!fOldHistogram && fDimension == 3) fAction = -3;
772 fVmin[2] = xmin;
773 fVmax[2] = xmax;
774 fVmin[1] = ymin;
775 fVmax[1] = ymax;
776 fVmin[0] = zmin;
777 fVmax[0] = zmax;
778 if (xmin < xmax && ymin < ymax && zmin < zmax) canExtend = false;
779 }
780 }
781 if ((fDimension == 3) && (profile || opt.Contains("prof"))) {
782 TProfile2D *hp;
783 if (fOldHistogram) {
784 fAction = 23;
786 } else {
787 if (fAction < 0) {
788 fAction = -23;
789 fVmin[2] = xmin;
790 fVmax[2] = xmax;
791 fVmin[1] = ymin;
792 fVmax[1] = ymax;
793 if (xmin < xmax && ymin < ymax) canExtend = false;
794 }
795 if (opt.Contains("profs")) {
796 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "s");
797 } else if (opt.Contains("profi")) {
798 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "i");
799 } else if (opt.Contains("profg")) {
800 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "g");
801 } else {
802 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "");
803 }
804 if (!hkeep) {
805 hp->SetBit(kCanDelete);
806 if (!opt.Contains("goff")) hp->SetDirectory(nullptr);
807 }
816 if (canExtend) hp->SetCanExtend(TH1::kAllAxes);
817 }
818 fVar[1]->SetAxis(hp->GetYaxis());
819 fVar[2]->SetAxis(hp->GetXaxis());
820 fObject = hp;
821 } else if (fDimension == 3 && opt.Contains("col")) {
822 TH2F *h2;
823 if (fOldHistogram) {
824 h2 = (TH2F*)fOldHistogram;
825 } else {
826 h2 = new TH2F(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
827 h2->SetLineColor(fTree->GetLineColor());
828 h2->SetLineWidth(fTree->GetLineWidth());
829 h2->SetLineStyle(fTree->GetLineStyle());
830 h2->SetFillColor(fTree->GetFillColor());
831 h2->SetFillStyle(fTree->GetFillStyle());
832 h2->SetMarkerStyle(fTree->GetMarkerStyle());
833 h2->SetMarkerColor(fTree->GetMarkerColor());
834 h2->SetMarkerSize(fTree->GetMarkerSize());
835 if (canExtend) h2->SetCanExtend(TH1::kAllAxes);
836 if (!hkeep) {
837 h2->GetXaxis()->SetTitle(fVar[1]->GetTitle());
838 h2->GetYaxis()->SetTitle(fVar[0]->GetTitle());
839 h2->GetZaxis()->SetTitle(fVar[2]->GetTitle());
840 h2->SetBit(TH1::kNoStats);
841 h2->SetBit(kCanDelete);
842 if (!opt.Contains("goff")) h2->SetDirectory(nullptr);
843 }
844 }
845 fVar[0]->SetAxis(h2->GetYaxis());
846 fVar[1]->SetAxis(h2->GetXaxis());
847 fObject = h2;
848 fAction = 33;
849 } else {
850 TH3 *h3;
851 if (fOldHistogram) {
852 h3 = (TH3F*)fOldHistogram;
853 } else {
854 TString precision = gEnv->GetValue("Hist.Precision.3D", "float");
855 if (precision.Contains("float")) {
856 h3 = new TH3F(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
857 } else {
858 h3 = new TH3D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
859 }
860 h3->SetLineColor(fTree->GetLineColor());
861 h3->SetLineWidth(fTree->GetLineWidth());
862 h3->SetLineStyle(fTree->GetLineStyle());
863 h3->SetFillColor(fTree->GetFillColor());
864 h3->SetFillStyle(fTree->GetFillStyle());
865 h3->SetMarkerStyle(fTree->GetMarkerStyle());
866 h3->SetMarkerColor(fTree->GetMarkerColor());
867 h3->SetMarkerSize(fTree->GetMarkerSize());
868 if (canExtend) h3->SetCanExtend(TH1::kAllAxes);
869 if (!hkeep) {
870 //small correction for the title offsets in x,y to take into account the angles
871 Double_t xoffset = h3->GetXaxis()->GetTitleOffset();
872 Double_t yoffset = h3->GetYaxis()->GetTitleOffset();
873 h3->GetXaxis()->SetTitleOffset(1.2 * xoffset);
874 h3->GetYaxis()->SetTitleOffset(1.2 * yoffset);
875 h3->GetXaxis()->SetTitle(fVar[2]->GetTitle());
876 h3->GetYaxis()->SetTitle(fVar[1]->GetTitle());
877 h3->GetZaxis()->SetTitle(fVar[0]->GetTitle());
878 h3->SetBit(kCanDelete);
879 h3->SetBit(TH1::kNoStats);
880 if (!opt.Contains("goff")) h3->SetDirectory(nullptr);
881 }
882 }
883 fVar[0]->SetAxis(h3->GetZaxis());
884 fVar[1]->SetAxis(h3->GetYaxis());
885 fVar[2]->SetAxis(h3->GetXaxis());
886 fObject = h3;
887 Int_t noscat = strlen(option);
888 if (optSame) noscat -= 4;
889 if (!noscat && fDimension == 3) {
890 fAction = 13;
891 if (!fOldHistogram && !optSame) fAction = -13;
892 }
893 }
894 // An Event List
895 } else if (enlist) {
896 fAction = 5;
898 fTree->SetEstimate(1);
899 fObject = enlist;
900 } else if (evlist) {
901 fAction = 5;
903 fTree->SetEstimate(1);
904 fObject = evlist;
905 } else if (optcandle || optpara || opt5d) {
906 if (optcandle) fAction = 7;
907 else if (opt5d) fAction = 8;
908 else fAction = 6;
909 }
910 if (varexp) delete[] varexp;
911 for (i = 0; i < fValSize; ++i)
912 fVarMultiple[i] = false;
913 fSelectMultiple = false;
914 for (i = 0; i < fDimension; ++i) {
915 if (fVar[i] && fVar[i]->GetMultiplicity()) fVarMultiple[i] = true;
916 }
917
919
922 fNfill = 0;
923
924 for (i = 0; i < fDimension; ++i) {
925 if (!fVal[i] && fVar[i]) {
926 fVal[i] = new Double_t[(Int_t)fTree->GetEstimate()];
927 }
928 }
929
930 if (!fW) fW = new Double_t[(Int_t)fTree->GetEstimate()];
931
932 for (i = 0; i < fValSize; ++i) {
933 fVmin[i] = DBL_MAX;
934 fVmax[i] = -DBL_MAX;
935 }
936}
937
938////////////////////////////////////////////////////////////////////////////////
939/// Delete internal buffers.
940
942{
944 for (Int_t i = 0; i < fValSize; ++i) {
945 delete fVar[i];
946 fVar[i] = nullptr;
947 }
948 delete fSelect; fSelect = nullptr;
949 fManager = nullptr;
950 fMultiplicity = 0;
951}
952
953////////////////////////////////////////////////////////////////////////////////
954/// Compile input variables and selection expression.
955///
956/// varexp is an expression of the general form e1:e2:e3
957/// where e1,etc is a formula referencing a combination of the columns
958///
959/// Example:
960///
961/// varexp = x simplest case: draw a 1-Dim distribution of column named x
962/// = sqrt(x) : draw distribution of sqrt(x)
963/// = x*y/z
964/// = y:sqrt(x) 2-Dim distribution of y versus sqrt(x)
965///
966/// selection is an expression with a combination of the columns
967///
968/// Example:
969///
970/// selection = "x<y && sqrt(z)>3.2"
971///
972/// in a selection all the C++ operators are authorized
973///
974/// Return false if any of the variable is not compilable.
975
976bool TSelectorDraw::CompileVariables(const char *varexp, const char *selection)
977{
978 Int_t i, nch, ncols;
979
980 // Compile selection expression if there is one
981 fDimension = 0;
982 ClearFormula();
983 fMultiplicity = 0;
984 fObjEval = false;
985
986 if (strlen(selection)) {
987 fSelect = new TTreeFormula("Selection", selection, fTree);
988 fSelect->SetQuickLoad(true);
989 if (!fSelect->GetNdim()) {
990 delete fSelect;
991 fSelect = nullptr;
992 return false;
993 }
994 }
995
996 // if varexp is empty, take first column by default
997 nch = strlen(varexp);
998 if (nch == 0) {
999 fDimension = 0;
1000 if (fSelect) {
1002 }
1004
1005 if (fManager) {
1006 fManager->Sync();
1007
1010 }
1011
1012 return true;
1013 }
1014
1015 // otherwise select only the specified columns
1016 std::vector<TString> varnames;
1017 ncols = SplitNames(varexp, varnames);
1018
1019 InitArrays(ncols);
1020
1022 if (fSelect) fManager->Add(fSelect);
1024 for (i = 0; i < ncols; ++i) {
1025 fVar[i] = new TTreeFormula(TString::Format("Var%i", i + 1), varnames[i].Data(), fTree);
1026 fVar[i]->SetQuickLoad(true);
1027 if(!fVar[i]->GetNdim()) { ClearFormula(); return false; }
1028 fManager->Add(fVar[i]);
1029 }
1030 fManager->Sync();
1031
1034
1035 fDimension = ncols;
1036
1037 if (ncols == 1) {
1038 TClass *cl = fVar[0]->EvalClass();
1039 if (cl) {
1040 fObjEval = true;
1041 }
1042 }
1043 return true;
1044}
1045
1046////////////////////////////////////////////////////////////////////////////////
1047/// Return the last values corresponding to the i-th component
1048/// of the formula being processed (where the component are ':' separated).
1049/// The actual number of entries is:
1050///
1051/// GetSelectedRows() % tree->GetEstimate()
1052///
1053/// Note GetSelectedRows currently returns the actual number of values plotted
1054/// and thus if the formula contains arrays, this number might be greater than
1055/// the number of entries in the trees.
1056///
1057/// By default TTree::Draw creates the arrays obtained
1058/// with all GetVal and GetW with a length corresponding to the
1059/// parameter fEstimate. By default fEstimate=10000 and can be modified
1060/// via TTree::SetEstimate. A possible recipe is to do
1061///
1062/// tree->SetEstimate(tree->GetEntries());
1063///
1064/// You must call SetEstimate if the expected number of selected rows
1065/// is greater than 10000.
1066///
1067/// See TTree::Draw for additional details.
1068
1070{
1071 if (i < 0 || i >= fDimension)
1072 return nullptr;
1073 else
1074 return fVal[i];
1075}
1076
1077////////////////////////////////////////////////////////////////////////////////
1078/// Return the TTreeFormula corresponding to the i-th component
1079/// of the request formula (where the component are ':' separated).
1080
1082{
1083 if (i < 0 || i >= fDimension)
1084 return nullptr;
1085 else
1086 return fVar[i];
1087}
1088
1089////////////////////////////////////////////////////////////////////////////////
1090/// Initialization of the primitive type arrays if the new size is bigger than the available space.
1091
1093{
1094 if (newsize > fValSize) {
1095 Int_t oldsize = fValSize;
1096 while (fValSize < newsize)
1097 fValSize *= 2; // Double the available space until it matches the new size.
1098 if (fNbins) delete [] fNbins;
1099 if (fVmin) delete [] fVmin;
1100 if (fVmax) delete [] fVmax;
1101 if (fVarMultiple) delete [] fVarMultiple;
1102
1103 fNbins = new Int_t[fValSize];
1104 fVmin = new Double_t[fValSize];
1105 fVmax = new Double_t[fValSize];
1106 fVarMultiple = new bool[fValSize];
1107
1108 for (Int_t i = 0; i < oldsize; ++i)
1109 delete [] fVal[i];
1110 delete [] fVal;
1111 delete [] fVar;
1112 fVal = new Double_t*[fValSize];
1113 fVar = new TTreeFormula*[fValSize];
1114 for (Int_t i = 0; i < fValSize; ++i) {
1115 fVal[i] = nullptr;
1116 fVar[i] = nullptr;
1117 }
1118 }
1119}
1120
1121////////////////////////////////////////////////////////////////////////////////
1122/// Build Index array for names in varexp.
1123/// This will allocated a C style array of TString and Ints
1124
1125UInt_t TSelectorDraw::SplitNames(const TString &varexp, std::vector<TString> &names)
1126{
1127 names.clear();
1128
1129 bool ternary = false;
1130 Int_t prev = 0;
1131 for (Int_t i = 0; i < varexp.Length(); i++) {
1132 if (varexp[i] == ':'
1133 && !((i > 0 && varexp[i-1] == ':') || varexp[i+1] == ':')
1134 ) {
1135 if (ternary) {
1136 ternary = false;
1137 } else {
1138 names.push_back(varexp(prev, i - prev));
1139 prev = i + 1;
1140 }
1141 }
1142 if (varexp[i] == '?') {
1143 ternary = true;
1144 }
1145 }
1146 names.push_back(varexp(prev, varexp.Length() - prev));
1147 return names.size();
1148}
1149
1150
1151////////////////////////////////////////////////////////////////////////////////
1152/// This function is called at the first entry of a new tree in a chain.
1153
1155{
1156 if (fTree) fWeight = fTree->GetWeight();
1157 if (fVar) {
1158 for (Int_t i = 0; i < fDimension; ++i) {
1159 if (fVar[i]) fVar[i]->UpdateFormulaLeaves();
1160 }
1161 }
1163 return true;
1164}
1165
1166////////////////////////////////////////////////////////////////////////////////
1167/// Called in the entry loop for all entries accepted by Select.
1168
1170{
1171 if (fObjEval) {
1172 ProcessFillObject(entry);
1173 return;
1174 }
1175
1176 if (fMultiplicity) {
1177 ProcessFillMultiple(entry);
1178 return;
1179 }
1180
1181 if (fNfill >= fTree->GetEstimate())
1182 fNfill = 0;
1183
1184 // simple case with no multiplicity
1185 if (fForceRead && fManager->GetNdata() <= 0) return;
1186
1187 if (fSelect) {
1189 if (!fW[fNfill]) return;
1190 } else fW[fNfill] = fWeight;
1191 if (fVal) {
1192 for (Int_t i = 0; i < fDimension; ++i) {
1193 if (fVar[i]) fVal[i][fNfill] = fVar[i]->EvalInstance(0);
1194 }
1195 }
1196 fNfill++;
1197 if (fNfill >= fTree->GetEstimate()) {
1198 TakeAction();
1199 }
1200}
1201
1202////////////////////////////////////////////////////////////////////////////////
1203/// Called in the entry loop for all entries accepted by Select.
1204/// Complex case with multiplicity.
1205
1207{
1208 if (fNfill >= fTree->GetEstimate())
1209 fNfill = 0;
1210
1211 // Grab the array size of the formulas for this entry
1212 Int_t ndata = fManager->GetNdata();
1213
1214 // No data at all, let's move on to the next entry.
1215 if (!ndata) return;
1216
1217 // If the entry list is a TEntryListArray, get the selected subentries for this entry
1218 TEntryList *subList = nullptr;
1219 if (fTreeElistArray) {
1220 subList = fTreeElistArray->GetSubListForEntry(entry, fTree->GetTree());
1221 }
1222
1223 Int_t nfill0 = fNfill;
1224
1225 // Calculate the first values
1226 if (fSelect) {
1227 // coverity[var_deref_model] fSelectMultiple==true => fSelect != 0
1229 if (!fW[fNfill] && !fSelectMultiple) return;
1230 } else fW[fNfill] = fWeight;
1231
1232 // Always call EvalInstance(0) to insure the loading
1233 // of the branches.
1234 if (fW[fNfill] && (!subList || subList->Contains(0))) {
1235 if (fDimension == 0 && fSelectMultiple) fCurrentSubEntry = (Long64_t) 0; // to fill TEntryListArray
1236 for (Int_t i = 0; i < fDimension; ++i) {
1237 if (fVar[i]) fVal[i][fNfill] = fVar[i]->EvalInstance(0);
1238 }
1239 fNfill++;
1240 if (fNfill >= fTree->GetEstimate()) {
1241 TakeAction();
1242 }
1243 } else {
1244 for (Int_t i = 0; i < fDimension; ++i) {
1245 if (fVar[i]) fVar[i]->ResetLoading();
1246 }
1247 }
1248 Double_t ww = fW[nfill0];
1249
1250 for (Int_t i = 1; i < ndata; i++) {
1251 if (fNfill >= fTree->GetEstimate())
1252 fNfill = 0;
1253 if (subList && !subList->Contains(i))
1254 continue;
1255 if (fSelectMultiple) {
1256 // coverity[var_deref_model] fSelectMultiple==true => fSelect != 0
1257 ww = fWeight * fSelect->EvalInstance(i);
1258 if (ww == 0) continue;
1259 if (fNfill == nfill0) {
1260 for (Int_t k = 0; k < fDimension; ++k) {
1261 if (!fVarMultiple[k]) fVal[k][fNfill] = fVar[k]->EvalInstance(0);
1262 }
1263 }
1264 if (fDimension == 0) fCurrentSubEntry = (Long64_t) i; // to fill TEntryListArray
1265 }
1266 for (Int_t k = 0; k < fDimension; ++k) {
1267 if (fVarMultiple[k]) fVal[k][fNfill] = fVar[k]->EvalInstance(i);
1268 else fVal[k][fNfill] = fVal[k][nfill0];
1269 }
1270 fW[fNfill] = ww;
1271
1272 fNfill++;
1273 if (fNfill >= fTree->GetEstimate()) {
1274 TakeAction();
1275 }
1276 }
1277}
1278
1279////////////////////////////////////////////////////////////////////////////////
1280/// Called in the entry loop for all entries accepted by Select.
1281/// Case where the only variable returns an object (or pointer to).
1282
1284{
1285 // Complex case with multiplicity.
1286
1287 if (fNfill >= fTree->GetEstimate())
1288 fNfill = 0;
1289
1290 // Grab the array size of the formulas for this entry
1291 Int_t ndata = fManager->GetNdata();
1292
1293 // No data at all, let's move on to the next entry.
1294 if (!ndata) return;
1295
1296 Int_t nfill0 = fNfill;
1297 Double_t ww = 0;
1298
1299 for (Int_t i = 0; i < ndata; i++) {
1300 if (i == 0) {
1301 if (fSelect) {
1303 if (!fW[fNfill] && !fSelectMultiple) return;
1304 } else fW[fNfill] = fWeight;
1305 ww = fW[nfill0];
1306 } else if (fSelectMultiple) {
1307 ww = fWeight * fSelect->EvalInstance(i);
1308 if (ww == 0) continue;
1309 }
1310 if (fDimension >= 1 && fVar[0]) {
1311 TClass *cl = fVar[0]->EvalClass();
1312 if (cl == TBits::Class()) {
1313
1314 void *obj = fVar[0]->EvalObject(i);
1315
1316 if (obj) {
1317 TBits *bits = (TBits*)obj;
1318 Int_t nbits = bits->GetNbits();
1319
1320 Int_t nextbit = -1;
1321 while (true) {
1322 nextbit = bits->FirstSetBit(nextbit + 1);
1323 if (nextbit >= nbits) break;
1324 fVal[0][fNfill] = nextbit;
1325 fW[fNfill] = ww;
1326 fNfill++;
1327 }
1328 }
1329
1330 } else {
1331
1332 if (!TestBit(kWarn)) {
1333 Warning("ProcessFillObject",
1334 "Not implemented for %s",
1335 cl ? cl->GetName() : "unknown class");
1336 SetBit(kWarn);
1337 }
1338
1339 }
1340 }
1341 if (fNfill >= fTree->GetEstimate()) {
1342 TakeAction();
1343 }
1344 }
1345
1346}
1347
1348////////////////////////////////////////////////////////////////////////////////
1349/// Set number of entries to estimate variable limits.
1350
1352{
1353 if (fVal) {
1354 for (Int_t i = 0; i < fValSize; ++i) {
1355 delete [] fVal[i];
1356 fVal[i] = nullptr;
1357 }
1358 }
1359 delete [] fW; fW = nullptr;
1360}
1361
1362////////////////////////////////////////////////////////////////////////////////
1363/// Execute action for object obj fNfill times.
1364
1366{
1367 Int_t i;
1368 //__________________________1-D histogram_______________________
1369 if (fAction == 1)((TH1*)fObject)->FillN(fNfill, fVal[0], fW);
1370 //__________________________2-D histogram_______________________
1371 else if (fAction == 2) {
1372 TH2 *h2 = (TH2*)fObject;
1373 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1374 }
1375 //__________________________Profile histogram_______________________
1376 else if (fAction == 4)((TProfile*)fObject)->FillN(fNfill, fVal[1], fVal[0], fW);
1377 //__________________________Event List______________________________
1378 else if (fAction == 5) {
1380 TEntryListArray *enlistarray = (TEntryListArray*)fObject;
1381 Long64_t enumb = fTree->GetTree()->GetReadEntry();
1382 enlistarray->Enter(enumb, nullptr, fCurrentSubEntry);
1383 } else if (fObject->InheritsFrom(TEntryList::Class())) {
1384 TEntryList *enlist = (TEntryList*)fObject;
1385 Long64_t enumb = fTree->GetTree()->GetReadEntry();
1386 enlist->Enter(enumb);
1387 } else {
1388 TEventList *evlist = (TEventList*)fObject;
1390 if (evlist->GetIndex(enumb) < 0) evlist->Enter(enumb);
1391 }
1392 }
1393 //__________________________2D scatter plot_______________________
1394 else if (fAction == 12) {
1395 TH2 *h2 = (TH2*)fObject;
1396 if (h2->CanExtendAllAxes() && h2->TestBit(kCanDelete)) {
1397 for (i = 0; i < fNfill; i++) {
1398 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1399 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1400 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1401 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1402 }
1404 }
1405 TGraph *pm = new TGraph(fNfill, fVal[1], fVal[0]);
1406 pm->SetEditable(false);
1407 pm->SetBit(kCanDelete);
1416
1417 if (!fDraw && !strstr(fOption.Data(), "goff")) {
1418 if (fOption.Length() == 0 || strcasecmp(fOption.Data(), "same") == 0) pm->Draw("p");
1419 else pm->Draw(fOption.Data());
1420 }
1421 if (!h2->TestBit(kCanDelete)) {
1422 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1423 }
1424 }
1425 //__________________________3D scatter plot_______________________
1426 else if (fAction == 3) {
1427 TH3 *h3 = (TH3*)fObject;
1428 if (!h3->TestBit(kCanDelete)) {
1429 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1430 }
1431 } else if (fAction == 13) {
1432 TPolyMarker3D *pm3d = new TPolyMarker3D(fNfill);
1436 for (i = 0; i < fNfill; i++) {
1437 pm3d->SetPoint(i, fVal[2][i], fVal[1][i], fVal[0][i]);
1438 }
1439 pm3d->Draw();
1440 TH3 *h3 = (TH3*)fObject;
1441 if (!h3->TestBit(kCanDelete)) {
1442 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1443 }
1444 }
1445 //__________________________3D scatter plot (3rd variable = col)__
1446 else if (fAction == 33) {
1447 TH2 *h2 = (TH2*)fObject;
1448 TakeEstimate();
1449 Int_t ncolors = gStyle->GetNumberOfColors();
1450 TObjArray *grs = (TObjArray*)h2->GetListOfFunctions()->FindObject("graphs");
1451 Int_t col;
1452 TGraph *gr;
1453 if (!grs) {
1454 grs = new TObjArray(ncolors);
1455 grs->SetOwner();
1456 grs->SetName("graphs");
1457 h2->GetListOfFunctions()->Add(grs, "P");
1458 for (col = 0; col < ncolors; col++) {
1459 gr = new TGraph();
1463 grs->AddAt(gr, col);
1464 }
1465 }
1466 h2->SetEntries(fNfill);
1467 h2->SetMinimum(fVmin[2]);
1468 h2->SetMaximum(fVmax[2]);
1469 // Fill the graphs according to the color
1470 for (i = 0; i < fNfill; i++) {
1471 col = Int_t(ncolors * ((fVal[2][i] - fVmin[2]) / (fVmax[2] - fVmin[2])));
1472 if (col < 0) col = 0;
1473 if (col > ncolors - 1) col = ncolors - 1;
1474 gr = (TGraph*)grs->UncheckedAt(col);
1475 if (gr) gr->SetPoint(gr->GetN(), fVal[1][i], fVal[0][i]);
1476 }
1477 // Remove potential empty graphs
1478 for (col = 0; col < ncolors; col++) {
1479 gr = (TGraph*)grs->At(col);
1480 if (gr && gr->GetN() <= 0) grs->Remove(gr);
1481 }
1482 }
1483 //__________________________2D Profile Histogram__________________
1484 else if (fAction == 23) {
1485 TProfile2D *hp2 = (TProfile2D*)fObject;
1486 for (i = 0; i < fNfill; i++) hp2->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1487 }
1488 //__________________________4D scatter plot_______________________
1489 else if (fAction == 40) {
1490 TakeEstimate();
1491 TH3 *h3 = (TH3*)fObject;
1492 Int_t ncolors = gStyle->GetNumberOfColors();
1493 if (ncolors == 0) {
1495 ncolors = gStyle->GetNumberOfColors();
1496 }
1497 TObjArray *pms = (TObjArray*)h3->GetListOfFunctions()->FindObject("polymarkers");
1498 Int_t col;
1499 TPolyMarker3D *pm3d;
1500 if (!pms) {
1501 pms = new TObjArray(ncolors);
1502 pms->SetOwner();
1503 pms->SetName("polymarkers");
1504 h3->GetListOfFunctions()->Add(pms);
1505 for (col = 0; col < ncolors; col++) {
1506 pm3d = new TPolyMarker3D();
1510 pms->AddAt(pm3d, col);
1511 }
1512 }
1513 h3->SetEntries(fNfill);
1514 h3->SetMinimum(fVmin[3]);
1515 h3->SetMaximum(fVmax[3]);
1516 for (i = 0; i < fNfill; i++) {
1517 col = Int_t(ncolors * ((fVal[3][i] - fVmin[3]) / (fVmax[3] - fVmin[3])));
1518 if (col > ncolors-1) col = ncolors-1;
1519 if (col < 0) col = 0;
1520 pm3d = (TPolyMarker3D*)pms->UncheckedAt(col);
1521 pm3d->SetPoint(pm3d->GetLastPoint() + 1, fVal[2][i], fVal[1][i], fVal[0][i]);
1522 }
1523 }
1524 //__________________________Parallel coordinates / candle chart_______________________
1525 else if (fAction == 6 || fAction == 7) {
1526 TakeEstimate();
1527 bool candle = (fAction == 7);
1528 // Using CINT to avoid a dependency in TParallelCoord
1529 if (!fOption.Contains("goff"))
1530 gROOT->ProcessLine(TString::Format("TParallelCoord::BuildParallelCoord((TSelectorDraw*)0x%zx,0x%zx)",
1531 (size_t)this, (size_t)candle));
1532 } else if (fAction == 8) {
1533 //gROOT->ProcessLineFast(TString::Format("(new TGL5DDataSet((TTree *)0x%1x))->Draw(\"%s\");", fTree, fOption.Data()));
1534 }
1535 //__________________________something else_______________________
1536 else if (fAction < 0) {
1537 fAction = -fAction;
1538 TakeEstimate();
1539 }
1540
1541 // Do we need to update screen?
1543 if (!fTree->GetUpdate()) return;
1544 if (fSelectedRows > fDraw + fTree->GetUpdate()) {
1545 if (fDraw) gPad->Modified();
1546 else fObject->Draw(fOption.Data());
1547 gPad->Update();
1549 }
1550}
1551
1552////////////////////////////////////////////////////////////////////////////////
1553/// Estimate limits for 1-D, 2-D or 3-D objects.
1554
1556{
1557 Int_t i;
1558 Double_t rmin[3], rmax[3];
1559 Double_t vminOld[4], vmaxOld[4];
1560 for (i = 0; i < fValSize && i < 4; i++) {
1561 vminOld[i] = fVmin[i];
1562 vmaxOld[i] = fVmax[i];
1563 }
1564 for (i = 0; i < fValSize; ++i) {
1565 fVmin[i] = DBL_MAX;
1566 fVmax[i] = - DBL_MAX;
1567 }
1568 //__________________________1-D histogram_______________________
1569 if (fAction == 1) {
1570 TH1 *h1 = (TH1*)fObject;
1571 if (h1->CanExtendAllAxes()) {
1572 for (i = 0; i < fNfill; i++) {
1573 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1574 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1575 }
1577 }
1578 h1->FillN(fNfill, fVal[0], fW);
1579 //__________________________2-D histogram_______________________
1580 } else if (fAction == 2) {
1581 TH2 *h2 = (TH2*)fObject;
1582 if (h2->CanExtendAllAxes()) {
1583 for (i = 0; i < fNfill; i++) {
1584 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1585 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1586 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1587 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1588 }
1590 }
1591 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1592 //__________________________Profile histogram_______________________
1593 } else if (fAction == 4) {
1594 TProfile *hp = (TProfile*)fObject;
1595 if (hp->CanExtendAllAxes()) {
1596 for (i = 0; i < fNfill; i++) {
1597 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1598 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1599 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1600 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1601 }
1603 }
1604 hp->FillN(fNfill, fVal[1], fVal[0], fW);
1605 //__________________________2D scatter plot_______________________
1606 } else if (fAction == 12) {
1607 TH2 *h2 = (TH2*)fObject;
1608 if (h2->CanExtendAllAxes()) {
1609 for (i = 0; i < fNfill; i++) {
1610 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1611 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1612 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1613 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1614 }
1616 // In case the new lower limits of h2 axis are 0, it is better to set them to the minimum of
1617 // the data set (which should be >0) to avoid data cut when plotting in log scale.
1618 TAxis *aX = h2->GetXaxis();
1619 TAxis *aY = h2->GetYaxis();
1620 Double_t xmin = aX->GetXmin();
1621 Double_t ymin = aY->GetXmin();
1622 if (xmin == 0 || ymin == 0) {
1623 if (aX->GetBinUpEdge(aX->FindFixBin(0.01*aX->GetBinWidth(aX->GetFirst()))) > fVmin[1]) xmin = fVmin[1];
1624 if (aY->GetBinUpEdge(aY->FindFixBin(0.01*aY->GetBinWidth(aY->GetFirst()))) > fVmin[0]) ymin = fVmin[0];
1625 h2->SetBins(aX->GetNbins(), xmin, aX->GetXmax(), aY->GetNbins(), ymin, aY->GetXmax());
1626 }
1627 }
1628
1629 if (!strstr(fOption.Data(), "same") && !strstr(fOption.Data(), "goff")) {
1630 if (!h2->TestBit(kCanDelete)) {
1631 // case like: T.Draw("y:x>>myhist")
1632 // we must draw a copy before filling the histogram h2=myhist
1633 // because h2 will be filled below and we do not want to show
1634 // the binned scatter-plot, the TGraph being better.
1635 TH1 *h2c = h2->DrawCopy(fOption.Data(),"");
1636 if (h2c) h2c->SetStats(false);
1637 } else {
1638 // case like: T.Draw("y:x")
1639 // h2 is a temporary histogram (htemp). This histogram
1640 // will be automatically deleted by TPad::Clear
1641 h2->Draw();
1642 }
1643 gPad->Update();
1644 }
1645 TGraph *pm = new TGraph(fNfill, fVal[1], fVal[0]);
1646 pm->SetEditable(false);
1647 pm->SetBit(kCanDelete);
1656 if (!fDraw && !strstr(fOption.Data(),"goff")) {
1657 if (fOption.Length() == 0 || strcasecmp(fOption.Data(),"same")==0) {
1658 pm->Draw("p");
1659 }
1660 else {
1661 TString opt = fOption;
1662 opt.ToLower();
1663 if (opt.Contains("a")) {
1664 TString temp(opt);
1665 temp.ReplaceAll("same","");
1666 if (temp.Contains("a")) {
1667 if (h2->TestBit(kCanDelete)) {
1668 // h2 will be deleted, the axis setting is delegated to only
1669 // the TGraph.
1670 h2 = nullptr;
1671 fObject = pm->GetHistogram();
1672 }
1673 }
1674 }
1675 pm->Draw(fOption.Data());
1676 }
1677 }
1678 if (h2 && !h2->TestBit(kCanDelete)) {
1679 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1680 }
1681 //__________________________3D scatter plot with option col_______________________
1682 } else if (fAction == 33) {
1683 TH2 *h2 = (TH2*)fObject;
1684 bool process2 = false;
1685 if (h2->CanExtendAllAxes()) {
1686 if (vminOld[2] == DBL_MAX)
1687 process2 = true;
1688 for (i = 0; i < fValSize && i < 4; i++) {
1689 fVmin[i] = vminOld[i];
1690 fVmax[i] = vmaxOld[i];
1691 }
1692 for (i = 0; i < fNfill; i++) {
1693 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1694 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1695 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1696 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1697 if (process2) {
1698 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1699 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1700 }
1701 }
1703 // In case the new lower limits of h2 axis are 0, it is better to set them to the minimum of
1704 // the data set (which should be >0) to avoid data cut when plotting in log scale.
1705 TAxis *aX = h2->GetXaxis();
1706 TAxis *aY = h2->GetYaxis();
1707 Double_t xmin = aX->GetXmin();
1708 Double_t ymin = aY->GetXmin();
1709 if (xmin == 0 || ymin == 0) {
1710 if (aX->GetBinUpEdge(aX->FindFixBin(0.01*aX->GetBinWidth(aX->GetFirst()))) > fVmin[1]) xmin = fVmin[1];
1711 if (aY->GetBinUpEdge(aY->FindFixBin(0.01*aY->GetBinWidth(aY->GetFirst()))) > fVmin[0]) ymin = fVmin[0];
1712 h2->SetBins(aX->GetNbins(), xmin, aX->GetXmax(), aY->GetNbins(), ymin, aY->GetXmax());
1713 }
1714 } else {
1715 for (i = 0; i < fNfill; i++) {
1716 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1717 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1718 }
1719 }
1720 //__________________________3D scatter plot_______________________
1721 } else if (fAction == 3 || fAction == 13) {
1722 TH3 *h3 = (TH3*)fObject;
1723 if (h3->CanExtendAllAxes()) {
1724 for (i = 0; i < fNfill; i++) {
1725 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1726 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1727 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1728 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1729 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1730 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1731 }
1733 }
1734 if (fAction == 3) {
1735 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1736 return;
1737 }
1738 if (!strstr(fOption.Data(), "same") && !strstr(fOption.Data(), "goff")) {
1739 if (!h3->TestBit(kCanDelete)) {
1740 // case like: T.Draw("y:x>>myhist")
1741 // we must draw a copy before filling the histogram h3=myhist
1742 // because h3 will be filled below and we do not want to show
1743 // the binned scatter-plot, the TGraph being better.
1744 TH1 *h3c = h3->DrawCopy(fOption.Data(),"");
1745 if (h3c) h3c->SetStats(false);
1746 } else {
1747 // case like: T.Draw("y:x")
1748 // h3 is a temporary histogram (htemp). This histogram
1749 // will be automatically deleted by TPad::Clear
1750 h3->Draw(fOption.Data());
1751 }
1752 gPad->Update();
1753 } else {
1754 rmin[0] = fVmin[2]; rmin[1] = fVmin[1]; rmin[2] = fVmin[0];
1755 rmax[0] = fVmax[2]; rmax[1] = fVmax[1]; rmax[2] = fVmax[0];
1756 gPad->Clear();
1757 gPad->Range(-1, -1, 1, 1);
1758 TView::CreateView(1, rmin, rmax);
1759 }
1760 TPolyMarker3D *pm3d = new TPolyMarker3D(fNfill);
1764 for (i = 0; i < fNfill; i++) {
1765 pm3d->SetPoint(i, fVal[2][i], fVal[1][i], fVal[0][i]);
1766 }
1767 if (!fDraw && !strstr(fOption.Data(), "goff")) pm3d->Draw();
1768 if (!h3->TestBit(kCanDelete)) {
1769 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1770 }
1771
1772 //__________________________2D Profile Histogram__________________
1773 } else if (fAction == 23) {
1775 if (hp->CanExtendAllAxes()) {
1776 for (i = 0; i < fNfill; i++) {
1777 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1778 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1779 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1780 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1781 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1782 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1783 }
1785 }
1786 for (i = 0; i < fNfill; i++) hp->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1787 //__________________________4D scatter plot_______________________
1788 } else if (fAction == 40) {
1789 TH3 *h3 = (TH3*)fObject;
1790 if (h3->CanExtendAllAxes()) {
1791 for (i = 0; i < fValSize && i < 4; i++) {
1792 fVmin[i] = vminOld[i];
1793 fVmax[i] = vmaxOld[i];
1794 }
1795 for (i = 0; i < fNfill; i++) {
1796 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1797 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1798 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1799 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1800 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1801 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1802 if (fVmin[3] > fVal[3][i]) fVmin[3] = fVal[3][i];
1803 if (fVmax[3] < fVal[3][i]) fVmax[3] = fVal[3][i];
1804 }
1806 } else {
1807 for (i = 0; i < fNfill; i++) {
1808 if (fVmin[3] > fVal[3][i]) fVmin[3] = fVal[3][i];
1809 if (fVmax[3] < fVal[3][i]) fVmax[3] = fVal[3][i];
1810 }
1811 }
1812 }
1813 //__________________________Parallel coordinates plot / candle chart_______________________
1814 else if (fAction == 6 || fAction == 7) {
1815 for (i = 0; i < fDimension; ++i) {
1816 for (Long64_t entry = 0; entry < fNfill; entry++) {
1817 if (fVmin[i] > fVal[i][entry]) fVmin[i] = fVal[i][entry];
1818 if (fVmax[i] < fVal[i][entry]) fVmax[i] = fVal[i][entry];
1819 }
1820 }
1821 }
1822}
1823
1824////////////////////////////////////////////////////////////////////////////////
1825/// Called at the end of a loop on a TTree.
1826
1828{
1829 // We take action (in Process) when we reach GetEstimate but
1830 // we reset at the beginning of Process, so
1831 // if fNfill == GetEstimate(), we just took action.
1832 if (fNfill && fNfill < fTree->GetEstimate())
1833 TakeAction();
1834
1835 if ((fSelectedRows == 0) && (TestBit(kCustomHistogram) == 0)) fDraw = 1; // do not draw
1836
1838}
int Int_t
Definition RtypesCore.h:45
long long Long64_t
Definition RtypesCore.h:80
#define BIT(n)
Definition Rtypes.h:85
#define ClassImp(name)
Definition Rtypes.h:377
#define gDirectory
Definition TDirectory.h:384
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
float xmin
float ymin
float xmax
float ymax
Int_t gDebug
Definition TROOT.cxx:595
#define gROOT
Definition TROOT.h:406
const Int_t kCustomHistogram
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
#define gPad
virtual Int_t GetNdim() const
Definition TFormula.h:237
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:30
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:31
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:39
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:35
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:34
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition TAttMarker.h:32
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:31
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition TAttMarker.h:33
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
Class to manage histogram axis.
Definition TAxis.h:31
Double_t GetXmax() const
Definition TAxis.h:140
virtual Int_t FindFixBin(Double_t x) const
Find bin number corresponding to abscissa x.
Definition TAxis.cxx:419
Double_t GetXmin() const
Definition TAxis.h:139
Int_t GetNbins() const
Definition TAxis.h:125
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition TAxis.cxx:540
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:528
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition TAxis.cxx:458
Container of bits.
Definition TBits.h:26
UInt_t GetNbits() const
Definition TBits.h:134
UInt_t FirstSetBit(UInt_t startBit=0) const
Return position of first non null bit (starting from position 0 and up)
Definition TBits.cxx:359
static TClass * Class()
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
void SetName(const char *name)
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
static void InitializeColors()
Initialize colors used by the TCanvas based graphics (via TColor objects).
Definition TColor.cxx:1148
A specialized string object used for TTree selections.
Definition TCut.h:25
A list of entries and subentries in a TTree or TChain.
virtual TEntryListArray * GetSubListForEntry(Long64_t entry, TTree *tree=nullptr)
Return the list holding the subentries for the given entry or 0.
static TClass * Class()
virtual bool Enter(Long64_t entry, TTree *tree, Long64_t subentry)
Add entry #entry (, #subentry) to the list.
A List of entry numbers in a TTree or TChain.
Definition TEntryList.h:26
virtual bool Enter(Long64_t entry, TTree *tree=nullptr)
Add entry #entry to the list.
static TClass * Class()
virtual Int_t Contains(Long64_t entry, TTree *tree=nullptr)
virtual void Reset()
Reset this list.
virtual bool GetReapplyCut() const
Definition TEntryList.h:82
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
virtual void SetValue(const char *name, const char *value, EEnvLevel level=kEnvChange, const char *type=nullptr)
Set the value of a resource or create a new resource.
Definition TEnv.cxx:736
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.
virtual Int_t GetIndex(Long64_t entry) const
Return index in the list of element with value entry array is supposed to be sorted prior to this cal...
virtual void Enter(Long64_t entry)
Enter element entry into the list.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2319
Int_t GetN() const
Definition TGraph.h:131
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:809
virtual TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1406
virtual void SetEditable(Bool_t editable=kTRUE)
if editable=kFALSE, the graph cannot be modified with the mouse by default a TGraph is editable
Definition TGraph.cxx:2278
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:669
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:621
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
virtual void SetDirectory(TDirectory *dir)
By default, when a histogram is created, it is added to the list of histogram objects in the current ...
Definition TH1.cxx:8905
TAxis * GetZaxis()
Definition TH1.h:326
@ kAllAxes
Definition TH1.h:76
static TClass * Class()
virtual Int_t GetDimension() const
Definition TH1.h:283
@ kNoStats
Don't draw stats box.
Definition TH1.h:165
virtual Bool_t CanExtendAllAxes() const
Returns true if all axes are extendable.
Definition TH1.cxx:6604
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
Definition TH1.cxx:7071
TAxis * GetXaxis()
Definition TH1.h:324
TObject * FindObject(const char *name) const override
Search object named name in the list of functions.
Definition TH1.cxx:3857
TAxis * GetYaxis()
Definition TH1.h:325
virtual UInt_t SetCanExtend(UInt_t extendBitMask)
Make the histogram axes extendable / not extendable according to the bit mask returns the previous bi...
Definition TH1.cxx:6617
virtual void FillN(Int_t ntimes, const Double_t *x, const Double_t *w, Int_t stride=1)
Fill this histogram with an array x and weights w.
Definition TH1.cxx:3447
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition TH1.cxx:8988
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:8958
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:357
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:307
Service class for 2-D histogram classes.
Definition TH2.h:30
3-D histogram with a double per channel (see TH1 documentation)
Definition TH3.h:363
3-D histogram with a float per channel (see TH1 documentation)
Definition TH3.h:317
The 3-D histogram classes derived from the 1-D histogram classes.
Definition TH3.h:31
static THLimitsFinder * GetLimitsFinder()
Return pointer to the current finder.
virtual Int_t FindGoodLimits(TH1 *h, Double_t xmin, Double_t xmax)
Compute the best axis limits for the X axis.
Iterator of linked list.
Definition TList.h:191
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:576
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
An array of TObjects.
Definition TObjArray.h:31
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
TObject * Remove(TObject *obj) override
Remove object from array.
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:973
virtual void Delete(Option_t *option="")
Delete this object.
Definition TObject.cxx:248
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
virtual const char * GetTitle() const
Returns title of object.
Definition TObject.cxx:483
virtual TClass * IsA() const
Definition TObject.h:245
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
void ResetBit(UInt_t f)
Definition TObject.h:200
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:62
A 3D polymarker.
void SetPoint(Int_t n, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
virtual Int_t GetLastPoint() const
void Draw(Option_t *option="") override
Draws 3-D polymarker with its current attributes.
Profile2D histograms are used to display the mean value of Z and its error for each cell in X,...
Definition TProfile2D.h:27
Int_t Fill(const Double_t *v)
Definition TProfile2D.h:51
static TClass * Class()
Profile Histogram.
Definition TProfile.h:32
static TClass * Class()
void FillN(Int_t, const Double_t *, const Double_t *, Int_t) override
Fill this histogram with an array x and weights w.
Definition TProfile.h:63
A specialized TSelector for TTree::Draw.
void ProcessFill(Long64_t entry) override
Called in the entry loop for all entries accepted by Select.
TEntryListArray * fTreeElistArray
! Pointer to Tree Event list array
virtual void SetEstimate(Long64_t n)
Set number of entries to estimate variable limits.
virtual void InitArrays(Int_t newsize)
Initialization of the primitive type arrays if the new size is bigger than the available space.
void Terminate() override
Called at the end of a loop on a TTree.
Int_t fAction
! Action type
Int_t GetMultiplicity() const
TTreeFormulaManager * fManager
Pointer to the formula manager.
TTreeFormula * GetVar(Int_t i) const
Return the TTreeFormula corresponding to the i-th component of the request formula (where the compone...
bool * fVarMultiple
![fDimension] True if fVar[i] has a variable index
~TSelectorDraw() override
Selector destructor.
TSelectorDraw()
Default selector constructor.
Double_t * fW
![fSelectedRows]Local buffer for weights
virtual UInt_t SplitNames(const TString &varexp, std::vector< TString > &names)
Build Index array for names in varexp.
Long64_t fCurrentSubEntry
Current subentry when fSelectMultiple is true. Used to fill TEntryListArray.
TTreeFormula * fSelect
Pointer to selection formula.
Long64_t fSelectedRows
Number of selected entries.
Int_t fForceRead
Force Read flag.
virtual void ClearFormula()
Delete internal buffers.
virtual void TakeAction()
Execute action for object obj fNfill times.
Long64_t fOldEstimate
Value of Tree fEstimate when selector is called.
Double_t fWeight
Tree weight (see TTree::SetWeight)
Int_t fNfill
! Total number of histogram fills
TH1 * fOldHistogram
! Pointer to previously used histogram
virtual bool CompileVariables(const char *varexp="", const char *selection="")
Compile input variables and selection expression.
Double_t * fVmax
![fDimension] Maxima of varexp columns
Int_t fMultiplicity
Indicator of the variability of the size of entries.
TTreeFormula ** fVar
![fDimension] Array of pointers to variables formula
virtual void ProcessFillMultiple(Long64_t entry)
Called in the entry loop for all entries accepted by Select.
bool fSelectMultiple
True if selection has a variable index.
TTree * fTree
Pointer to current Tree.
Int_t fDimension
Dimension of the current expression.
TObject * fTreeElist
Pointer to Tree Event list.
virtual Double_t * GetVal(Int_t i) const
Return the last values corresponding to the i-th component of the formula being processed (where the ...
Int_t * fNbins
![fDimension] Number of bins per dimension
Double_t ** fVal
![fSelectedRows][fDimension] Local buffer for the variables
Long64_t fDraw
! Last entry loop number when object was drawn
virtual void TakeEstimate()
Estimate limits for 1-D, 2-D or 3-D objects.
bool Notify() override
This function is called at the first entry of a new tree in a chain.
bool fCleanElist
True if original Tree elist must be saved.
bool fObjEval
True if fVar1 returns an object (or pointer to).
void Begin(TTree *tree) override
Called every time a loop on the tree(s) starts.
virtual void ProcessFillObject(Long64_t entry)
Called in the entry loop for all entries accepted by Select.
Double_t * fVmin
![fDimension] Minima of varexp columns
virtual void SetStatus(Long64_t status)
Definition TSelector.h:67
TList * fInput
List of objects available during processing.
Definition TSelector.h:41
TString fOption
Option given to TTree::Process.
Definition TSelector.h:39
virtual void Abort(const char *why, EAbort what=kAbortProcess)
Abort processing.
const char * GetOption() const override
Definition TSelector.h:57
TObject * fObject
! Current object if processing object (vs. TTree)
Definition TSelector.h:40
virtual void ResetAbort()
Definition TSelector.h:74
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:417
void ToLower()
Change string to lower-case.
Definition TString.cxx:1182
const char * Data() const
Definition TString.h:376
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2356
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
Int_t GetColorPalette(Int_t i) const
Return color number i in current palette.
Definition TStyle.cxx:1097
Int_t GetNumberOfColors() const
Return number of colors in the color palette.
Definition TStyle.cxx:1171
Used to coordinate one or more TTreeFormula objects.
virtual Int_t GetNdata(bool forceLoadDim=false)
Return number of available instances in the formulas.
virtual void Add(TTreeFormula *)
Add a new formula to the list of formulas managed The manager of the formula will be changed and the ...
virtual bool Sync()
Synchronize all the formulae.
virtual Int_t GetMultiplicity() const
Used to pass a selection expression to the Tree drawing routine.
virtual void ResetLoading()
Tell the formula that we are going to request a new entry.
TTreeFormulaManager * GetManager() const
virtual Int_t GetMultiplicity() const
T EvalInstance(Int_t i=0, const char *stringStack[]=nullptr)
Evaluate this treeformula.
void SetQuickLoad(bool quick)
virtual void UpdateFormulaLeaves()
This function is called TTreePlayer::UpdateFormulaLeaves, itself called by TChain::LoadTree when a ne...
virtual void * EvalObject(Int_t i=0)
Evaluate this treeformula.
virtual void SetAxis(TAxis *axis=nullptr)
Set the axis (in particular get the type).
virtual TClass * EvalClass(Int_t oper) const
Evaluate the class of the operation oper.
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Long64_t GetEstimate() const
Definition TTree.h:467
virtual Double_t GetWeight() const
Definition TTree.h:544
virtual TEntryList * GetEntryList()
Returns the entry list assigned to this tree.
Definition TTree.cxx:5854
virtual void SetEstimate(Long64_t nentries=1000000)
Set number of entries to estimate variable limits.
Definition TTree.cxx:9110
virtual Long64_t GetReadEntry() const
Definition TTree.h:509
virtual TTree * GetTree() const
Definition TTree.h:517
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition TTree.cxx:9046
TEventList * GetEventList() const
Definition TTree.h:473
@ kForceRead
Definition TTree.h:251
virtual Int_t GetUpdate() const
Definition TTree.h:521
virtual Long64_t GetChainOffset() const
Definition TTree.h:456
See TView3D.
Definition TView.h:25
virtual Double_t * GetRmax()=0
virtual Double_t * GetRmin()=0
static TView * CreateView(Int_t system=1, const Double_t *rmin=nullptr, const Double_t *rmax=nullptr)
Create a concrete default 3-d view via the plug-in manager.
Definition TView.cxx:27
TGraphErrors * gr
Definition legend1.C:25
TH1F * h1
Definition legend1.C:5
Definition graph.py:1
TLine l
Definition textangle.C:4