Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
TSpider.cxx
Go to the documentation of this file.
1// @(#)root/treeviewer:$Id: 2bb6def14de16b049d4c979e73ad2b08d0936520 $
2// Author: Bastien Dalla Piazza 20/07/07
3
4/*************************************************************************
5 * Copyright (C) 1995-2007, 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#include "TSpider.h"
13#include "TAttFill.h"
14#include "TAttText.h"
15#include "TAttLine.h"
16#include "TGraphPolargram.h"
17#include "TPolyLine.h"
18#include "TNtuple.h"
19#include "TBranch.h"
20#include "TTreeFormula.h"
21#include "TTreeFormulaManager.h"
22#include "TList.h"
23#include "TSelectorDraw.h"
24#include "TROOT.h"
25#include "TEntryList.h"
26#include "TLatex.h"
27#include "TVirtualPad.h"
28#include "TMath.h"
29#include "TCanvas.h"
30#include "TArc.h"
31#include <cfloat>
32#include "TEnv.h"
33
34
35/** \class TSpider
36Spider class.
37
38A spider view is a handy way to visualize a set of data stored in a TTree.
39It draws as many polar axes as selected data members. For each of them, it draws
40on the axis the position of the present event between the min and max of the
41data member. Two modes are available:
42
43 - The spider view: With each points on the axes is drawn a polyline.
44 - The segment view: For each data member is drawn an arc segment with the
45 radius corresponding to the event.
46
47The spider plot is available from the treeviewer called by
48"atree->StartViewer()", or simply by calling its constructor and defining the
49variables to display.
50
51Begin_Macro(source)
52{
53 TCanvas *c1 = new TCanvas("c1","TSpider example",200,10,700,700);
54 TFile *f = new TFile("$(ROOTSYS)/tutorials/hsimple.root");
55 if (!f || f->IsZombie()) {
56 printf("Please run <ROOT location>/tutorials/hsimple.C before.");
57 return;
58 }
59 TNtuple* ntuple = (TNtuple*)f->Get("ntuple");
60 TString varexp = "px:py:pz:random:sin(px):log(px/py):log(pz)";
61 TString selectStr = "px>0 && py>0 && pz>0";
62 TString options = "average";
63 TSpider *spider = new TSpider(ntuple,varexp.Data(),selectStr.Data(),options.Data());
64 spider->Draw();
65 c1->ToggleEditor();
66 c1->Selected(c1,spider,1);
67 return c1;
68}
69End_Macro
70*/
71
72////////////////////////////////////////////////////////////////////////////////
73/// Default constructor.
74
76{
77 fDisplayAverage=false;
78 fForceDim=false;
79 fPolargram=nullptr;
80 fInput=nullptr;
81 fManager=nullptr;
82 fNcols=0;
83 fNx=3;
84 fNy=4;
85 fPolyList=nullptr;
86 fSelect=nullptr;
87 fSelector=nullptr;
88 fTree=nullptr;
89 fMax=nullptr;
90 fMin=nullptr;
91 fAve=nullptr;
92 fCanvas=nullptr;
93 fAveragePoly=nullptr;
94 fEntry=0;
95 fSuperposed=nullptr;
96 fShowRange=false;
97 fAngularLabels=false;
98 fAverageSlices=nullptr;
99 fSegmentDisplay=false;
100 fNentries=0;
101 fFirstEntry=0;
102 fArraySize=0;
103 fCurrentEntries = nullptr;
104 fFormulas = nullptr;
105}
106
107////////////////////////////////////////////////////////////////////////////////
108/// Normal constructor. Options are:
109/// - "average"
110/// - "showrange"
111/// - "segment"
112
113TSpider::TSpider(TTree* tree ,const char *varexp, const char *selection,
114 Option_t *option, Long64_t nentries, Long64_t firstentry)
115 : TAttFill(2,3003), TAttLine(1,1,1)
116{
117 UInt_t ui=0;
118
119 fArraySize = 16;
120 fTree=tree;
122 fFormulas= new TList();
123 fInput= new TList();
124 fInput->Add(new TNamed("varexp",""));
125 fInput->Add(new TNamed("selection",""));
126 fSelector->SetInputList(fInput);
127 gROOT->GetListOfCleanups()->Add(this);
128 fNx=2;
129 fNy=2;
130 fDisplayAverage=false;
131 fSelect=nullptr;
132 fManager=nullptr;
133 fCanvas=nullptr;
134 fAveragePoly=nullptr;
136 fSuperposed=nullptr;
137 fShowRange=false;
138 fAngularLabels=true;
139 fForceDim=false;
140 fAverageSlices=nullptr;
141 fSegmentDisplay=false;
142 if (firstentry < 0 || firstentry > tree->GetEstimate()) firstentry = 0;
143 fFirstEntry = firstentry;
144 if (nentries>0) fNentries = nentries;
145 else fNentries = nentries = tree->GetEstimate()-firstentry;
146
148
149 fPolargram=nullptr;
150 fPolyList=nullptr;
151
152 fTree->SetScanField(fNx*fNy);
154 for(ui=0;ui<fNx*fNy;++ui) fCurrentEntries[ui]=0;
155
156 TString opt = option;
157
158 if (opt.Contains("average")) fDisplayAverage=true;
159 if (opt.Contains("showrange")) fShowRange=true;
160 if (opt.Contains("segment")) fSegmentDisplay=true;
161
162 fNcols=8;
163
165 SetSelectionExpression(selection);
166 SyncFormulas();
167 InitVariables(firstentry,nentries);
168}
169
170////////////////////////////////////////////////////////////////////////////////
171/// Destructor.
172
174{
175 delete [] fCurrentEntries;
176 if(fPolyList){
177 fPolyList->Delete();
178 delete fPolyList;
179 }
181 delete [] fAverageSlices;
182 if(fFormulas){
183 fFormulas->Delete();
184 delete fFormulas;
185 }
186 if(fSelect) delete fSelect;
187 if(fSelector) delete fSelector;
188 if(fInput){
189 fInput->Delete();
190 delete fInput;
191 }
192 if(fMax) delete [] fMax;
193 if(fMin) delete [] fMin;
194 if(fAve) delete [] fAve;
195 if(fSuperposed){
196 fSuperposed->Delete();
197 delete fSuperposed;
198 }
199 if (fCanvas) fCanvas->cd(0);
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Allow to superpose several spider views.
204
206{
207 if(!fSuperposed) fSuperposed=new TList();
208 fSuperposed->Add(sp);
209}
210
211////////////////////////////////////////////////////////////////////////////////
212/// Add a variable to the plot from its expression.
213
214void TSpider::AddVariable(const char* varexp)
215{
216 if(!varexp[0]) return;
217 TTreeFormula *fvar = new TTreeFormula("Var1",varexp,fTree);
218 if(fvar->GetNdim() <= 0) return;
219
220 fFormulas->AddAfter(fFormulas->At(fNcols-1),fvar);
221
222 InitArrays(fNcols + 1);
223 ++fNcols;
224 SyncFormulas();
225
226 UInt_t ui=0;
227 Long64_t notSkipped=0;
228 Int_t tnumber=-1;
229 Long64_t entryNumber;
230 Long64_t entry = fFirstEntry;
231 Int_t entriesToDisplay = fNentries;
232 while(entriesToDisplay!=0){
233 entryNumber = fTree->GetEntryNumber(entry);
234 if(entryNumber < 0) break;
235 Long64_t localEntry = fTree->LoadTree(entryNumber);
236 if(localEntry < 0) break;
237 if(tnumber != fTree->GetTreeNumber()) {
238 tnumber = fTree->GetTreeNumber();
239 if(fManager) fManager->UpdateFormulaLeaves();
240 else {
241 for(Int_t i=0;i<=fFormulas->LastIndex();++i)
242 ((TTreeFormula*)fFormulas->At(i))->UpdateFormulaLeaves();
243 }
244 }
245 Int_t ndata=1;
246 if(fForceDim){
247 if(fManager)
248 ndata = fManager->GetNdata(true);
249 else {
250 for(ui=0;ui<fNcols;++ui){
251 if(ndata<((TTreeFormula*)fFormulas->At(ui))->GetNdata())
252 ndata = ((TTreeFormula*)fFormulas->At(ui))->GetNdata();
253 }
254 if(fSelect && fSelect->GetNdata() == 0)
255 ndata = 0;
256 }
257 }
258
259 bool loaded = false;
260 bool skip = false;
261 // Loop over the instances of the selection condition
262 for(Int_t inst=0;inst<ndata;++inst){
263 if(fSelect){
264 if(fSelect->EvalInstance(inst) == 0){
265 skip = true;
266 ++entry;
267 }
268 }
269 if (!loaded) {
270 // EvalInstance(0) always needs to be called so that
271 // the proper branches are loaded.
272 ((TTreeFormula*)fFormulas->At(fNcols-1))->EvalInstance(0);
273 loaded = true;
274 } else if (inst == 0) {
275 loaded = true;
276 }
277 }
278 if(!skip){
279 fTree->LoadTree(entryNumber);
281 if(var->EvalInstance()>fMax[fNcols-1]) fMax[fNcols-1]=var->EvalInstance();
282 if(var->EvalInstance()<fMin[fNcols-1]) fMin[fNcols-1]=var->EvalInstance();
283 fAve[fNcols-1]+=var->EvalInstance();
284 ++notSkipped;
285 --entriesToDisplay;
286 ++entry;
287 }
288 }
289 if (notSkipped) fAve[fNcols-1]/=notSkipped;
290
291 Color_t lc;
292 Style_t lt;
293 Width_t lw;
294 Color_t fc;
295 Style_t fs;
296
297 if(fAverageSlices){
298 lc = fAverageSlices[0]->GetLineColor();
299 lt = fAverageSlices[0]->GetLineStyle();
300 lw = fAverageSlices[0]->GetLineWidth();
301 fc = fAverageSlices[0]->GetFillColor();
302 fs = fAverageSlices[0]->GetFillStyle();
303 } else {
304 lc = fAveragePoly->GetLineColor();
305 lt = fAveragePoly->GetLineStyle();
306 lw = fAveragePoly->GetLineWidth();
307 fc = fAveragePoly->GetFillColor();
308 fs = fAveragePoly->GetFillStyle();
309 }
310
311 delete fPolargram;
312 fPolargram = nullptr;
313
314 if(fSegmentDisplay){
315 for(ui=0;ui<fNx*fNy;++ui) ((TList*)fPolyList->At(ui))->Delete();
316 if (fAverageSlices) for(ui=0;ui<fNcols-1;++ui) delete fAverageSlices[ui];
317 }
318 fPolyList->Delete();
319 delete fPolyList;
320 fPolyList = nullptr;
321 delete [] fAverageSlices;
322 fAverageSlices = nullptr;
323 delete fAveragePoly;
324 fAveragePoly = nullptr;
325
326 if (fCanvas) {
327 fCanvas->Clear();
328 fCanvas->Divide(fNx,fNy);
329 }
330 Draw("");
331
332 if(fAverageSlices){
333 for(ui = 0;ui<fNcols;++ui){
334 fAverageSlices[ui]->SetLineColor(lc);
335 fAverageSlices[ui]->SetLineStyle(lt);
336 fAverageSlices[ui]->SetLineWidth(lw);
337 fAverageSlices[ui]->SetFillColor(fc);
338 fAverageSlices[ui]->SetFillStyle(fs);
339 }
340 } else {
341 fAveragePoly->SetLineColor(lc);
342 fAveragePoly->SetLineStyle(lt);
343 fAveragePoly->SetLineWidth(lw);
344 fAveragePoly->SetFillColor(fc);
345 fAveragePoly->SetFillStyle(fs);
346 }
347}
348
349////////////////////////////////////////////////////////////////////////////////
350/// Delete a variable from its expression.
351
352void TSpider::DeleteVariable(const char* varexp)
353{
354 Int_t var=-1;
355 UInt_t ui=0;
356
357 if(fNcols == 2) return;
358 for(ui=0; ui<fNcols;++ui){
359 if(!strcmp(varexp,((TTreeFormula*)fFormulas->At(ui))->GetTitle())) var = ui;
360 }
361 if(var<0) return;
362
363 fFormulas->Remove(fFormulas->At(var));
364 SyncFormulas();
365
366 for(ui=var+1;ui<fNcols;++ui){
367 fMin[ui-1] = fMin[ui];
368 fMax[ui-1] = fMax[ui];
369 fAve[ui-1] = fAve[ui];
370 }
371 fMin[fNcols-1] = DBL_MAX;
372 fMax[fNcols-1] = -DBL_MAX;
373 fAve[fNcols-1] = 0;
374 --fNcols;
375
376 Color_t lc;
377 Style_t lt;
378 Width_t lw;
379 Color_t fc;
380 Style_t fs;
381
382 if(fAverageSlices){
383 lc = fAverageSlices[0]->GetLineColor();
384 lt = fAverageSlices[0]->GetLineStyle();
385 lw = fAverageSlices[0]->GetLineWidth();
386 fc = fAverageSlices[0]->GetFillColor();
387 fs = fAverageSlices[0]->GetFillStyle();
388 } else {
389 lc = fAveragePoly->GetLineColor();
390 lt = fAveragePoly->GetLineStyle();
391 lw = fAveragePoly->GetLineWidth();
392 fc = fAveragePoly->GetFillColor();
393 fs = fAveragePoly->GetFillStyle();
394 }
395
396 delete fPolargram;
397 fPolargram = nullptr;
398
399 if(fSegmentDisplay){
400 for(ui=0;ui<fNx*fNy;++ui) ((TList*)fPolyList->At(ui))->Delete();
401 if (fAverageSlices) for(ui=0;ui<=fNcols;++ui) delete fAverageSlices[ui];
402 }
403 fPolyList->Delete();
404 delete fPolyList;
405 fPolyList = nullptr;
406 delete [] fAverageSlices;
407 fAverageSlices = nullptr;
408 delete fAveragePoly;
409 fAveragePoly = nullptr;
410
411 if (fCanvas) {
412 fCanvas->Clear();
413 fCanvas->Divide(fNx,fNy);
414 }
415 Draw("");
416 if(fNcols == 2) SetSegmentDisplay(true);
417
418 if(fAverageSlices){
419 for(ui = 0;ui<fNcols;++ui){
420 fAverageSlices[ui]->SetLineColor(lc);
421 fAverageSlices[ui]->SetLineStyle(lt);
422 fAverageSlices[ui]->SetLineWidth(lw);
423 fAverageSlices[ui]->SetFillColor(fc);
424 fAverageSlices[ui]->SetFillStyle(fs);
425 }
426 } else {
427 fAveragePoly->SetLineColor(lc);
428 fAveragePoly->SetLineStyle(lt);
429 fAveragePoly->SetLineWidth(lw);
430 fAveragePoly->SetFillColor(fc);
431 fAveragePoly->SetFillStyle(fs);
432 }
433}
434
435////////////////////////////////////////////////////////////////////////////////
436/// Compute the distance to the spider.
437
439{
440 if(!gPad) return 9999;
441 Double_t xx,yy,r2;
442 xx=gPad->AbsPixeltoX(px);
443 yy=gPad->AbsPixeltoY(py);
444 r2 = xx*xx + yy*yy;
445 if(r2>1 && r2<1.5)
446 return 0;
447 else return 9999;
448}
449
450////////////////////////////////////////////////////////////////////////////////
451/// Draw the spider.
452
454{
455 UInt_t ui=0;
456
457 gEnv->SetValue("Canvas.ShowEditor",1);
458 if(!gPad && !fCanvas){
459 fCanvas = new TCanvas("screen","Spider Plot",fNx*256,fNy*256);
460 if (fCanvas) fCanvas->Divide(fNx,fNy);
461 } else if(!fCanvas){
462 fCanvas = (TCanvas*)gPad;
463 if (fCanvas) fCanvas->Divide(fNx,fNy);
464 }
465 if(fPolargram) delete fPolargram;
466 fPolargram=new TGraphPolargram("fPolargram");
467 fPolargram->SetNdivPolar(fNcols);
468 fPolargram->SetNdivRadial(0);
469 if (fCanvas) fCanvas->cd();
471 AppendPad(options);
472 for(ui=0;ui<fNx*fNy;++ui){
473 if (fCanvas) fCanvas->cd(ui+1);
474 fPolargram->Draw("pn");
475 fTree->LoadTree(fCurrentEntries[ui]);
476 if(fSegmentDisplay){
478 DrawSlices("");
479 } else {
481 DrawPoly("");
482 }
483 AppendPad();
484 }
485 if (fCanvas) fCanvas->Selected(fCanvas,this,1);
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// Paint the Polygon representing the average value of the variables.
490
492{
493 Int_t linecolor=4;
494 Int_t fillstyle=0;
495 Int_t fillcolor=linecolor;
496 Int_t linewidth=1;
497 Int_t linestyle=1;
498
499 UInt_t ui=0;
500 Double_t slice = 2*TMath::Pi()/fNcols;
501 Double_t *x = new Double_t[fNcols+1];
502 Double_t *y = new Double_t[fNcols+1];
503
504 for(ui=0;ui<fNcols;++ui){
505 x[ui]=(fAve[ui]-fMin[ui])/(fMax[ui]-fMin[ui])*TMath::Cos(ui*slice);
506 y[ui]=(fAve[ui]-fMin[ui])/(fMax[ui]-fMin[ui])*TMath::Sin(ui*slice);
507 }
508 x[fNcols]=(fAve[0]-fMin[0])/(fMax[0]-fMin[0]);
509 y[fNcols]=0;
510
511 if(!fAveragePoly){
513 fAveragePoly->SetLineColor(linecolor);
514 fAveragePoly->SetLineWidth(linewidth);
515 fAveragePoly->SetLineStyle(linestyle);
516 fAveragePoly->SetFillStyle(fillstyle);
517 fAveragePoly->SetFillColor(fillcolor);
518 }
519 fAveragePoly->Draw();
520 fAveragePoly->Draw("f");
521
522 delete [] x;
523 delete [] y;
524}
525
526////////////////////////////////////////////////////////////////////////////////
527/// Paint the polygon representing the current entry.
528
529void TSpider::DrawPoly(Option_t* /*options*/)
530{
531 if(!fPolyList) fPolyList = new TList();
532 Double_t *x = new Double_t[fNcols+1];
533 Double_t *y = new Double_t[fNcols+1];
534
535 Double_t slice = 2*TMath::Pi()/fNcols;
536 for(UInt_t i=0;i<fNcols;++i){
537 x[i]=(((TTreeFormula*)fFormulas->At(i))->EvalInstance()-fMin[i])/(fMax[i]-fMin[i])*TMath::Cos(i*slice);
538 y[i]=(((TTreeFormula*)fFormulas->At(i))->EvalInstance()-fMin[i])/(fMax[i]-fMin[i])*TMath::Sin(i*slice);
539 }
540 x[fNcols]=(((TTreeFormula*)fFormulas->At(0))->EvalInstance()-fMin[0])/(fMax[0]-fMin[0]);
541 y[fNcols]=0;
542
543 TPolyLine* poly= new TPolyLine(fNcols+1,x,y);
544 poly->SetFillColor(GetFillColor());
545 poly->SetFillStyle(GetFillStyle());
546 poly->SetLineWidth(GetLineWidth());
547 poly->SetLineColor(GetLineColor());
548 poly->SetLineStyle(GetLineStyle());
549 poly->Draw("f");
550 poly->Draw();
551 fPolyList->Add(poly);
552 delete [] x;
553 delete [] y;
554}
555
556////////////////////////////////////////////////////////////////////////////////
557/// Draw the slices of the segment plot.
558
560{
561 UInt_t ui=0;
562
563 Double_t angle = 2*TMath::Pi()/fNcols;
564 Double_t conv = 180.0/TMath::Pi();
565
566 if(!fPolyList) fPolyList = new TList;
567 TList* li = new TList();
568 for(ui=0;ui<fNcols;++ui){
569 Double_t r = (((TTreeFormula*)fFormulas->At(ui))->EvalInstance()-fMin[ui])/(fMax[ui]-fMin[ui]);
570 TArc* slice = new TArc(0,0,r,(ui-0.25)*angle*conv,(ui+0.25)*angle*conv);
571 slice->SetFillColor(GetFillColor());
572 slice->SetFillStyle(GetFillStyle());
573 slice->SetLineWidth(GetLineWidth());
574 slice->SetLineColor(GetLineColor());
575 slice->SetLineStyle(GetLineStyle());
576 li->Add(slice);
577 slice->Draw(options);
578 }
579 fPolyList->Add(li);
580}
581
582////////////////////////////////////////////////////////////////////////////////
583/// Draw the slices representing the average for the segment plot.
584
586{
587 UInt_t ui=0;
588
589 Int_t fillstyle=3002;
590 Int_t linecolor=4;
591 Int_t fillcolor=linecolor;
592 Int_t linewidth=1;
593 Int_t linestyle=1;
594
595 Double_t angle = 2*TMath::Pi()/fNcols;
596 Double_t conv = 180.0/TMath::Pi();
597
598 if(!fAverageSlices){
599 fAverageSlices = new TArc*[fNcols];
600 for(ui=0;ui<fNcols;++ui){
601 Double_t r = (fAve[ui]-fMin[ui])/(fMax[ui]-fMin[ui]);
602 fAverageSlices[ui] = new TArc(0,0,r,(ui-0.5)*angle*conv,(ui+0.5)*angle*conv);
603 fAverageSlices[ui]->SetFillColor(fillcolor);
604 fAverageSlices[ui]->SetFillStyle(fillstyle);
605 fAverageSlices[ui]->SetLineWidth(linewidth);
606 fAverageSlices[ui]->SetLineColor(linecolor);
607 fAverageSlices[ui]->SetLineStyle(linestyle);
608 }
609 }
610 for(ui=0;ui<fNcols;++ui) fAverageSlices[ui]->Draw();
611}
612
613////////////////////////////////////////////////////////////////////////////////
614/// Get the LineStyle of the average.
615
617{
618 if(fAverageSlices) return fAverageSlices[0]->GetLineStyle();
619 else if(fAveragePoly) return fAveragePoly->GetLineStyle();
620 else return 0;
621}
622
623////////////////////////////////////////////////////////////////////////////////
624/// Get the LineColor of the average.
625
627{
628 if(fAverageSlices) return fAverageSlices[0]->GetLineColor();
629 else if(fAveragePoly) return fAveragePoly->GetLineColor();
630 else return 0;
631}
632
633////////////////////////////////////////////////////////////////////////////////
634/// Get the LineWidth of the average.
635
637{
638 if(fAverageSlices) return fAverageSlices[0]->GetLineWidth();
639 else if(fAveragePoly) return fAveragePoly->GetLineWidth();
640 else return 0;
641}
642
643////////////////////////////////////////////////////////////////////////////////
644/// Get the Fill Color of the average.
645
647{
648 if(fAverageSlices) return fAverageSlices[0]->GetFillColor();
649 else if(fAveragePoly) return fAveragePoly->GetFillColor();
650 else return 0;
651}
652
653////////////////////////////////////////////////////////////////////////////////
654/// Get the FillStyle of the average.
655
657{
658 if(fAverageSlices) return fAverageSlices[0]->GetFillStyle();
659 else if(fAveragePoly) return fAveragePoly->GetFillStyle();
660 else return 0;
661}
662
663////////////////////////////////////////////////////////////////////////////////
664/// Execute the corresponding event.
665
666void TSpider::ExecuteEvent(Int_t /*event*/,Int_t /*px*/, Int_t /*py*/)
667{
668 if (!gPad) return;
669 gPad->SetCursor(kHand);
670}
671
672////////////////////////////////////////////////////////////////////////////////
673/// Find the alignement rule to apply for TText::SetTextAlign(Short_t).
674
676{
677 Double_t pi = TMath::Pi();
678
679 while(angle < 0 || angle > 2*pi){
680 if(angle < 0) angle+=2*pi;
681 if(angle > 2*pi) angle-=2*pi;
682 }
683 if(!fAngularLabels){
684 if(angle > 0 && angle < pi/2) return 11;
685 else if(angle > pi/2 && angle < pi) return 31;
686 else if(angle > pi && angle < 3*pi/2) return 33;
687 else if(angle > 3*pi/2 && angle < 2*pi) return 13;
688 else if(angle == 0 || angle == 2*pi) return 12;
689 else if(angle == pi/2) return 21;
690 else if(angle == pi) return 32;
691 else if(angle == 3*pi/2) return 23;
692 else return 0;
693 }
694 else{
695 if(angle >= 0 && angle < pi) return 21;
696 else if(angle >=pi && angle <= 2*pi) return 23;
697 else return 0;
698 }
699}
700
701////////////////////////////////////////////////////////////////////////////////
702/// Determine the orientation of the polar labels according to their angle.
703
705{
706 Double_t pi = TMath::Pi();
707 Double_t convraddeg = 180.0/pi;
708
709 while(angle < 0 || angle > 2*pi){
710 if(angle < 0) angle+=2*pi;
711 if(angle > 2*pi) angle-=2*pi;
712 }
713
714 if(angle >= 0 && angle <= pi/2) return angle*convraddeg - 90;
715 else if(angle > pi/2 && angle <= pi) return (angle + pi)*convraddeg + 90;
716 else if(angle > pi && angle <= 3*pi/2) return (angle - pi)*convraddeg - 90;
717 else if(angle > 3*pi/2 && angle <= 2*pi) return angle*convraddeg + 90;
718 else return 0;
719}
720
721////////////////////////////////////////////////////////////////////////////////
722/// return the number of entries to be processed
723/// this function checks that nentries is not bigger than the number
724/// of entries in the Tree or in the associated TEventlist
725
727{
728 Long64_t lastentry = firstentry + nentries - 1;
729 if (lastentry > fTree->GetEntriesFriend()-1) {
730 lastentry = fTree->GetEntriesFriend() - 1;
731 nentries = lastentry - firstentry + 1;
732 }
733 //TEventList *elist = fTree->GetEventList();
734 //if (elist && elist->GetN() < nentries) nentries = elist->GetN();
735 TEntryList *elist = fTree->GetEntryList();
736 if (elist && elist->GetN() < nentries) nentries = elist->GetN();
737 return nentries;
738}
739
740////////////////////////////////////////////////////////////////////////////////
741/// Go to a specified entry.
742
744{
745 if(e<fFirstEntry || e+fTree->GetScanField()>=fFirstEntry + fNentries) return;
746 fEntry = e;
748}
749
750////////////////////////////////////////////////////////////////////////////////
751/// Go to the next entries.
752
754{
755 if(fEntry + 2*fTree->GetScanField() -1 >= fFirstEntry + fNentries) fEntry = fFirstEntry;
756 else fEntry=fCurrentEntries[fTree->GetScanField()-1]+1;
758}
759
760////////////////////////////////////////////////////////////////////////////////
761/// Go to the previous entries.
762
764{
765 if(fEntry-fTree->GetScanField() < fFirstEntry) fEntry = fFirstEntry + fNentries -1 - fTree->GetScanField();
766 else fEntry -= fTree->GetScanField();
768}
769
770////////////////////////////////////////////////////////////////////////////////
771/// Go to the next entry.
772
774{
775 if(fEntry + fTree->GetScanField() >= fFirstEntry + fNentries) return;
776 ++fEntry;
778}
779
780////////////////////////////////////////////////////////////////////////////////
781/// Go to the last entry.
782
784{
785 if(fEntry - 1 < fFirstEntry) return;
786 --fEntry;
788}
789
790////////////////////////////////////////////////////////////////////////////////
791/// Check if the arrays size is enough and reallocate them if not.
792
794{
795 if(newsize>fArraySize){
796
797 Int_t i;
798 Int_t old = fArraySize;
799
800 while(fArraySize<newsize) fArraySize*=2;
801
802 Double_t *memmax = new Double_t[fArraySize];
803 Double_t *memmin = new Double_t[fArraySize];
804 Double_t *memave = new Double_t[fArraySize];
805
806 for(i=0;i<fArraySize;++i){
807 if(i<old){
808 memmax[i] = fMax[i];
809 memmin[i] = fMin[i];
810 memave[i] = fAve[i];
811 } else {
812 memmax[i] = -DBL_MAX;
813 memmin[i] = DBL_MAX;
814 memave[i] = 0;
815 }
816 }
817
818 delete [] fMax;
819 delete [] fMin;
820 delete [] fAve;
821
822 fMax = memmax;
823 fMin = memmin;
824 fAve = memave;
825 }
826}
827
828////////////////////////////////////////////////////////////////////////////////
829/// Browse the tree to set the min, max and average value of each variable of fVar.
830
832{
833 UInt_t ui=0;
834 Int_t i;
835
836 fMax = new Double_t [fArraySize];
837 fMin= new Double_t [fArraySize];
838 fAve= new Double_t [fArraySize];
839
840 for(i=0;i<fArraySize;++i){
841 fMax[i]= -DBL_MAX;
842 fMin[i]= DBL_MAX;
843 fAve[i]=0;
844 }
845
846 Long64_t notSkipped=0;
847 Int_t tnumber=-1;
848 Long64_t entryNumber;
849 Long64_t entry = firstentry;
850 Int_t entriesToDisplay = nentries;
851 while(entriesToDisplay!=0){
852 entryNumber = fTree->GetEntryNumber(entry);
853 if(entryNumber < 0) break;
854 Long64_t localEntry = fTree->LoadTree(entryNumber);
855 if(localEntry < 0) break;
856 if(tnumber != fTree->GetTreeNumber()) {
857 tnumber = fTree->GetTreeNumber();
858 if(fManager) fManager->UpdateFormulaLeaves();
859 else {
860 for(i=0;i<=fFormulas->LastIndex();++i)
861 ((TTreeFormula*)fFormulas->At(i))->UpdateFormulaLeaves();
862 }
863 }
864 Int_t ndata=1;
865 if(fForceDim){
866 if(fManager)
867 ndata = fManager->GetNdata(true);
868 else {
869 for(ui=0;ui<fNcols;++ui){
870 if(ndata<((TTreeFormula*)fFormulas->At(ui))->GetNdata())
871 ndata = ((TTreeFormula*)fFormulas->At(ui))->GetNdata();
872 }
873 if(fSelect && fSelect->GetNdata() == 0)
874 ndata = 0;
875 }
876 }
877 bool loaded = false;
878 bool skip = false;
879 // Loop over the instances of the selection condition
880 for(Int_t inst=0;inst<ndata;++inst){
881 if(fSelect){
882 if(fSelect->EvalInstance(inst) == 0){
883 skip = true;
884 ++entry;
885 }
886 }
887 if (!loaded) {
888 // EvalInstance(0) always needs to be called so that
889 // the proper branches are loaded.
890 for (ui=0;ui<fNcols;ui++) {
891 ((TTreeFormula*)fFormulas->At(ui))->EvalInstance(0);
892 }
893 loaded = true;
894 } else if (inst == 0) {
895 loaded = true;
896 }
897 }
898 if(!skip){
899 fTree->LoadTree(entryNumber);
900 for(ui=0;ui<fNcols;++ui){
901 Double_t inst = ((TTreeFormula*)fFormulas->At(ui))->EvalInstance();
902 if(inst > fMax[ui]) fMax[ui] = inst;
903 if(inst < fMin[ui]) fMin[ui] = inst;
904 fAve[ui] += inst;
905 }
906 ++notSkipped;
907 --entriesToDisplay;
908 ++entry;
909 }
910 }
911 if (notSkipped) {for(ui=0;ui<fNcols;++ui) fAve[ui]/=notSkipped;}
912}
913
914////////////////////////////////////////////////////////////////////////////////
915/// Paint the spider.
916
918{
919 UInt_t ui=0;
920 TString opt = options;
921
922 if(opt.Contains("n")) return;
923
924 Double_t slice = 2*TMath::Pi()/fNcols;
925 Double_t offset(1.0);
926 if (!fCanvas) {
927 if (gPad) fCanvas = gPad->GetCanvas();
928 else return;
929 }
930
931 TLatex *txt = new TLatex();
932 for(ui=0;ui<fNx*fNy;++ui){
933 txt->SetTextAlign(13);
934 if (fCanvas) fCanvas->cd(ui+1);
935 if (fCurrentEntries) {
936 txt->PaintLatex(-1.2,1.2,0,0.08,Form("#%d",(int)fCurrentEntries[ui]));
937 }
938 txt->SetTextSize(0.035);
939 for(UInt_t var=0;var<fNcols;++var){ // Print labels.
940 if(ui==0){
941 txt->SetTextAlign(FindTextAlign(var*slice));
942 offset = 1.09 + txt->GetTextSize();
943 txt->PaintLatex(offset*TMath::Cos(var*slice),offset*TMath::Sin(var*slice),
944 FindTextAngle(var*slice),0.035,fFormulas->At(var)->GetTitle());
945 offset= 1.03;
946 txt->PaintLatex(offset*TMath::Cos(var*slice),offset*TMath::Sin(var*slice),
947 FindTextAngle(var*slice),0.035,Form("[%5.3f,%5.3f]",fMin[var],fMax[var]));
948 }
949 else {
950 txt->SetTextAlign(FindTextAlign(var*slice));
951 if(var*slice >=0 && var*slice <= TMath::Pi()) offset =1.13 + txt->GetTextSize();
952 else offset = 1.09 + txt->GetTextSize();
953 txt->PaintLatex(offset*TMath::Cos(var*slice),offset*TMath::Sin(var*slice),
954 FindTextAngle(var*slice),0.035,fFormulas->At(var)->GetTitle());
955 }
956 }
957 }
958 delete txt;
959}
960
961////////////////////////////////////////////////////////////////////////////////
962/// Set the LineStyle of the average.
963
965{
966 UInt_t ui=0;
967
968 if(fAverageSlices){
969 for(ui=0;ui<fNcols;++ui) fAverageSlices[ui]->SetLineStyle(sty);
970 } else if(fAveragePoly) fAveragePoly->SetLineStyle(sty);
971}
972
973////////////////////////////////////////////////////////////////////////////////
974/// Set the LineColor of the average.
975
977{
978 UInt_t ui=0;
979
980 if(fAverageSlices){
981 for(ui=0;ui<fNcols;++ui) fAverageSlices[ui]->SetLineColor(col);
982 } else if(fAveragePoly) fAveragePoly->SetLineColor(col);
983}
984
985////////////////////////////////////////////////////////////////////////////////
986/// Set the LineWidth of the average.
987
989{
990 UInt_t ui=0;
991
992 if(fAverageSlices){
993 for(ui=0;ui<fNcols;++ui) fAverageSlices[ui]->SetLineWidth(wid);
994 } else if(fAveragePoly) fAveragePoly->SetLineWidth(wid);
995}
996
997////////////////////////////////////////////////////////////////////////////////
998/// Set the Fill Color of the average.
999
1001{
1002 UInt_t ui=0;
1003
1004 if(fAverageSlices){
1005 for(ui=0;ui<fNcols;++ui) fAverageSlices[ui]->SetFillColor(col);
1006 } else if(fAveragePoly) fAveragePoly->SetFillColor(col);
1007}
1008
1009////////////////////////////////////////////////////////////////////////////////
1010/// Set the FillStyle of the average.
1011
1013{
1014 UInt_t ui=0;
1015
1016 if(fAverageSlices){
1017 for(ui=0;ui<fNcols;++ui) fAverageSlices[ui]->SetFillStyle(sty);
1018 } else if(fAveragePoly) fAveragePoly->SetFillStyle(sty);
1019}
1020
1021////////////////////////////////////////////////////////////////////////////////
1022/// Display or not the average.
1023
1025{
1026 if(disp == fDisplayAverage) return;
1027
1028 UInt_t ui=0;
1029
1030 fDisplayAverage = disp;
1031 delete fAveragePoly;
1032 fAveragePoly = nullptr;
1033 if(fAverageSlices){
1034 for(ui = 0;ui<fNcols;++ui) delete fAverageSlices[ui];
1035 }
1036 delete [] fAverageSlices;
1037 fAverageSlices = nullptr;
1038
1039 for(ui=0;ui<fNx*fNy;++ui){
1040 if (fCanvas) fCanvas->cd(ui+1);
1041 gPad->Clear();
1042 }
1043
1044 for(ui = 0; ui < fNx*fNy; ++ui){
1045 if (fCanvas) fCanvas->cd(ui+1);
1046 fPolargram->Draw("pn");
1047 fTree->LoadTree(fCurrentEntries[ui]);
1048 if(fSegmentDisplay){
1049 if(disp) DrawSlicesAverage("");
1050 DrawSlices("");
1051 } else {
1052 if(disp) DrawPolyAverage("");
1053 DrawPoly("");
1054 }
1055 AppendPad();
1056 }
1057 if (fCanvas) {
1058 fCanvas->Modified();
1059 fCanvas->Update();
1060 }
1061}
1062
1063////////////////////////////////////////////////////////////////////////////////
1064/// Set the current selected entries.
1065
1067{
1068 Int_t i;
1069 UInt_t ui=0;
1070 Int_t tnumber=-1;
1071 Long64_t entryNumber;
1072 Long64_t entry = fEntry;
1073 Int_t entriesToDisplay = fTree->GetScanField();
1074
1075 if(!fCurrentEntries) fCurrentEntries = new Long64_t[fTree->GetScanField()];
1076
1077 while(entriesToDisplay!=0){
1078 entryNumber = fTree->GetEntryNumber(entry);
1079 if(entryNumber < 0) break;
1080 Long64_t localEntry = fTree->LoadTree(entryNumber);
1081 if(localEntry < 0) break;
1082 if(tnumber != fTree->GetTreeNumber()) {
1083 tnumber = fTree->GetTreeNumber();
1084 if(fManager) fManager->UpdateFormulaLeaves();
1085 else {
1086 for(i=0;i<=fFormulas->LastIndex();++i)
1087 ((TTreeFormula*)fFormulas->At(i))->UpdateFormulaLeaves();
1088 }
1089 }
1090 Int_t ndata=1;
1091 if(fForceDim){
1092 if(fManager)
1093 ndata = fManager->GetNdata(true);
1094 else {
1095 for(ui=0;ui<fNcols;++ui){
1096 if(ndata < ((TTreeFormula*)fFormulas->At(ui))->GetNdata())
1097 ndata = ((TTreeFormula*)fFormulas->At(ui))->GetNdata();
1098 }
1099 if(fSelect && fSelect->GetNdata() == 0)
1100 ndata = 0;
1101 }
1102 }
1103 bool loaded = false;
1104 bool skip = false;
1105 // Loop over the instances of the selection condition
1106 for(Int_t inst=0;inst<ndata;++inst){
1107 if(fSelect){
1108 if(fSelect->EvalInstance(inst) == 0){
1109 skip = true;
1110 ++entry;
1111 }
1112 }
1113 if (!loaded) {
1114 // EvalInstance(0) always needs to be called so that
1115 // the proper branches are loaded.
1116 for (ui=0;ui<fNcols;ui++) {
1117 ((TTreeFormula*)fFormulas->At(ui))->EvalInstance(0);
1118 }
1119 loaded = true;
1120 } else if (inst == 0) {
1121 loaded = true;
1122 }
1123 }
1124 if(!skip){
1125 fCurrentEntries[fTree->GetScanField()-entriesToDisplay] = entryNumber;
1126 --entriesToDisplay;
1127 ++entry;
1128 }
1129 }
1130 if(fPolyList) UpdateView();
1131}
1132
1133////////////////////////////////////////////////////////////////////////////////
1134/// Set line style.
1135
1137{
1138 UInt_t ui=0;
1139
1141 for(ui=0; ui<fNx*fNy;++ui){
1142 if(fSegmentDisplay){
1143 TList *li = (TList*)fPolyList->At(ui);
1144 for(UInt_t var=0;var<fNcols;++var) ((TArc*)li->At(var))->SetLineStyle(sty);
1145 } else ((TPolyLine*)fPolyList->At(ui))->SetLineStyle(sty);
1146 }
1147}
1148
1149////////////////////////////////////////////////////////////////////////////////
1150/// Set lin color.
1151
1153{
1154 UInt_t ui=0;
1155
1157 for(ui=0; ui<fNx*fNy;++ui){
1158 if(fSegmentDisplay){
1159 TList *li = (TList*)fPolyList->At(ui);
1160 for(UInt_t var=0;var<fNcols;++var) ((TArc*)li->At(var))->SetLineColor(col);
1161 } else ((TPolyLine*)fPolyList->At(ui))->SetLineColor(col);
1162 }
1163}
1164
1165////////////////////////////////////////////////////////////////////////////////
1166///Set line width.
1167
1169{
1170 UInt_t ui=0;
1171
1173 for(ui=0; ui<fNx*fNy;++ui){
1174 if(fSegmentDisplay){
1175 TList *li = (TList*)fPolyList->At(ui);
1176 for(UInt_t var=0;var<fNcols;++var) ((TArc*)li->At(var))->SetLineWidth(wid);
1177 } else ((TPolyLine*)fPolyList->At(ui))->SetLineWidth(wid);
1178 }
1179}
1180
1181////////////////////////////////////////////////////////////////////////////////
1182/// Set fill color.
1183
1185{
1186 UInt_t ui=0;
1187
1189 for(ui=0; ui<fNx*fNy;++ui){
1190 if(fSegmentDisplay){
1191 TList *li = (TList*)fPolyList->At(ui);
1192 for(UInt_t var=0;var<fNcols;++var) ((TArc*)li->At(var))->SetFillColor(col);
1193 } else ((TPolyLine*)fPolyList->At(ui))->SetFillColor(col);
1194 }
1195}
1196
1197////////////////////////////////////////////////////////////////////////////////
1198/// Set fill style.
1199
1201{
1202 UInt_t ui=0;
1203
1205 for(ui=0; ui<fNx*fNy;++ui){
1206 if(fSegmentDisplay){
1207 TList *li = (TList*)fPolyList->At(ui);
1208 for(UInt_t var=0;var<fNcols;++var) ((TArc*)li->At(var))->SetFillStyle(sty);
1209 } else ((TPolyLine*)fPolyList->At(ui))->SetFillStyle(sty);
1210 }
1211}
1212
1213////////////////////////////////////////////////////////////////////////////////
1214/// Set number of radial divisions.
1215
1217{
1218 if(fPolargram->GetNdivRadial() == ndiv) return;
1219 fPolargram->SetNdivRadial(ndiv);
1220}
1221
1222////////////////////////////////////////////////////////////////////////////////
1223/// Set the X number of sub pads.
1224
1226{
1227 if(fNx == nx || nx <= 0) return;
1229
1230 UInt_t ui=0;
1231 Color_t lc;
1232 Style_t lt;
1233 Width_t lw;
1234 Color_t fc;
1235 Style_t fs;
1236 if(fAverageSlices){
1237 lc = fAverageSlices[0]->GetLineColor();
1238 lt = fAverageSlices[0]->GetLineStyle();
1239 lw = fAverageSlices[0]->GetLineWidth();
1240 fc = fAverageSlices[0]->GetFillColor();
1241 fs = fAverageSlices[0]->GetFillStyle();
1242 } else {
1243 lc = fAveragePoly->GetLineColor();
1244 lt = fAveragePoly->GetLineStyle();
1245 lw = fAveragePoly->GetLineWidth();
1246 fc = fAveragePoly->GetFillColor();
1247 fs = fAveragePoly->GetFillStyle();
1248 }
1249
1250 if(fSegmentDisplay){
1251 for(ui=0; ui<fNx*fNy;++ui) ((TList*)fPolyList->At(ui))->Delete();
1252 }
1253 fPolyList->Delete();
1254 delete fPolyList;
1255 fPolyList = nullptr;
1256 delete [] fCurrentEntries;
1257 fCurrentEntries = nullptr;
1258
1259 fNx = nx;
1260
1261 fTree->SetScanField(fNx*fNy);
1263 if (fCanvas) {
1264 fCanvas->Clear();
1265 fCanvas->Divide(fNx,fNy);
1266 }
1267
1268 for(ui=0; ui < fNx*fNy;++ui){
1269 if (fCanvas) fCanvas->cd(ui+1);
1270 fPolargram->Draw("pn");
1271 fTree->LoadTree(fCurrentEntries[ui]);
1272 if(fSegmentDisplay){
1274 DrawSlices("");
1275 } else {
1277 DrawPoly("");
1278 }
1279 AppendPad();
1280 }
1281
1282 if(fAverageSlices){
1283 for(ui = 0;ui<fNcols;++ui){
1284 fAverageSlices[ui]->SetLineColor(lc);
1285 fAverageSlices[ui]->SetLineStyle(lt);
1286 fAverageSlices[ui]->SetLineWidth(lw);
1287 fAverageSlices[ui]->SetFillColor(fc);
1288 fAverageSlices[ui]->SetFillStyle(fs);
1289 }
1290 } else {
1291 fAveragePoly->SetLineColor(lc);
1292 fAveragePoly->SetLineStyle(lt);
1293 fAveragePoly->SetLineWidth(lw);
1294 fAveragePoly->SetFillColor(fc);
1295 fAveragePoly->SetFillStyle(fs);
1296 }
1297}
1298
1299////////////////////////////////////////////////////////////////////////////////
1300/// Set the Y number of sub pads.
1301
1303{
1304 if(fNy == ny || ny <= 0) return;
1306
1307 UInt_t ui=0;
1308 Color_t lc;
1309 Style_t lt;
1310 Width_t lw;
1311 Color_t fc;
1312 Style_t fs;
1313 if(fAverageSlices){
1314 lc = fAverageSlices[0]->GetLineColor();
1315 lt = fAverageSlices[0]->GetLineStyle();
1316 lw = fAverageSlices[0]->GetLineWidth();
1317 fc = fAverageSlices[0]->GetFillColor();
1318 fs = fAverageSlices[0]->GetFillStyle();
1319 } else {
1320 lc = fAveragePoly->GetLineColor();
1321 lt = fAveragePoly->GetLineStyle();
1322 lw = fAveragePoly->GetLineWidth();
1323 fc = fAveragePoly->GetFillColor();
1324 fs = fAveragePoly->GetFillStyle();
1325 }
1326
1327 if(fSegmentDisplay){
1328 for(ui=0; ui<fNx*fNy;++ui) ((TList*)fPolyList->At(ui))->Delete();
1329 }
1330 fPolyList->Delete();
1331 delete fPolyList;
1332 fPolyList = nullptr;
1333 delete [] fCurrentEntries;
1334 fCurrentEntries = nullptr;
1335
1336 fNy = ny;
1337
1338 fTree->SetScanField(fNx*fNy);
1340 if (fCanvas) {
1341 fCanvas->Clear();
1342 fCanvas->Divide(fNx,fNy);
1343 }
1344
1345 for(ui=0; ui < fNx*fNy;++ui){
1346 if (fCanvas) fCanvas->cd(ui+1);
1347 fPolargram->Draw("pn");
1348 fTree->LoadTree(fCurrentEntries[ui]);
1349 if(fSegmentDisplay){
1351 DrawSlices("");
1352 } else {
1354 DrawPoly("");
1355 }
1356 AppendPad();
1357 }
1358
1359 if(fAverageSlices){
1360 for(ui = 0;ui<fNcols;++ui){
1361 fAverageSlices[ui]->SetLineColor(lc);
1362 fAverageSlices[ui]->SetLineStyle(lt);
1363 fAverageSlices[ui]->SetLineWidth(lw);
1364 fAverageSlices[ui]->SetFillColor(fc);
1365 fAverageSlices[ui]->SetFillStyle(fs);
1366 }
1367 } else {
1368 fAveragePoly->SetLineColor(lc);
1369 fAveragePoly->SetLineStyle(lt);
1370 fAveragePoly->SetLineWidth(lw);
1371 fAveragePoly->SetFillColor(fc);
1372 fAveragePoly->SetFillStyle(fs);
1373 }
1374}
1375
1376////////////////////////////////////////////////////////////////////////////////
1377/// Set the segment display or not.
1378
1380{
1381 if(seg == fSegmentDisplay) return;
1382
1383 UInt_t ui=0;
1384
1385 if(fSegmentDisplay){
1386 for(ui=0;ui<fNx*fNy;++ui){
1387 ((TList*)fPolyList->At(ui))->Delete();
1388 }
1389 }
1390 fPolyList->Delete();
1391
1392 Color_t lc = 1;
1393 Style_t lt = 1;
1394 Width_t lw = 1;
1395 Color_t fc = 1;
1396 Style_t fs = 0;
1397 if(fAverageSlices){
1398 lc = fAverageSlices[0]->GetLineColor();
1399 lt = fAverageSlices[0]->GetLineStyle();
1400 lw = fAverageSlices[0]->GetLineWidth();
1401 fc = fAverageSlices[0]->GetFillColor();
1402 fs = fAverageSlices[0]->GetFillStyle();
1403 } else {
1404 if (fAveragePoly) {
1405 lc = fAveragePoly->GetLineColor();
1406 lt = fAveragePoly->GetLineStyle();
1407 lw = fAveragePoly->GetLineWidth();
1408 fc = fAveragePoly->GetFillColor();
1409 fs = fAveragePoly->GetFillStyle();
1410 }
1411 }
1412
1413 delete fPolyList;
1414 fPolyList = nullptr;
1415 if(fAverageSlices){
1416 for(ui=0;ui<fNcols;++ui) delete fAverageSlices[ui];
1417 }
1418 delete [] fAverageSlices;
1419 fAverageSlices = nullptr;
1420 delete fAveragePoly;
1421 fAveragePoly = nullptr;
1422
1423 for(ui=0;ui<fNx*fNy;++ui){
1424 if (fCanvas) fCanvas->cd(ui+1);
1425 gPad->Clear();
1426 }
1427
1428 fSegmentDisplay = seg;
1429
1430 for(ui=0; ui < fNx*fNy;++ui){
1431 if (fCanvas) fCanvas->cd(ui+1);
1432 fPolargram->Draw("pn");
1433 fTree->LoadTree(fCurrentEntries[ui]);
1434 if(fSegmentDisplay){
1436 DrawSlices("");
1437 } else {
1439 DrawPoly("");
1440 }
1441 AppendPad();
1442 }
1443
1444 if(fAverageSlices){
1445 for(ui=0;ui<fNcols;++ui){
1446 fAverageSlices[ui]->SetLineColor(lc);
1447 fAverageSlices[ui]->SetLineStyle(lt);
1448 fAverageSlices[ui]->SetLineWidth(lw);
1449 fAverageSlices[ui]->SetFillColor(fc);
1450 fAverageSlices[ui]->SetFillStyle(fs);
1451 }
1452 } else {
1453 if (fAveragePoly) {
1454 fAveragePoly->SetLineColor(lc);
1455 fAveragePoly->SetLineStyle(lt);
1456 fAveragePoly->SetLineWidth(lw);
1457 fAveragePoly->SetFillColor(fc);
1458 fAveragePoly->SetFillStyle(fs);
1459 }
1460 }
1461 if (fCanvas) {
1462 fCanvas->Modified();
1463 fCanvas->Update();
1464 }
1465}
1466
1467////////////////////////////////////////////////////////////////////////////////
1468/// Compile selection expression if there is one.
1469
1470void TSpider::SetSelectionExpression(const char* selection)
1471{
1472 if (selection && strlen(selection)) {
1473 fSelect = new TTreeFormula("Selection",selection,fTree);
1474 // if (!fSelect) return -1;
1475 // if (!fSelect->GetNdim()) { delete fSelect; return -1; }
1476 fFormulas->Add(fSelect);
1477 }
1478}
1479
1480////////////////////////////////////////////////////////////////////////////////
1481/// Compile the variables expression from the given string varexp.
1482
1483void TSpider::SetVariablesExpression(const char* varexp)
1484{
1485 Int_t nch;
1486 fNcols=8;
1487
1488 if (!varexp) return;
1489 TObjArray *leaves = fTree->GetListOfLeaves();
1490 UInt_t nleaves = leaves->GetEntriesFast();
1491 if (nleaves < fNcols) fNcols = nleaves;
1492 nch = strlen(varexp);
1493
1494 // if varexp is empty, take first 8 columns by default
1495 Int_t allvar = 0;
1496 std::vector<TString> cnames;
1497 if (!strcmp(varexp, "*")) { fNcols = nleaves; allvar = 1; }
1498 if (nch == 0 || allvar) {
1499 UInt_t ncs = fNcols;
1500 fNcols = 0;
1501 for (UInt_t ui=0;ui<ncs;++ui) {
1502 TLeaf *lf = (TLeaf*)leaves->At(ui);
1503 if (lf->GetBranch()->GetListOfBranches()->GetEntries() > 0) continue;
1504 cnames.push_back(lf->GetName());
1505 fNcols++;
1506 }
1507 // otherwise select only the specified columns
1508 } else {
1509 fNcols = fSelector->SplitNames(varexp,cnames);
1510 }
1511
1512 // Create the TreeFormula objects corresponding to each column
1513 for (UInt_t ui=0;ui<fNcols;ui++) {
1514 fFormulas->Add(new TTreeFormula("Var1",cnames[ui].Data(),fTree));
1515 }
1516}
1517
1518////////////////////////////////////////////////////////////////////////////////
1519/// Create a TreeFormulaManager to coordinate the formulas.
1520
1522{
1523 Int_t i;
1524 if (fFormulas->LastIndex()>=0) {
1525 if (fSelect) {
1526 if (fSelect->GetManager()->GetMultiplicity() > 0 ) {
1528 for(i=0;i<=fFormulas->LastIndex();i++) {
1529 fManager->Add((TTreeFormula*)fFormulas->At(i));
1530 }
1531 fManager->Sync();
1532 }
1533 }
1534 for(i=0;i<=fFormulas->LastIndex();i++) {
1535 TTreeFormula *form = ((TTreeFormula*)fFormulas->At(i));
1536 switch( form->GetManager()->GetMultiplicity() ) {
1537 case 1:
1538 case 2:
1539 case -1:
1540 fForceDim = true;
1541 break;
1542 case 0:
1543 break;
1544 }
1545
1546 }
1547 }
1548}
1549
1550////////////////////////////////////////////////////////////////////////////////
1551/// Update the polylines or the arcs for the current entries.
1552
1554{
1555 Double_t slice = 2*TMath::Pi()/fNcols;
1556
1557 Double_t x,y,r;
1558
1559 for(UInt_t pad=1;pad <= fNx*fNy;++pad){
1560 fTree->LoadTree(fCurrentEntries[pad-1]);
1561 for(UInt_t i=0;i<fNcols;++i){
1562 r = (((TTreeFormula*)fFormulas->At(i))->EvalInstance()-fMin[i])/(fMax[i]-fMin[i]);
1563 x=r*TMath::Cos(i*slice);
1564 y=r*TMath::Sin(i*slice);
1565 if(!fSegmentDisplay) ((TPolyLine*)fPolyList->At(pad-1))->SetPoint(i,x,y);
1566 else {
1567 ((TArc*)((TList*)fPolyList->At(pad-1))->At(i))->SetR1(r);
1568 ((TArc*)((TList*)fPolyList->At(pad-1))->At(i))->SetR2(r);
1569 }
1570 }
1571 x=(((TTreeFormula*)fFormulas->At(0))->EvalInstance()-fMin[0])/(fMax[0]-fMin[0]);
1572 y=0;
1573 if(!fSegmentDisplay) ((TPolyLine*)fPolyList->At(pad-1))->SetPoint(fNcols,x,y);
1574 }
1575}
@ kHand
Definition GuiTypes.h:375
ROOT::R::TRInterface & r
Definition Object.C:4
#define e(i)
Definition RSha256.hxx:103
int Int_t
Signed integer 4 bytes (int).
Definition RtypesCore.h:59
unsigned int UInt_t
Unsigned integer 4 bytes (unsigned int).
Definition RtypesCore.h:60
short Width_t
Line width (short).
Definition RtypesCore.h:98
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
short Color_t
Color number (short).
Definition RtypesCore.h:99
long long Long64_t
Portable signed long integer 8 bytes.
Definition RtypesCore.h:83
short Style_t
Style number (short).
Definition RtypesCore.h:96
const char Option_t
Option string (const char).
Definition RtypesCore.h:80
externTEnv * gEnv
Definition TEnv.h:170
int nentries
#define gROOT
Definition TROOT.h:417
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2496
#define gPad
virtual Int_t GetNdim() const
Definition TFormula.h:238
Create an Arc.
Definition TArc.h:26
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 Float_t GetTextSize() const
Return the text size.
Definition TAttText.h:39
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition TAttText.h:48
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:53
TObjArray * GetListOfBranches()
Definition TBranch.h:255
The Canvas class.
Definition TCanvas.h:23
void Draw(Option_t *option="") override
Default Draw method for all objects.
virtual Long64_t GetN() const
Definition TEntryList.h:78
virtual void PaintLatex(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition TLeaf.h:57
TBranch * GetBranch() const
Definition TLeaf.h:119
A doubly linked list.
Definition TList.h:38
void Add(TObject *obj) override
Definition TList.h:81
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:487
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
Int_t GetEntries() const override
Return the number of objects in array (i.e.
TObject * At(Int_t idx) const override
Definition TObjArray.h:170
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:204
void Draw(Option_t *option="") override
Default Draw method for all objects.
A specialized TSelector for TTree::Draw.
void SetAverageLineStyle(Style_t sty)
Set the LineStyle of the average.
Definition TSpider.cxx:964
TPolyLine * fAveragePoly
Polygon representing the average variables value.
Definition TSpider.h:55
void SetLineColor(Color_t col) override
Set lin color.
Definition TSpider.cxx:1152
void SetAverageFillColor(Color_t col)
Set the Fill Color of the average.
Definition TSpider.cxx:1000
~TSpider() override
Destructor.
Definition TSpider.cxx:173
UInt_t fNy
Number of vertical spider plots.
Definition TSpider.h:43
void Draw(Option_t *options="") override
Draw the spider.
Definition TSpider.cxx:453
TList * fInput
Used for fSelector.
Definition TSpider.h:59
void SetNx(UInt_t nx)
Set the X number of sub pads.
Definition TSpider.cxx:1225
Double_t * fAve
[fNcols] Average value of each variable.
Definition TSpider.h:50
Style_t GetAverageLineStyle() const
Get the LineStyle of the average.
Definition TSpider.cxx:616
void GotoPrevious()
Go to the previous entries.
Definition TSpider.cxx:763
void SetFillStyle(Style_t sty) override
Set fill style.
Definition TSpider.cxx:1200
TList * fSuperposed
Superposed spider plots.
Definition TSpider.h:53
void DrawSlices(Option_t *options)
Draw the slices of the segment plot.
Definition TSpider.cxx:559
Long64_t fFirstEntry
First entry.
Definition TSpider.h:48
void SetDisplayAverage(bool disp)
Display or not the average.
Definition TSpider.cxx:1024
void UpdateView()
Update the polylines or the arcs for the current entries.
Definition TSpider.cxx:1553
TTree * fTree
Pointer to the TTree to represent.
Definition TSpider.h:54
bool fSegmentDisplay
True if displaying a segment plot.
Definition TSpider.h:68
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute the distance to the spider.
Definition TSpider.cxx:438
void SetSegmentDisplay(bool seg)
Set the segment display or not.
Definition TSpider.cxx:1379
void SetNdivRadial(Int_t div)
Set number of radial divisions.
Definition TSpider.cxx:1216
void SetCurrentEntries()
Set the current selected entries.
Definition TSpider.cxx:1066
void AddVariable(const char *varexp)
Add a variable to the plot from its expression.
Definition TSpider.cxx:214
UInt_t fNcols
Number of variables.
Definition TSpider.h:44
bool fShowRange
Show range of variables or not.
Definition TSpider.h:69
void DrawPoly(Option_t *options)
Paint the polygon representing the current entry.
Definition TSpider.cxx:529
bool fAngularLabels
True if the labels are oriented according to their axis.
Definition TSpider.h:65
void DrawSlicesAverage(Option_t *options)
Draw the slices representing the average for the segment plot.
Definition TSpider.cxx:585
void GotoEntry(Long64_t e)
Go to a specified entry.
Definition TSpider.cxx:743
TCanvas * fCanvas
! Pointer to the mother pad.
Definition TSpider.h:57
Int_t fArraySize
Actual size of the arrays.
Definition TSpider.h:45
void SetAverageLineColor(Color_t col)
Set the LineColor of the average.
Definition TSpider.cxx:976
TArc ** fAverageSlices
! Average slices.
Definition TSpider.h:56
void SetSelectionExpression(const char *selexp)
Compile selection expression if there is one.
Definition TSpider.cxx:1470
Long64_t fEntry
Present entry number in fTree.
Definition TSpider.h:46
void DeleteVariable(const char *varexp)
Delete a variable from its expression.
Definition TSpider.cxx:352
void SetAverageLineWidth(Width_t wid)
Set the LineWidth of the average.
Definition TSpider.cxx:988
TTreeFormulaManager * fManager
Coordinator for the formulas.
Definition TSpider.h:60
Width_t GetAverageLineWidth() const
Get the LineWidth of the average.
Definition TSpider.cxx:636
void SetLineWidth(Width_t wid) override
Set line width.
Definition TSpider.cxx:1168
TList * fFormulas
List of all formulas to represent.
Definition TSpider.h:58
void SetLineStyle(Style_t sty) override
Set line style.
Definition TSpider.cxx:1136
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute the corresponding event.
Definition TSpider.cxx:666
Color_t GetAverageLineColor() const
Get the LineColor of the average.
Definition TSpider.cxx:626
TTreeFormula * fSelect
Selection condition.
Definition TSpider.h:63
void InitVariables(Long64_t firstentry, Long64_t nentries)
Browse the tree to set the min, max and average value of each variable of fVar.
Definition TSpider.cxx:831
Int_t FindTextAlign(Double_t theta)
Find the alignement rule to apply for TText::SetTextAlign(Short_t).
Definition TSpider.cxx:675
void SetVariablesExpression(const char *varexp)
Compile the variables expression from the given string varexp.
Definition TSpider.cxx:1483
void SetAverageFillStyle(Style_t sty)
Set the FillStyle of the average.
Definition TSpider.cxx:1012
TSelectorDraw * fSelector
! Selector.
Definition TSpider.h:64
TSpider()
Default constructor.
Definition TSpider.cxx:75
Double_t FindTextAngle(Double_t theta)
Determine the orientation of the polar labels according to their angle.
Definition TSpider.cxx:704
Color_t GetAverageFillColor() const
Get the Fill Color of the average.
Definition TSpider.cxx:646
void SyncFormulas()
Create a TreeFormulaManager to coordinate the formulas.
Definition TSpider.cxx:1521
void AddSuperposed(TSpider *sp)
Allow to superpose several spider views.
Definition TSpider.cxx:205
void SetFillColor(Color_t col) override
Set fill color.
Definition TSpider.cxx:1184
Double_t * fMax
[fNcols] Maximum value of the variables.
Definition TSpider.h:51
Double_t * fMin
[fNcols] Minimum value of the variables.
Definition TSpider.h:52
UInt_t fNx
Number of horizontal spider plots.
Definition TSpider.h:42
void SetNy(UInt_t ny)
Set the Y number of sub pads.
Definition TSpider.cxx:1302
bool fForceDim
Force dimension.
Definition TSpider.h:67
void Paint(Option_t *options) override
Paint the spider.
Definition TSpider.cxx:917
void InitArrays(Int_t newsize)
Check if the arrays size is enough and reallocate them if not.
Definition TSpider.cxx:793
Long64_t fNentries
Number of entries.
Definition TSpider.h:47
Style_t GetAverageFillStyle() const
Get the FillStyle of the average.
Definition TSpider.cxx:656
TGraphPolargram * fPolargram
Polar graph.
Definition TSpider.h:61
void GotoNext()
Go to the next entries.
Definition TSpider.cxx:753
void DrawPolyAverage(Option_t *options)
Paint the Polygon representing the average value of the variables.
Definition TSpider.cxx:491
Long64_t * fCurrentEntries
![fNx*fNy] current selected entries;
Definition TSpider.h:49
TList * fPolyList
Polygons representing the variables.
Definition TSpider.h:62
void GotoPreceding()
Go to the last entry.
Definition TSpider.cxx:783
bool fDisplayAverage
Display or not the average.
Definition TSpider.h:66
Long64_t GetEntriesToProcess(Long64_t firstentry, Long64_t nentries) const
return the number of entries to be processed this function checks that nentries is not bigger than th...
Definition TSpider.cxx:726
void GotoFollowing()
Go to the next entry.
Definition TSpider.cxx:773
Basic string class.
Definition TString.h:138
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:641
Used to coordinate one or more TTreeFormula objects.
virtual Int_t GetMultiplicity() const
Used to pass a selection expression to the Tree drawing routine.
TTreeFormulaManager * GetManager() const
T EvalInstance(Int_t i=0, const char *stringStack[]=nullptr)
Evaluate this treeformula.
A TTree represents a columnar dataset.
Definition TTree.h:89
virtual Long64_t GetEstimate() const
Definition TTree.h:554
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
constexpr Double_t Pi()
Definition TMath.h:40
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
th1 Draw()