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