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 }
590 if (!hkeep) {
591 hist->GetXaxis()->SetTitle(fVar[0]->GetTitle());
592 hist->SetBit(kCanDelete);
593 if (!opt.Contains("goff")) hist->SetDirectory(nullptr);
594 }
595 if (opt.Length() && opt.Contains("e")) hist->Sumw2();
596 }
597 fVar[0]->SetAxis(hist->GetXaxis());
598 fObject = hist;
599
600 // 2-D distribution
601 } else if (fDimension == 2 && !(optpara || optcandle)) {
602 fAction = 2;
603 if (!fOldHistogram || !optSame) {
604 fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y", 40);
605 fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x", 40);
606 if (opt.Contains("prof")) fNbins[1] = gEnv->GetValue("Hist.Binning.2D.Prof", 100);
607 if (optSame) {
608 TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
609 if (oldhtemp) {
610 fNbins[1] = oldhtemp->GetXaxis()->GetNbins();
611 fVmin[1] = oldhtemp->GetXaxis()->GetXmin();
612 fVmax[1] = oldhtemp->GetXaxis()->GetXmax();
613 fNbins[0] = oldhtemp->GetYaxis()->GetNbins();
614 fVmin[0] = oldhtemp->GetYaxis()->GetXmin();
615 fVmax[0] = oldhtemp->GetYaxis()->GetXmax();
616 } else {
617 fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x", 40);
618 fVmin[1] = gPad->GetUxmin();
619 fVmax[1] = gPad->GetUxmax();
620 fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y", 40);
621 fVmin[0] = gPad->GetUymin();
622 fVmax[0] = gPad->GetUymax();
623 }
624 } else {
625 if (!fOldHistogram) fAction = -2;
626 fVmin[1] = xmin;
627 fVmax[1] = xmax;
628 fVmin[0] = ymin;
629 fVmax[0] = ymax;
630 if (xmin < xmax && ymin < ymax) canExtend = false;
631 }
632 }
633 if (profile || opt.Contains("prof")) {
634 TProfile *hp;
635 if (fOldHistogram) {
636 fAction = 4;
638 } else {
639 if (fAction < 0) {
640 fAction = -4;
641 fVmin[1] = xmin;
642 fVmax[1] = xmax;
643 if (xmin < xmax) canExtend = false;
644 }
645 if (fAction == 2) {
646 //we come here when option = "same prof"
647 fAction = -4;
648 TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
649 if (oldhtemp) {
650 fNbins[1] = oldhtemp->GetXaxis()->GetNbins();
651 fVmin[1] = oldhtemp->GetXaxis()->GetXmin();
652 fVmax[1] = oldhtemp->GetXaxis()->GetXmax();
653 }
654 }
655 if (opt.Contains("profs")) {
656 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "s");
657 } else if (opt.Contains("profi")) {
658 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "i");
659 } else if (opt.Contains("profg")) {
660 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "g");
661 } else {
662 hp = new TProfile(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], "");
663 }
664 if (!hkeep) {
665 hp->SetBit(kCanDelete);
666 if (!opt.Contains("goff")) hp->SetDirectory(nullptr);
667 }
668 hp->SetLineColor(fTree->GetLineColor());
669 hp->SetLineWidth(fTree->GetLineWidth());
670 hp->SetLineStyle(fTree->GetLineStyle());
671 hp->SetFillColor(fTree->GetFillColor());
672 hp->SetFillStyle(fTree->GetFillStyle());
673 hp->SetMarkerStyle(fTree->GetMarkerStyle());
674 hp->SetMarkerColor(fTree->GetMarkerColor());
675 hp->SetMarkerSize(fTree->GetMarkerSize());
676 if (canExtend) hp->SetCanExtend(TH1::kAllAxes);
677 }
678 fVar[1]->SetAxis(hp->GetXaxis());
679 fObject = hp;
680
681 } else {
682 TH2 *h2;
683 if (fOldHistogram) {
684 h2 = (TH2F*)fOldHistogram;
685 } else {
686 TString precision = gEnv->GetValue("Hist.Precision.2D", "float");
687 if (precision.Contains("float")) {
688 h2 = new TH2F(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
689 } else {
690 h2 = new TH2D(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
691 }
692 h2->SetLineColor(fTree->GetLineColor());
693 h2->SetLineWidth(fTree->GetLineWidth());
694 h2->SetLineStyle(fTree->GetLineStyle());
695 h2->SetFillColor(fTree->GetFillColor());
696 h2->SetFillStyle(fTree->GetFillStyle());
697 h2->SetMarkerStyle(fTree->GetMarkerStyle());
698 h2->SetMarkerColor(fTree->GetMarkerColor());
699 h2->SetMarkerSize(fTree->GetMarkerSize());
700 if (canExtend) h2->SetCanExtend(TH1::kAllAxes);
701 if (!hkeep) {
702 h2->GetXaxis()->SetTitle(fVar[1]->GetTitle());
703 h2->GetYaxis()->SetTitle(fVar[0]->GetTitle());
704 h2->SetBit(TH1::kNoStats);
705 h2->SetBit(kCanDelete);
706 if (!opt.Contains("goff")) h2->SetDirectory(nullptr);
707 }
708 }
709 fVar[0]->SetAxis(h2->GetYaxis());
710 fVar[1]->SetAxis(h2->GetXaxis());
711 bool graph = false;
712 Int_t l = opt.Length();
713 if (l == 0 || optSame) graph = true;
714 if (opt.Contains("p") || opt.Contains("*") || opt.Contains("l")) graph = true;
715 if (opt.Contains("surf") || opt.Contains("lego") || opt.Contains("cont")) graph = false;
716 if (opt.Contains("col") || opt.Contains("hist") || opt.Contains("scat")) graph = false;
717 if (opt.Contains("box")) graph = false;
718 fObject = h2;
719 if (graph) {
720 fAction = 12;
721 if (!fOldHistogram && !optSame) fAction = -12;
722 }
723 }
724
725 // 3-D distribution
726 } else if ((fDimension == 3 || fDimension == 4) && !(optpara || optcandle)) {
727 fAction = 3;
728 if (fDimension == 4) fAction = 40;
729 if (!fOldHistogram || !optSame) {
730 fNbins[0] = gEnv->GetValue("Hist.Binning.3D.z", 20);
731 fNbins[1] = gEnv->GetValue("Hist.Binning.3D.y", 20);
732 fNbins[2] = gEnv->GetValue("Hist.Binning.3D.x", 20);
733 if (fDimension == 3 && opt.Contains("prof")) {
734 fNbins[1] = gEnv->GetValue("Hist.Binning.3D.Profy", 20);
735 fNbins[2] = gEnv->GetValue("Hist.Binning.3D.Profx", 20);
736 } else if (fDimension == 3 && opt.Contains("col")) {
737 fNbins[0] = gEnv->GetValue("Hist.Binning.2D.y", 40);
738 fNbins[1] = gEnv->GetValue("Hist.Binning.2D.x", 40);
739 }
740 if (optSame) {
741 TH1 *oldhtemp = (TH1*)gPad->FindObject(hdefault);
742 if (oldhtemp) {
743 fNbins[2] = oldhtemp->GetXaxis()->GetNbins();
744 fVmin[2] = oldhtemp->GetXaxis()->GetXmin();
745 fVmax[2] = oldhtemp->GetXaxis()->GetXmax();
746 fNbins[1] = oldhtemp->GetYaxis()->GetNbins();
747 fVmin[1] = oldhtemp->GetYaxis()->GetXmin();
748 fVmax[1] = oldhtemp->GetYaxis()->GetXmax();
749 fNbins[0] = oldhtemp->GetZaxis()->GetNbins();
750 fVmin[0] = oldhtemp->GetZaxis()->GetXmin();
751 fVmax[0] = oldhtemp->GetZaxis()->GetXmax();
752 } else {
753 TView *view = gPad->GetView();
754 if (!view) {
755 Error("Begin", "You cannot use option same when no 3D view exists");
756 fVmin[0] = fVmin[1] = fVmin[2] = -1;
757 fVmax[0] = fVmax[1] = fVmax[2] = 1;
758 view = TView::CreateView(1, fVmin, fVmax);
759 }
760 Double_t *rmin = view->GetRmin();
761 Double_t *rmax = view->GetRmax();
762 fNbins[2] = gEnv->GetValue("Hist.Binning.3D.z", 20);
763 fVmin[2] = rmin[0];
764 fVmax[2] = rmax[0];
765 fNbins[1] = gEnv->GetValue("Hist.Binning.3D.y", 20);
766 fVmin[1] = rmin[1];
767 fVmax[1] = rmax[1];
768 fNbins[0] = gEnv->GetValue("Hist.Binning.3D.x", 20);
769 fVmin[0] = rmin[2];
770 fVmax[0] = rmax[2];
771 }
772 } else {
773 if (!fOldHistogram && fDimension == 3) fAction = -3;
774 fVmin[2] = xmin;
775 fVmax[2] = xmax;
776 fVmin[1] = ymin;
777 fVmax[1] = ymax;
778 fVmin[0] = zmin;
779 fVmax[0] = zmax;
780 if (xmin < xmax && ymin < ymax && zmin < zmax) canExtend = false;
781 }
782 }
783 if ((fDimension == 3) && (profile || opt.Contains("prof"))) {
784 TProfile2D *hp;
785 if (fOldHistogram) {
786 fAction = 23;
788 } else {
789 if (fAction < 0) {
790 fAction = -23;
791 fVmin[2] = xmin;
792 fVmax[2] = xmax;
793 fVmin[1] = ymin;
794 fVmax[1] = ymax;
795 if (xmin < xmax && ymin < ymax) canExtend = false;
796 }
797 if (opt.Contains("profs")) {
798 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "s");
799 } else if (opt.Contains("profi")) {
800 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "i");
801 } else if (opt.Contains("profg")) {
802 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "g");
803 } else {
804 hp = new TProfile2D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], "");
805 }
806 if (!hkeep) {
807 hp->SetBit(kCanDelete);
808 if (!opt.Contains("goff")) hp->SetDirectory(nullptr);
809 }
810 hp->SetLineColor(fTree->GetLineColor());
811 hp->SetLineWidth(fTree->GetLineWidth());
812 hp->SetLineStyle(fTree->GetLineStyle());
813 hp->SetFillColor(fTree->GetFillColor());
814 hp->SetFillStyle(fTree->GetFillStyle());
815 hp->SetMarkerStyle(fTree->GetMarkerStyle());
816 hp->SetMarkerColor(fTree->GetMarkerColor());
817 hp->SetMarkerSize(fTree->GetMarkerSize());
818 if (canExtend) hp->SetCanExtend(TH1::kAllAxes);
819 }
820 fVar[1]->SetAxis(hp->GetYaxis());
821 fVar[2]->SetAxis(hp->GetXaxis());
822 fObject = hp;
823 } else if (fDimension == 3 && opt.Contains("col")) {
824 TH2F *h2;
825 if (fOldHistogram) {
826 h2 = (TH2F*)fOldHistogram;
827 } else {
828 h2 = new TH2F(hname, htitle.Data(), fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
829 h2->SetLineColor(fTree->GetLineColor());
830 h2->SetLineWidth(fTree->GetLineWidth());
831 h2->SetLineStyle(fTree->GetLineStyle());
832 h2->SetFillColor(fTree->GetFillColor());
833 h2->SetFillStyle(fTree->GetFillStyle());
834 h2->SetMarkerStyle(fTree->GetMarkerStyle());
835 h2->SetMarkerColor(fTree->GetMarkerColor());
836 h2->SetMarkerSize(fTree->GetMarkerSize());
837 if (canExtend) h2->SetCanExtend(TH1::kAllAxes);
838 if (!hkeep) {
839 h2->GetXaxis()->SetTitle(fVar[1]->GetTitle());
840 h2->GetYaxis()->SetTitle(fVar[0]->GetTitle());
841 h2->GetZaxis()->SetTitle(fVar[2]->GetTitle());
842 h2->SetBit(TH1::kNoStats);
843 h2->SetBit(kCanDelete);
844 if (!opt.Contains("goff")) h2->SetDirectory(nullptr);
845 }
846 }
847 fVar[0]->SetAxis(h2->GetYaxis());
848 fVar[1]->SetAxis(h2->GetXaxis());
849 fObject = h2;
850 fAction = 33;
851 } else {
852 TH3 *h3;
853 if (fOldHistogram) {
854 h3 = (TH3F*)fOldHistogram;
855 } else {
856 TString precision = gEnv->GetValue("Hist.Precision.3D", "float");
857 if (precision.Contains("float")) {
858 h3 = new TH3F(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
859 } else {
860 h3 = new TH3D(hname, htitle.Data(), fNbins[2], fVmin[2], fVmax[2], fNbins[1], fVmin[1], fVmax[1], fNbins[0], fVmin[0], fVmax[0]);
861 }
862 h3->SetLineColor(fTree->GetLineColor());
863 h3->SetLineWidth(fTree->GetLineWidth());
864 h3->SetLineStyle(fTree->GetLineStyle());
865 h3->SetFillColor(fTree->GetFillColor());
866 h3->SetFillStyle(fTree->GetFillStyle());
867 h3->SetMarkerStyle(fTree->GetMarkerStyle());
868 h3->SetMarkerColor(fTree->GetMarkerColor());
869 h3->SetMarkerSize(fTree->GetMarkerSize());
870 if (canExtend) h3->SetCanExtend(TH1::kAllAxes);
871 if (!hkeep) {
872 //small correction for the title offsets in x,y to take into account the angles
873 Double_t xoffset = h3->GetXaxis()->GetTitleOffset();
874 Double_t yoffset = h3->GetYaxis()->GetTitleOffset();
875 h3->GetXaxis()->SetTitleOffset(1.2 * xoffset);
876 h3->GetYaxis()->SetTitleOffset(1.2 * yoffset);
877 h3->GetXaxis()->SetTitle(fVar[2]->GetTitle());
878 h3->GetYaxis()->SetTitle(fVar[1]->GetTitle());
879 h3->GetZaxis()->SetTitle(fVar[0]->GetTitle());
880 h3->SetBit(kCanDelete);
881 h3->SetBit(TH1::kNoStats);
882 if (!opt.Contains("goff")) h3->SetDirectory(nullptr);
883 }
884 }
885 fVar[0]->SetAxis(h3->GetZaxis());
886 fVar[1]->SetAxis(h3->GetYaxis());
887 fVar[2]->SetAxis(h3->GetXaxis());
888 fObject = h3;
890 if (optSame) noscat -= 4;
891 if (!noscat && fDimension == 3) {
892 fAction = 13;
893 if (!fOldHistogram && !optSame) fAction = -13;
894 }
895 }
896 // An Event List
897 } else if (enlist) {
898 fAction = 5;
900 fTree->SetEstimate(1);
901 fObject = enlist;
902 } else if (evlist) {
903 fAction = 5;
905 fTree->SetEstimate(1);
906 fObject = evlist;
907 } else if (optcandle || optpara || opt5d) {
908 if (optcandle) fAction = 7;
909 else if (opt5d) fAction = 8;
910 else fAction = 6;
911 }
912 if (varexp) delete[] varexp;
913 for (i = 0; i < fValSize; ++i)
914 fVarMultiple[i] = false;
915 fSelectMultiple = false;
916 for (i = 0; i < fDimension; ++i) {
917 if (fVar[i] && fVar[i]->GetMultiplicity()) fVarMultiple[i] = true;
918 }
919
921
924 fNfill = 0;
925
926 for (i = 0; i < fDimension; ++i) {
927 if (!fVal[i] && fVar[i]) {
928 fVal[i] = new Double_t[(Int_t)fTree->GetEstimate()];
929 }
930 }
931
932 if (!fW) fW = new Double_t[(Int_t)fTree->GetEstimate()];
933
934 for (i = 0; i < fValSize; ++i) {
935 fVmin[i] = DBL_MAX;
936 fVmax[i] = -DBL_MAX;
937 }
938}
939
940////////////////////////////////////////////////////////////////////////////////
941/// Delete internal buffers.
942
944{
946 for (Int_t i = 0; i < fValSize; ++i) {
947 delete fVar[i];
948 fVar[i] = nullptr;
949 }
950 delete fSelect; fSelect = nullptr;
951 fManager = nullptr;
952 fMultiplicity = 0;
953}
954
955////////////////////////////////////////////////////////////////////////////////
956/// Compile input variables and selection expression.
957///
958/// varexp is an expression of the general form e1:e2:e3
959/// where e1,etc is a formula referencing a combination of the columns
960///
961/// Example:
962///
963/// varexp = x simplest case: draw a 1-Dim distribution of column named x
964/// = sqrt(x) : draw distribution of sqrt(x)
965/// = x*y/z
966/// = y:sqrt(x) 2-Dim distribution of y versus sqrt(x)
967///
968/// selection is an expression with a combination of the columns
969///
970/// Example:
971///
972/// selection = "x<y && sqrt(z)>3.2"
973///
974/// in a selection all the C++ operators are authorized
975///
976/// Return false if any of the variable is not compilable.
977
979{
980 Int_t i, nch, ncols;
981
982 // Compile selection expression if there is one
983 fDimension = 0;
984 ClearFormula();
985 fMultiplicity = 0;
986 fObjEval = false;
987
988 if (strlen(selection)) {
989 fSelect = new TTreeFormula("Selection", selection, fTree);
990 fSelect->SetQuickLoad(true);
991 if (!fSelect->GetNdim()) {
992 delete fSelect;
993 fSelect = nullptr;
994 return false;
995 }
996 }
997
998 // if varexp is empty, take first column by default
999 nch = strlen(varexp);
1000 if (nch == 0) {
1001 fDimension = 0;
1002 if (fSelect) {
1004 }
1006
1007 if (fManager) {
1008 fManager->Sync();
1009
1012 }
1013
1014 return true;
1015 }
1016
1017 // otherwise select only the specified columns
1018 std::vector<TString> varnames;
1020
1022
1024 if (fSelect) fManager->Add(fSelect);
1026 for (i = 0; i < ncols; ++i) {
1027 fVar[i] = new TTreeFormula(TString::Format("Var%i", i + 1), varnames[i].Data(), fTree);
1028 fVar[i]->SetQuickLoad(true);
1029 if(!fVar[i]->GetNdim()) { ClearFormula(); return false; }
1030 fManager->Add(fVar[i]);
1031 }
1032 fManager->Sync();
1033
1036
1037 fDimension = ncols;
1038
1039 if (ncols == 1) {
1040 TClass *cl = fVar[0]->EvalClass();
1041 if (cl) {
1042 fObjEval = true;
1043 }
1044 }
1045 return true;
1046}
1047
1048////////////////////////////////////////////////////////////////////////////////
1049/// Return the last values corresponding to the i-th component
1050/// of the formula being processed (where the component are ':' separated).
1051/// The actual number of entries is:
1052///
1053/// GetSelectedRows() % tree->GetEstimate()
1054///
1055/// Note GetSelectedRows currently returns the actual number of values plotted
1056/// and thus if the formula contains arrays, this number might be greater than
1057/// the number of entries in the trees.
1058///
1059/// By default TTree::Draw creates the arrays obtained
1060/// with all GetVal and GetW with a length corresponding to the
1061/// parameter fEstimate. By default fEstimate=10000 and can be modified
1062/// via TTree::SetEstimate. A possible recipe is to do
1063///
1064/// tree->SetEstimate(tree->GetEntries());
1065///
1066/// You must call SetEstimate if the expected number of selected rows
1067/// is greater than 10000.
1068///
1069/// See TTree::Draw for additional details.
1070
1072{
1074 return nullptr;
1075 else
1076 return fVal[i];
1077}
1078
1079////////////////////////////////////////////////////////////////////////////////
1080/// Return the TTreeFormula corresponding to the i-th component
1081/// of the request formula (where the component are ':' separated).
1082
1084{
1086 return nullptr;
1087 else
1088 return fVar[i];
1089}
1090
1091////////////////////////////////////////////////////////////////////////////////
1092/// Initialization of the primitive type arrays if the new size is bigger than the available space.
1093
1095{
1096 if (newsize > fValSize) {
1098 while (fValSize < newsize)
1099 fValSize *= 2; // Double the available space until it matches the new size.
1100 if (fNbins) delete [] fNbins;
1101 if (fVmin) delete [] fVmin;
1102 if (fVmax) delete [] fVmax;
1103 if (fVarMultiple) delete [] fVarMultiple;
1104
1105 fNbins = new Int_t[fValSize];
1106 fVmin = new Double_t[fValSize];
1107 fVmax = new Double_t[fValSize];
1108 fVarMultiple = new bool[fValSize];
1109
1110 for (Int_t i = 0; i < oldsize; ++i)
1111 delete [] fVal[i];
1112 delete [] fVal;
1113 delete [] fVar;
1114 fVal = new Double_t*[fValSize];
1115 fVar = new TTreeFormula*[fValSize];
1116 for (Int_t i = 0; i < fValSize; ++i) {
1117 fVal[i] = nullptr;
1118 fVar[i] = nullptr;
1119 }
1120 }
1121}
1122
1123////////////////////////////////////////////////////////////////////////////////
1124/// Build Index array for names in varexp.
1125/// This will allocated a C style array of TString and Ints
1126
1127UInt_t TSelectorDraw::SplitNames(const TString &varexp, std::vector<TString> &names)
1128{
1129 names.clear();
1130
1131 bool ternary = false;
1132 Int_t prev = 0;
1133 for (Int_t i = 0; i < varexp.Length(); i++) {
1134 if (varexp[i] == ':'
1135 && !((i > 0 && varexp[i-1] == ':') || varexp[i+1] == ':')
1136 ) {
1137 if (ternary) {
1138 ternary = false;
1139 } else {
1140 names.push_back(varexp(prev, i - prev));
1141 prev = i + 1;
1142 }
1143 }
1144 if (varexp[i] == '?') {
1145 ternary = true;
1146 }
1147 }
1148 names.push_back(varexp(prev, varexp.Length() - prev));
1149 return names.size();
1150}
1151
1152
1153////////////////////////////////////////////////////////////////////////////////
1154/// This function is called at the first entry of a new tree in a chain.
1155
1157{
1158 if (fTree) fWeight = fTree->GetWeight();
1159 if (fVar) {
1160 for (Int_t i = 0; i < fDimension; ++i) {
1161 if (fVar[i]) fVar[i]->UpdateFormulaLeaves();
1162 }
1163 }
1165 return true;
1166}
1167
1168////////////////////////////////////////////////////////////////////////////////
1169/// Called in the entry loop for all entries accepted by Select.
1170
1172{
1173 if (fObjEval) {
1175 return;
1176 }
1177
1178 if (fMultiplicity) {
1180 return;
1181 }
1182
1183 if (fNfill >= fTree->GetEstimate())
1184 fNfill = 0;
1185
1186 // simple case with no multiplicity
1187 if (fForceRead && fManager->GetNdata() <= 0) return;
1188
1189 if (fSelect) {
1191 if (!fW[fNfill]) return;
1192 } else fW[fNfill] = fWeight;
1193 if (fVal) {
1194 for (Int_t i = 0; i < fDimension; ++i) {
1195 if (fVar[i]) fVal[i][fNfill] = fVar[i]->EvalInstance(0);
1196 }
1197 }
1198 fNfill++;
1199 if (fNfill >= fTree->GetEstimate()) {
1200 TakeAction();
1201 }
1202}
1203
1204////////////////////////////////////////////////////////////////////////////////
1205/// Called in the entry loop for all entries accepted by Select.
1206/// Complex case with multiplicity.
1207
1209{
1210 if (fNfill >= fTree->GetEstimate())
1211 fNfill = 0;
1212
1213 // Grab the array size of the formulas for this entry
1215
1216 // No data at all, let's move on to the next entry.
1217 if (!ndata) return;
1218
1219 // If the entry list is a TEntryListArray, get the selected subentries for this entry
1220 TEntryList *subList = nullptr;
1221 if (fTreeElistArray) {
1223 }
1224
1226
1227 // Calculate the first values
1228 if (fSelect) {
1229 // coverity[var_deref_model] fSelectMultiple==true => fSelect != 0
1231 if (!fW[fNfill] && !fSelectMultiple) return;
1232 } else fW[fNfill] = fWeight;
1233
1234 // Always call EvalInstance(0) to insure the loading
1235 // of the branches.
1236 if (fW[fNfill] && (!subList || subList->Contains(0))) {
1237 if (fDimension == 0 && fSelectMultiple) fCurrentSubEntry = (Long64_t) 0; // to fill TEntryListArray
1238 for (Int_t i = 0; i < fDimension; ++i) {
1239 if (fVar[i]) fVal[i][fNfill] = fVar[i]->EvalInstance(0);
1240 }
1241 fNfill++;
1242 if (fNfill >= fTree->GetEstimate()) {
1243 TakeAction();
1244 }
1245 } else {
1246 for (Int_t i = 0; i < fDimension; ++i) {
1247 if (fVar[i]) fVar[i]->ResetLoading();
1248 }
1249 }
1250 Double_t ww = fW[nfill0];
1251
1252 for (Int_t i = 1; i < ndata; i++) {
1253 if (fNfill >= fTree->GetEstimate())
1254 fNfill = 0;
1255 if (subList && !subList->Contains(i))
1256 continue;
1257 if (fSelectMultiple) {
1258 // coverity[var_deref_model] fSelectMultiple==true => fSelect != 0
1259 ww = fWeight * fSelect->EvalInstance(i);
1260 if (ww == 0) continue;
1261 if (fNfill == nfill0) {
1262 for (Int_t k = 0; k < fDimension; ++k) {
1263 if (!fVarMultiple[k]) fVal[k][fNfill] = fVar[k]->EvalInstance(0);
1264 }
1265 }
1266 if (fDimension == 0) fCurrentSubEntry = (Long64_t) i; // to fill TEntryListArray
1267 }
1268 for (Int_t k = 0; k < fDimension; ++k) {
1269 if (fVarMultiple[k]) fVal[k][fNfill] = fVar[k]->EvalInstance(i);
1270 else fVal[k][fNfill] = fVal[k][nfill0];
1271 }
1272 fW[fNfill] = ww;
1273
1274 fNfill++;
1275 if (fNfill >= fTree->GetEstimate()) {
1276 TakeAction();
1277 }
1278 }
1279}
1280
1281////////////////////////////////////////////////////////////////////////////////
1282/// Called in the entry loop for all entries accepted by Select.
1283/// Case where the only variable returns an object (or pointer to).
1284
1286{
1287 // Complex case with multiplicity.
1288
1289 if (fNfill >= fTree->GetEstimate())
1290 fNfill = 0;
1291
1292 // Grab the array size of the formulas for this entry
1294
1295 // No data at all, let's move on to the next entry.
1296 if (!ndata) return;
1297
1299 Double_t ww = 0;
1300
1301 for (Int_t i = 0; i < ndata; i++) {
1302 if (i == 0) {
1303 if (fSelect) {
1305 if (!fW[fNfill] && !fSelectMultiple) return;
1306 } else fW[fNfill] = fWeight;
1307 ww = fW[nfill0];
1308 } else if (fSelectMultiple) {
1309 ww = fWeight * fSelect->EvalInstance(i);
1310 if (ww == 0) continue;
1311 }
1312 if (fDimension >= 1 && fVar[0]) {
1313 TClass *cl = fVar[0]->EvalClass();
1314 if (cl == TBits::Class()) {
1315
1316 void *obj = fVar[0]->EvalObject(i);
1317
1318 if (obj) {
1319 TBits *bits = (TBits*)obj;
1320 Int_t nbits = bits->GetNbits();
1321
1322 Int_t nextbit = -1;
1323 while (true) {
1324 nextbit = bits->FirstSetBit(nextbit + 1);
1325 if (nextbit >= nbits) break;
1326 fVal[0][fNfill] = nextbit;
1327 fW[fNfill] = ww;
1328 fNfill++;
1329 }
1330 }
1331
1332 } else {
1333
1334 if (!TestBit(kWarn)) {
1335 Warning("ProcessFillObject",
1336 "Not implemented for %s",
1337 cl ? cl->GetName() : "unknown class");
1338 SetBit(kWarn);
1339 }
1340
1341 }
1342 }
1343 if (fNfill >= fTree->GetEstimate()) {
1344 TakeAction();
1345 }
1346 }
1347
1348}
1349
1350////////////////////////////////////////////////////////////////////////////////
1351/// Set number of entries to estimate variable limits.
1352
1354{
1355 if (fVal) {
1356 for (Int_t i = 0; i < fValSize; ++i) {
1357 delete [] fVal[i];
1358 fVal[i] = nullptr;
1359 }
1360 }
1361 delete [] fW; fW = nullptr;
1362}
1363
1364////////////////////////////////////////////////////////////////////////////////
1365/// Execute action for object obj fNfill times.
1366
1368{
1369 Int_t i;
1370 //__________________________1-D histogram_______________________
1371 if (fAction == 1)((TH1*)fObject)->FillN(fNfill, fVal[0], fW);
1372 //__________________________2-D histogram_______________________
1373 else if (fAction == 2) {
1374 TH2 *h2 = (TH2*)fObject;
1375 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1376 }
1377 //__________________________Profile histogram_______________________
1378 else if (fAction == 4)((TProfile*)fObject)->FillN(fNfill, fVal[1], fVal[0], fW);
1379 //__________________________Event List______________________________
1380 else if (fAction == 5) {
1383 Long64_t enumb = fTree->GetTree()->GetReadEntry();
1384 enlistarray->Enter(enumb, nullptr, fCurrentSubEntry);
1385 } else if (fObject->InheritsFrom(TEntryList::Class())) {
1387 Long64_t enumb = fTree->GetTree()->GetReadEntry();
1388 enlist->Enter(enumb);
1389 } else {
1391 Long64_t enumb = fTree->GetChainOffset() + fTree->GetTree()->GetReadEntry();
1392 if (evlist->GetIndex(enumb) < 0) evlist->Enter(enumb);
1393 }
1394 }
1395 //__________________________2D scatter plot_______________________
1396 else if (fAction == 12) {
1397 TH2 *h2 = (TH2*)fObject;
1398 if (h2->CanExtendAllAxes() && h2->TestBit(kCanDelete)) {
1399 for (i = 0; i < fNfill; i++) {
1400 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1401 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1402 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1403 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1404 }
1405 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2, fVmin[1], fVmax[1], fVmin[0], fVmax[0]);
1406 }
1407 TGraph *pm = new TGraph(fNfill, fVal[1], fVal[0]);
1408 pm->SetEditable(false);
1409 pm->SetBit(kCanDelete);
1410 pm->SetMarkerStyle(fTree->GetMarkerStyle());
1411 pm->SetMarkerColor(fTree->GetMarkerColor());
1412 pm->SetMarkerSize(fTree->GetMarkerSize());
1413 pm->SetLineColor(fTree->GetLineColor());
1414 pm->SetLineWidth(fTree->GetLineWidth());
1415 pm->SetLineStyle(fTree->GetLineStyle());
1416 pm->SetFillColor(fTree->GetFillColor());
1417 pm->SetFillStyle(fTree->GetFillStyle());
1418
1419 if (!fDraw && !strstr(fOption.Data(), "goff")) {
1420 if (fOption.Length() == 0 || strcasecmp(fOption.Data(), "same") == 0) pm->Draw("p");
1421 else pm->Draw(fOption.Data());
1422 }
1423 if (!h2->TestBit(kCanDelete)) {
1424 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1425 }
1426 }
1427 //__________________________3D scatter plot_______________________
1428 else if (fAction == 3) {
1429 TH3 *h3 = (TH3*)fObject;
1430 if (!h3->TestBit(kCanDelete)) {
1431 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1432 }
1433 } else if (fAction == 13) {
1435 pm3d->SetMarkerStyle(fTree->GetMarkerStyle());
1436 pm3d->SetMarkerColor(fTree->GetMarkerColor());
1437 pm3d->SetMarkerSize(fTree->GetMarkerSize());
1438 for (i = 0; i < fNfill; i++) {
1439 pm3d->SetPoint(i, fVal[2][i], fVal[1][i], fVal[0][i]);
1440 }
1441 pm3d->Draw();
1442 TH3 *h3 = (TH3*)fObject;
1443 if (!h3->TestBit(kCanDelete)) {
1444 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1445 }
1446 }
1447 //__________________________3D scatter plot (3rd variable = col)__
1448 else if (fAction == 33) {
1449 TH2 *h2 = (TH2*)fObject;
1450 TakeEstimate();
1451 Int_t ncolors = gStyle->GetNumberOfColors();
1452 TObjArray *grs = (TObjArray*)h2->GetListOfFunctions()->FindObject("graphs");
1453 Int_t col;
1454 TGraph *gr;
1455 if (!grs) {
1456 grs = new TObjArray(ncolors);
1457 grs->SetOwner();
1458 grs->SetName("graphs");
1459 h2->GetListOfFunctions()->Add(grs, "P");
1460 for (col = 0; col < ncolors; col++) {
1461 gr = new TGraph();
1465 grs->AddAt(gr, col);
1466 }
1467 }
1468 h2->SetEntries(fNfill);
1469 h2->SetMinimum(fVmin[2]);
1470 h2->SetMaximum(fVmax[2]);
1471 // Fill the graphs according to the color
1472 for (i = 0; i < fNfill; i++) {
1473 col = Int_t(ncolors * ((fVal[2][i] - fVmin[2]) / (fVmax[2] - fVmin[2])));
1474 if (col < 0) col = 0;
1475 if (col > ncolors - 1) col = ncolors - 1;
1476 gr = (TGraph*)grs->UncheckedAt(col);
1477 if (gr) gr->SetPoint(gr->GetN(), fVal[1][i], fVal[0][i]);
1478 }
1479 // Remove potential empty graphs
1480 for (col = 0; col < ncolors; col++) {
1481 gr = (TGraph*)grs->At(col);
1482 if (gr && gr->GetN() <= 0) grs->Remove(gr);
1483 }
1484 }
1485 //__________________________2D Profile Histogram__________________
1486 else if (fAction == 23) {
1488 for (i = 0; i < fNfill; i++) hp2->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1489 }
1490 //__________________________4D scatter plot_______________________
1491 else if (fAction == 40) {
1492 TakeEstimate();
1493 TH3 *h3 = (TH3*)fObject;
1494 Int_t ncolors = gStyle->GetNumberOfColors();
1495 if (ncolors == 0) {
1497 ncolors = gStyle->GetNumberOfColors();
1498 }
1499 TObjArray *pms = (TObjArray*)h3->GetListOfFunctions()->FindObject("polymarkers");
1500 Int_t col;
1502 if (!pms) {
1503 pms = new TObjArray(ncolors);
1504 pms->SetOwner();
1505 pms->SetName("polymarkers");
1506 h3->GetListOfFunctions()->Add(pms);
1507 for (col = 0; col < ncolors; col++) {
1508 pm3d = new TPolyMarker3D();
1509 pm3d->SetMarkerColor(gStyle->GetColorPalette(col));
1510 pm3d->SetMarkerStyle(fTree->GetMarkerStyle());
1511 pm3d->SetMarkerSize(fTree->GetMarkerSize());
1512 pms->AddAt(pm3d, col);
1513 }
1514 }
1515 h3->SetEntries(fNfill);
1516 h3->SetMinimum(fVmin[3]);
1517 h3->SetMaximum(fVmax[3]);
1518 for (i = 0; i < fNfill; i++) {
1519 col = Int_t(ncolors * ((fVal[3][i] - fVmin[3]) / (fVmax[3] - fVmin[3])));
1520 if (col > ncolors-1) col = ncolors-1;
1521 if (col < 0) col = 0;
1522 pm3d = (TPolyMarker3D*)pms->UncheckedAt(col);
1523 pm3d->SetPoint(pm3d->GetLastPoint() + 1, fVal[2][i], fVal[1][i], fVal[0][i]);
1524 }
1525 }
1526 //__________________________Parallel coordinates / candle chart_______________________
1527 else if (fAction == 6 || fAction == 7) {
1528 TakeEstimate();
1529 bool candle = (fAction == 7);
1530 // Using CINT to avoid a dependency in TParallelCoord
1531 if (!fOption.Contains("goff"))
1532 gROOT->ProcessLine(TString::Format("TParallelCoord::BuildParallelCoord((TSelectorDraw*)0x%zx,0x%zx)",
1533 (size_t)this, (size_t)candle));
1534 } else if (fAction == 8) {
1535 //gROOT->ProcessLineFast(TString::Format("(new TGL5DDataSet((TTree *)0x%1x))->Draw(\"%s\");", fTree, fOption.Data()));
1536 }
1537 //__________________________something else_______________________
1538 else if (fAction < 0) {
1539 fAction = -fAction;
1540 TakeEstimate();
1541 }
1542
1543 // Do we need to update screen?
1545 if (!fTree->GetUpdate()) return;
1546 if (fSelectedRows > fDraw + fTree->GetUpdate()) {
1547 if (fDraw) gPad->Modified();
1548 else fObject->Draw(fOption.Data());
1549 gPad->Update();
1551 }
1552}
1553
1554////////////////////////////////////////////////////////////////////////////////
1555/// Estimate limits for 1-D, 2-D or 3-D objects.
1556
1558{
1559 Int_t i;
1560 Double_t rmin[3], rmax[3];
1561 Double_t vminOld[4], vmaxOld[4];
1562 for (i = 0; i < fValSize && i < 4; i++) {
1563 vminOld[i] = fVmin[i];
1564 vmaxOld[i] = fVmax[i];
1565 }
1566 for (i = 0; i < fValSize; ++i) {
1567 fVmin[i] = DBL_MAX;
1568 fVmax[i] = - DBL_MAX;
1569 }
1570 //__________________________1-D histogram_______________________
1571 if (fAction == 1) {
1572 TH1 *h1 = (TH1*)fObject;
1573 if (h1->CanExtendAllAxes()) {
1574 for (i = 0; i < fNfill; i++) {
1575 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1576 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1577 }
1578 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h1, fVmin[0], fVmax[0]);
1579 }
1580 h1->FillN(fNfill, fVal[0], fW);
1581 //__________________________2-D histogram_______________________
1582 } else if (fAction == 2) {
1583 TH2 *h2 = (TH2*)fObject;
1584 if (h2->CanExtendAllAxes()) {
1585 for (i = 0; i < fNfill; i++) {
1586 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1587 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1588 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1589 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1590 }
1591 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2, fVmin[1], fVmax[1], fVmin[0], fVmax[0]);
1592 }
1593 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1594 //__________________________Profile histogram_______________________
1595 } else if (fAction == 4) {
1597 if (hp->CanExtendAllAxes()) {
1598 for (i = 0; i < fNfill; i++) {
1599 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1600 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1601 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1602 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1603 }
1604 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(hp, fVmin[1], fVmax[1]);
1605 }
1606 hp->FillN(fNfill, fVal[1], fVal[0], fW);
1607 //__________________________2D scatter plot_______________________
1608 } else if (fAction == 12) {
1609 TH2 *h2 = (TH2*)fObject;
1610 if (h2->CanExtendAllAxes()) {
1611 for (i = 0; i < fNfill; i++) {
1612 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1613 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1614 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1615 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1616 }
1617 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2, fVmin[1], fVmax[1], fVmin[0], fVmax[0]);
1618 // In case the new lower limits of h2 axis are 0, it is better to set them to the minimum of
1619 // the data set (which should be >0) to avoid data cut when plotting in log scale.
1620 TAxis *aX = h2->GetXaxis();
1621 TAxis *aY = h2->GetYaxis();
1622 Double_t xmin = aX->GetXmin();
1623 Double_t ymin = aY->GetXmin();
1624 if (xmin == 0 || ymin == 0) {
1625 if (aX->GetBinUpEdge(aX->FindFixBin(0.01*aX->GetBinWidth(aX->GetFirst()))) > fVmin[1]) xmin = fVmin[1];
1626 if (aY->GetBinUpEdge(aY->FindFixBin(0.01*aY->GetBinWidth(aY->GetFirst()))) > fVmin[0]) ymin = fVmin[0];
1627 h2->SetBins(aX->GetNbins(), xmin, aX->GetXmax(), aY->GetNbins(), ymin, aY->GetXmax());
1628 }
1629 }
1630
1631 if (!strstr(fOption.Data(), "same") && !strstr(fOption.Data(), "goff")) {
1632 if (!h2->TestBit(kCanDelete)) {
1633 // case like: T.Draw("y:x>>myhist")
1634 // we must draw a copy before filling the histogram h2=myhist
1635 // because h2 will be filled below and we do not want to show
1636 // the binned scatter-plot, the TGraph being better.
1637 TH1 *h2c = h2->DrawCopy(fOption.Data(),"");
1638 if (h2c) h2c->SetStats(false);
1639 } else {
1640 // case like: T.Draw("y:x")
1641 // h2 is a temporary histogram (htemp). This histogram
1642 // will be automatically deleted by TPad::Clear
1643 h2->Draw();
1644 }
1645 gPad->Update();
1646 }
1647 TGraph *pm = new TGraph(fNfill, fVal[1], fVal[0]);
1648 pm->SetEditable(false);
1649 pm->SetBit(kCanDelete);
1650 pm->SetMarkerStyle(fTree->GetMarkerStyle());
1651 pm->SetMarkerColor(fTree->GetMarkerColor());
1652 pm->SetMarkerSize(fTree->GetMarkerSize());
1653 pm->SetLineColor(fTree->GetLineColor());
1654 pm->SetLineWidth(fTree->GetLineWidth());
1655 pm->SetLineStyle(fTree->GetLineStyle());
1656 pm->SetFillColor(fTree->GetFillColor());
1657 pm->SetFillStyle(fTree->GetFillStyle());
1658 if (!fDraw && !strstr(fOption.Data(),"goff")) {
1659 if (fOption.Length() == 0 || strcasecmp(fOption.Data(),"same")==0) {
1660 pm->Draw("p");
1661 }
1662 else {
1663 TString opt = fOption;
1664 opt.ToLower();
1665 if (opt.Contains("a")) {
1666 TString temp(opt);
1667 temp.ReplaceAll("same","");
1668 if (temp.Contains("a")) {
1669 if (h2->TestBit(kCanDelete)) {
1670 // h2 will be deleted, the axis setting is delegated to only
1671 // the TGraph.
1672 h2 = nullptr;
1673 fObject = pm->GetHistogram();
1674 }
1675 }
1676 }
1677 pm->Draw(fOption.Data());
1678 }
1679 }
1680 if (h2 && !h2->TestBit(kCanDelete)) {
1681 for (i = 0; i < fNfill; i++) h2->Fill(fVal[1][i], fVal[0][i], fW[i]);
1682 }
1683 //__________________________3D scatter plot with option col_______________________
1684 } else if (fAction == 33) {
1685 TH2 *h2 = (TH2*)fObject;
1686 bool process2 = false;
1687 if (h2->CanExtendAllAxes()) {
1688 if (vminOld[2] == DBL_MAX)
1689 process2 = true;
1690 for (i = 0; i < fValSize && i < 4; i++) {
1691 fVmin[i] = vminOld[i];
1692 fVmax[i] = vmaxOld[i];
1693 }
1694 for (i = 0; i < fNfill; i++) {
1695 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1696 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1697 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1698 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1699 if (process2) {
1700 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1701 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1702 }
1703 }
1704 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h2, fVmin[1], fVmax[1], fVmin[0], fVmax[0]);
1705 // In case the new lower limits of h2 axis are 0, it is better to set them to the minimum of
1706 // the data set (which should be >0) to avoid data cut when plotting in log scale.
1707 TAxis *aX = h2->GetXaxis();
1708 TAxis *aY = h2->GetYaxis();
1709 Double_t xmin = aX->GetXmin();
1710 Double_t ymin = aY->GetXmin();
1711 if (xmin == 0 || ymin == 0) {
1712 if (aX->GetBinUpEdge(aX->FindFixBin(0.01*aX->GetBinWidth(aX->GetFirst()))) > fVmin[1]) xmin = fVmin[1];
1713 if (aY->GetBinUpEdge(aY->FindFixBin(0.01*aY->GetBinWidth(aY->GetFirst()))) > fVmin[0]) ymin = fVmin[0];
1714 h2->SetBins(aX->GetNbins(), xmin, aX->GetXmax(), aY->GetNbins(), ymin, aY->GetXmax());
1715 }
1716 } else {
1717 for (i = 0; i < fNfill; i++) {
1718 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1719 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1720 }
1721 }
1722 //__________________________3D scatter plot_______________________
1723 } else if (fAction == 3 || fAction == 13) {
1724 TH3 *h3 = (TH3*)fObject;
1725 if (h3->CanExtendAllAxes()) {
1726 for (i = 0; i < fNfill; i++) {
1727 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1728 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1729 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1730 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1731 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1732 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1733 }
1734 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h3, fVmin[2], fVmax[2], fVmin[1], fVmax[1], fVmin[0], fVmax[0]);
1735 }
1736 if (fAction == 3) {
1737 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1738 return;
1739 }
1740 if (!strstr(fOption.Data(), "same") && !strstr(fOption.Data(), "goff")) {
1741 if (!h3->TestBit(kCanDelete)) {
1742 // case like: T.Draw("y:x>>myhist")
1743 // we must draw a copy before filling the histogram h3=myhist
1744 // because h3 will be filled below and we do not want to show
1745 // the binned scatter-plot, the TGraph being better.
1746 TH1 *h3c = h3->DrawCopy(fOption.Data(),"");
1747 if (h3c) h3c->SetStats(false);
1748 } else {
1749 // case like: T.Draw("y:x")
1750 // h3 is a temporary histogram (htemp). This histogram
1751 // will be automatically deleted by TPad::Clear
1752 h3->Draw(fOption.Data());
1753 }
1754 gPad->Update();
1755 } else {
1756 rmin[0] = fVmin[2]; rmin[1] = fVmin[1]; rmin[2] = fVmin[0];
1757 rmax[0] = fVmax[2]; rmax[1] = fVmax[1]; rmax[2] = fVmax[0];
1758 gPad->Clear();
1759 gPad->Range(-1, -1, 1, 1);
1761 }
1763 pm3d->SetMarkerStyle(fTree->GetMarkerStyle());
1764 pm3d->SetMarkerColor(fTree->GetMarkerColor());
1765 pm3d->SetMarkerSize(fTree->GetMarkerSize());
1766 for (i = 0; i < fNfill; i++) {
1767 pm3d->SetPoint(i, fVal[2][i], fVal[1][i], fVal[0][i]);
1768 }
1769 if (!fDraw && !strstr(fOption.Data(), "goff")) pm3d->Draw();
1770 if (!h3->TestBit(kCanDelete)) {
1771 for (i = 0; i < fNfill; i++) h3->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1772 }
1773
1774 //__________________________2D Profile Histogram__________________
1775 } else if (fAction == 23) {
1777 if (hp->CanExtendAllAxes()) {
1778 for (i = 0; i < fNfill; i++) {
1779 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1780 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1781 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1782 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1783 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1784 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1785 }
1786 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(hp, fVmin[2], fVmax[2], fVmin[1], fVmax[1]);
1787 }
1788 for (i = 0; i < fNfill; i++) hp->Fill(fVal[2][i], fVal[1][i], fVal[0][i], fW[i]);
1789 //__________________________4D scatter plot_______________________
1790 } else if (fAction == 40) {
1791 TH3 *h3 = (TH3*)fObject;
1792 if (h3->CanExtendAllAxes()) {
1793 for (i = 0; i < fValSize && i < 4; i++) {
1794 fVmin[i] = vminOld[i];
1795 fVmax[i] = vmaxOld[i];
1796 }
1797 for (i = 0; i < fNfill; i++) {
1798 if (fVmin[0] > fVal[0][i]) fVmin[0] = fVal[0][i];
1799 if (fVmax[0] < fVal[0][i]) fVmax[0] = fVal[0][i];
1800 if (fVmin[1] > fVal[1][i]) fVmin[1] = fVal[1][i];
1801 if (fVmax[1] < fVal[1][i]) fVmax[1] = fVal[1][i];
1802 if (fVmin[2] > fVal[2][i]) fVmin[2] = fVal[2][i];
1803 if (fVmax[2] < fVal[2][i]) fVmax[2] = fVal[2][i];
1804 if (fVmin[3] > fVal[3][i]) fVmin[3] = fVal[3][i];
1805 if (fVmax[3] < fVal[3][i]) fVmax[3] = fVal[3][i];
1806 }
1807 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(h3, fVmin[2], fVmax[2], fVmin[1], fVmax[1], fVmin[0], fVmax[0]);
1808 } else {
1809 for (i = 0; i < fNfill; i++) {
1810 if (fVmin[3] > fVal[3][i]) fVmin[3] = fVal[3][i];
1811 if (fVmax[3] < fVal[3][i]) fVmax[3] = fVal[3][i];
1812 }
1813 }
1814 }
1815 //__________________________Parallel coordinates plot / candle chart_______________________
1816 else if (fAction == 6 || fAction == 7) {
1817 for (i = 0; i < fDimension; ++i) {
1818 for (Long64_t entry = 0; entry < fNfill; entry++) {
1819 if (fVmin[i] > fVal[i][entry]) fVmin[i] = fVal[i][entry];
1820 if (fVmax[i] < fVal[i][entry]) fVmax[i] = fVal[i][entry];
1821 }
1822 }
1823 }
1824}
1825
1826////////////////////////////////////////////////////////////////////////////////
1827/// Called at the end of a loop on a TTree.
1828
1830{
1831 // We take action (in Process) when we reach GetEstimate but
1832 // we reset at the beginning of Process, so
1833 // if fNfill == GetEstimate(), we just took action.
1834 if (fNfill && fNfill < fTree->GetEstimate())
1835 TakeAction();
1836
1837 if ((fSelectedRows == 0) && (TestBit(kCustomHistogram) == 0)) fDraw = 1; // do not draw
1838
1840}
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:411
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:31
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:32
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:40
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:37
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:36
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition TAttMarker.h:33
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:32
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition TAttMarker.h:34
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:46
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:490
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:735
<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:2290
Int_t GetN() const
Definition TGraph.h:131
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:927
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:879
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:8965
@ kAllAxes
Definition TH1.h:126
static TClass * Class()
virtual Int_t GetDimension() const
Definition TH1.h:528
@ kNoStats
Don't draw stats box.
Definition TH1.h:404
virtual Bool_t CanExtendAllAxes() const
Returns true if all axes are extendable.
Definition TH1.cxx:6671
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
Definition TH1.cxx:7138
TAxis * GetXaxis()
Definition TH1.h:572
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:6684
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:3418
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition TH1.cxx:9048
2-D histogram with a double per channel (see TH1 documentation)
Definition TH2.h:356
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:371
3-D histogram with a float per channel (see TH1 documentation)
Definition TH3.h:326
The 3-D histogram classes derived from the 1-D histogram classes.
Definition TH3.h:40
static THLimitsFinder * GetLimitsFinder()
Return pointer to the current finder.
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:575
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:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:202
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1057
virtual void Delete(Option_t *option="")
Delete this object.
Definition TObject.cxx:267
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:864
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:543
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
virtual const char * GetTitle() const
Returns title of object.
Definition TObject.cxx:501
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
void ResetBit(UInt_t f)
Definition TObject.h:201
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:68
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:712
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:640
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:546
virtual Double_t GetWeight() const
Definition TTree.h:623
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:9317
virtual TTree * GetTree() const
Definition TTree.h:596
virtual void SetEntryList(TEntryList *list, Option_t *opt="")
Set an EntryList.
Definition TTree.cxx:9253
TEventList * GetEventList() const
Definition TTree.h:552
@ kForceRead
Definition TTree.h:287
virtual Int_t GetUpdate() const
Definition TTree.h:600
virtual Long64_t GetChainOffset() const
Definition TTree.h:495
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