Logo ROOT  
Reference Guide
TPad.cxx
Go to the documentation of this file.
1// @(#)root/gpad:$Id$
2// Author: Rene Brun 12/12/94
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#include <cstring>
13#include <cstdlib>
14#include <iostream>
15
16#include "TROOT.h"
17#include "TBuffer.h"
18#include "TError.h"
19#include "TMath.h"
20#include "TSystem.h"
21#include "TStyle.h"
22#include "TH1.h"
23#include "TH2.h"
24#include "TH3.h"
25#include "TClass.h"
26#include "TBaseClass.h"
27#include "TClassTable.h"
28#include "TVirtualPS.h"
29#include "TVirtualX.h"
30#include "TVirtualViewer3D.h"
31#include "TView.h"
32#include "TPoint.h"
33#include "TGraph.h"
34#include "TMultiGraph.h"
35#include "THStack.h"
36#include "TPaveText.h"
37#include "TPaveStats.h"
38#include "TGroupButton.h"
39#include "TBrowser.h"
40#include "TVirtualGL.h"
41#include "TString.h"
42#include "TDataMember.h"
43#include "TMethod.h"
44#include "TDataType.h"
45#include "TFrame.h"
46#include "TExec.h"
47#include "TDatime.h"
48#include "TColor.h"
49#include "TCanvas.h"
50#include "TPluginManager.h"
51#include "TEnv.h"
52#include "TImage.h"
53#include "TViewer3DPad.h"
54#include "TCreatePrimitives.h"
55#include "TLegend.h"
56#include "TAtt3D.h"
57#include "TVirtualPadPainter.h"
58#include "strlcpy.h"
59#include "snprintf.h"
60
61#include "TVirtualMutex.h"
62
63static Int_t gReadLevel = 0;
64
66
68
69/** \class TPad
70\ingroup gpad
71
72The most important graphics class in the ROOT system.
73
74A Pad is contained in a Canvas.
75
76A Pad may contain other pads (unlimited pad hierarchy).
77
78A pad is a linked list of primitives of any type (graphics objects,
79histograms, detectors, tracks, etc.).
80
81Adding a new element into a pad is in general performed by the Draw
82member function of the object classes.
83
84It is important to realize that the pad is a linked list of references
85to the original object.
86For example, in case of a histogram, the histogram.Draw() operation
87only stores a reference to the histogram object and not a graphical
88representation of this histogram.
89When the mouse is used to change (say the bin content), the bin content
90of the original histogram is changed.
91
92The convention used in ROOT is that a Draw operation only adds
93a reference to the object. The effective drawing is performed
94when the canvas receives a signal to be painted.
95
96\image html gpad_pad1.png
97
98This signal is generally sent when typing carriage return in the
99command input or when a graphical operation has been performed on one
100of the pads of this canvas.
101When a Canvas/Pad is repainted, the member function Paint for all
102objects in the Pad linked list is invoked.
103
104\image html gpad_pad2.png
105
106When the mouse is moved on the Pad, The member function DistancetoPrimitive
107is called for all the elements in the pad. DistancetoPrimitive returns
108the distance in pixels to this object.
109
110When the object is within the distance window, the member function
111ExecuteEvent is called for this object.
112
113In ExecuteEvent, move, changes can be performed on the object.
114
115For examples of DistancetoPrimitive and ExecuteEvent functions,
116see classes
117~~~ {.cpp}
118 TLine::DistancetoPrimitive, TLine::ExecuteEvent
119 TBox::DistancetoPrimitive, TBox::ExecuteEvent
120 TH1::DistancetoPrimitive, TH1::ExecuteEvent
121~~~
122A Pad supports linear and log scales coordinate systems.
123The transformation coefficients are explained in TPad::ResizePad.
124*/
125
126////////////////////////////////////////////////////////////////////////////////
127/// Pad default constructor.
128
130{
131 fModified = kTRUE;
132 fTip = nullptr;
133 fPadPointer = nullptr;
134 fPrimitives = nullptr;
135 fExecs = nullptr;
136 fCanvas = nullptr;
137 fPadPaint = 0;
138 fPixmapID = -1;
139 fGLDevice = -1;
140 fCopyGLDevice = kFALSE;
141 fEmbeddedGL = kFALSE;
142 fTheta = 30;
143 fPhi = 30;
144 fNumber = 0;
145 fAbsCoord = kFALSE;
146 fEditable = kTRUE;
147 fCrosshair = 0;
148 fCrosshairPos = 0;
149 fPadView3D = nullptr;
150 fMother = (TPad*)gPad;
151
152 fAbsHNDC = 0.;
153 fAbsPixeltoXk = 0.;
154 fAbsPixeltoYk = 0.;
155 fAbsWNDC = 0.;
156 fAbsXlowNDC = 0.;
157 fAbsYlowNDC = 0.;
158 fBorderMode = 0;
159 fBorderSize = 0;
160 fPixeltoX = 0;
161 fPixeltoXk = 0.;
162 fPixeltoY = 0.;
163 fPixeltoYk = 0.;
164 fUtoAbsPixelk = 0.;
165 fUtoPixel = 0.;
166 fUtoPixelk = 0.;
167 fVtoAbsPixelk = 0.;
168 fVtoPixel = 0.;
169 fVtoPixelk = 0.;
170 fXtoAbsPixelk = 0.;
171 fXtoPixel = 0.;
172 fXtoPixelk = 0.;
173 fYtoAbsPixelk = 0.;
174 fYtoPixel = 0.;
175 fYtoPixelk = 0.;
176 fXUpNDC = 0.;
177 fYUpNDC = 0.;
178
179 fFixedAspectRatio = kFALSE;
180 fAspectRatio = 0.;
181
182 fNumPaletteColor = 0;
183 fNextPaletteColor = 0;
184 fCollideGrid = nullptr;
185 fCGnx = 0;
186 fCGny = 0;
187
188 fLogx = 0;
189 fLogy = 0;
190 fLogz = 0;
191 fGridx = 0;
192 fGridy = 0;
193 fTickx = 0;
194 fTicky = 0;
195 fFrame = nullptr;
196 fView = nullptr;
197
198 fUxmin = fUymin = fUxmax = fUymax = 0;
199
200 // Set default world coordinates to NDC [0,1]
201 fX1 = 0;
202 fX2 = 1;
203 fY1 = 0;
204 fY2 = 1;
205
206 // Set default pad range
207 fXlowNDC = 0;
208 fYlowNDC = 0;
209 fWNDC = 1;
210 fHNDC = 1;
211
212 fViewer3D = nullptr;
213 SetBit(kMustCleanup);
214
215 // the following line is temporarily disabled. It has side effects
216 // when the pad is a TDrawPanelHist or a TFitPanel.
217 // the line was supposed to fix a problem with DrawClonePad
218 // gROOT->SetSelectedPad(this);
219}
220
221////////////////////////////////////////////////////////////////////////////////
222/// Pad constructor.
223///
224/// A pad is a linked list of primitives.
225/// A pad is contained in a canvas. It may contain other pads.
226/// A pad has attributes. When a pad is created, the attributes
227/// defined in the current style are copied to the pad attributes.
228///
229/// \param[in] name pad name
230/// \param[in] title pad title
231/// \param[in] xlow [0,1] is the position of the bottom left point of the pad
232/// expressed in the mother pad reference system
233/// \param[in] ylow [0,1] is the Y position of this point.
234/// \param[in] xup [0,1] is the x position of the top right point of the pad
235/// expressed in the mother pad reference system
236/// \param[in] yup [0,1] is the Y position of this point.
237/// \param[in] color pad color
238/// \param[in] bordersize border size in pixels
239/// \param[in] bordermode border mode
240/// - bordermode = -1 box looks as it is behind the screen
241/// - bordermode = 0 no special effects
242/// - bordermode = 1 box looks as it is in front of the screen
243
244TPad::TPad(const char *name, const char *title, Double_t xlow,
245 Double_t ylow, Double_t xup, Double_t yup,
246 Color_t color, Short_t bordersize, Short_t bordermode)
247 : TVirtualPad(name,title,xlow,ylow,xup,yup,color,bordersize,bordermode)
248{
250 fTip = nullptr;
251 fBorderSize = bordersize;
252 fBorderMode = bordermode;
253 if (gPad) fCanvas = gPad->GetCanvas();
254 else fCanvas = (TCanvas*)this;
255 fMother = (TPad*)gPad;
256 fPrimitives = new TList;
257 fExecs = new TList;
258 fPadPointer = nullptr;
259 fTheta = 30;
260 fPhi = 30;
265 fFrame = nullptr;
266 fView = nullptr;
267 fPadPaint = 0;
268 fPadView3D = nullptr;
269 fPixmapID = -1; // -1 means pixmap will be created by ResizePad()
272 fNumber = 0;
275 fCrosshair = 0;
276 fCrosshairPos = 0;
277
278 fVtoAbsPixelk = 0.;
279 fVtoPixelk = 0.;
280 fVtoPixel = 0.;
281 fAbsPixeltoXk = 0.;
282 fPixeltoXk = 0.;
283 fPixeltoX = 0;
284 fAbsPixeltoYk = 0.;
285 fPixeltoYk = 0.;
286 fPixeltoY = 0.;
287 fXlowNDC = 0;
288 fYlowNDC = 0;
289 fWNDC = 1;
290 fHNDC = 1;
291 fXUpNDC = 0.;
292 fYUpNDC = 0.;
293 fAbsXlowNDC = 0.;
294 fAbsYlowNDC = 0.;
295 fAbsWNDC = 0.;
296 fAbsHNDC = 0.;
297 fXtoAbsPixelk = 0.;
298 fXtoPixelk = 0.;
299 fXtoPixel = 0.;
300 fYtoAbsPixelk = 0.;
301 fYtoPixelk = 0.;
302 fYtoPixel = 0.;
303 fUtoAbsPixelk = 0.;
304 fUtoPixelk = 0.;
305 fUtoPixel = 0.;
306
307 fUxmin = fUymin = fUxmax = fUymax = 0;
311
313 fAspectRatio = 0.;
314
317 fCollideGrid = nullptr;
318 fCGnx = 0;
319 fCGny = 0;
320
321 fViewer3D = nullptr;
322
324 // Set default world coordinates to NDC [0,1]
325 fX1 = 0;
326 fX2 = 1;
327 fY1 = 0;
328 fY2 = 1;
329
330 if (!gPad) {
331 Error("TPad", "You must create a TCanvas before creating a TPad");
332 MakeZombie();
333 return;
334 }
335
336 TPad *padsav = (TPad*)gPad;
337
338 if ((xlow < 0) || (xlow > 1) || (ylow < 0) || (ylow > 1)) {
339 Error("TPad", "illegal bottom left position: x=%f, y=%f", xlow, ylow);
340 goto zombie;
341 }
342 if ((xup < 0) || (xup > 1) || (yup < 0) || (yup > 1)) {
343 Error("TPad", "illegal top right position: x=%f, y=%f", xup, yup);
344 goto zombie;
345 }
346 if (xup-xlow <= 0) {
347 Error("TPad", "illegal width: %f", xup-xlow);
348 goto zombie;
349 }
350 if (yup-ylow <= 0) {
351 Error("TPad", "illegal height: %f", yup-ylow);
352 goto zombie;
353 }
354
358
359 fUxmin = fUymin = fUxmax = fUymax = 0;
360
361 // Set pad parameters and Compute conversion coefficients
362 SetPad(name, title, xlow, ylow, xup, yup, color, bordersize, bordermode);
363 Range(0, 0, 1, 1);
366
367 padsav->cd();
368 return;
369
370zombie:
371 // error in creating pad occurred, make this pad a zombie
372 MakeZombie();
373 padsav->cd();
374}
375
376
377////////////////////////////////////////////////////////////////////////////////
378/// Pad destructor.
379
381{
382 if (!TestBit(kNotDeleted)) return;
383 Close();
386 auto primitives = fPrimitives;
387 // In some cases, fPrimitives has the kMustCleanup bit set which will lead
388 // its destructor to call RecursiveRemove and since this pad is still
389 // likely to be (indirectly) in the list of cleanups, we must set
390 // fPrimitives to nullptr to avoid TPad::RecursiveRemove from calling
391 // a member function of a partially destructed object.
392 fPrimitives = nullptr;
393 delete primitives;
395 delete fViewer3D;
396 if (fCollideGrid) delete [] fCollideGrid;
397
398 // Required since we overload TObject::Hash.
400 if (this == gPad) gPad=nullptr;
401}
402
403////////////////////////////////////////////////////////////////////////////////
404/// Add a new TExec object to the list of Execs.
405///
406/// When an event occurs in the pad (mouse click, etc) the list of C++ commands
407/// in the list of Execs are executed via TPad::AutoExec.
408///
409/// When a pad event occurs (mouse move, click, etc) all the commands
410/// contained in the fExecs list are executed in the order found in the list.
411///
412/// This facility is activated by default. It can be deactivated by using
413/// the canvas "Option" menu.
414///
415/// The following examples of TExec commands are provided in the tutorials:
416/// macros exec1.C and exec2.C.
417///
418/// ### Example1 of use of exec1.C
419///
420/// ~~~ {.cpp}
421/// Root > TFile f("hsimple.root")
422/// Root > hpx.Draw()
423/// Root > c1.AddExec("ex1",".x exec1.C")
424/// ~~~
425///
426/// At this point you can use the mouse to click on the contour of
427/// the histogram hpx. When the mouse is clicked, the bin number and its
428/// contents are printed.
429///
430/// ### Example2 of use of exec1.C
431///
432/// ~~~ {.cpp}
433/// Root > TFile f("hsimple.root")
434/// Root > hpxpy.Draw()
435/// Root > c1.AddExec("ex2",".x exec2.C")
436/// ~~~
437///
438/// When moving the mouse in the canvas, a second canvas shows the
439/// projection along X of the bin corresponding to the Y position
440/// of the mouse. The resulting histogram is fitted with a gaussian.
441/// A "dynamic" line shows the current bin position in Y.
442/// This more elaborated example can be used as a starting point
443/// to develop more powerful interactive applications exploiting the C++
444/// interpreter as a development engine.
445
446void TPad::AddExec(const char *name, const char*command)
447{
448 if (!fExecs) fExecs = new TList;
449 TExec *ex = new TExec(name,command);
450 fExecs->Add(ex);
451}
452
453////////////////////////////////////////////////////////////////////////////////
454/// Execute the list of Execs when a pad event occurs.
455
457{
459
460 if (!fExecs) fExecs = new TList;
461 TIter next(fExecs);
462 TExec *exec;
463 while ((exec = (TExec*)next())) {
464 exec->Exec();
465 }
466}
467
468////////////////////////////////////////////////////////////////////////////////
469/// Browse pad.
470
472{
473 cd();
475}
476
477////////////////////////////////////////////////////////////////////////////////
478/// Build a legend from the graphical objects in the pad.
479///
480/// A simple method to build automatically a TLegend from the primitives in a TPad.
481///
482/// Only those deriving from TAttLine, TAttMarker and TAttFill are added, excluding
483/// TPave and TFrame derived classes.
484///
485/// \return The built TLegend
486///
487/// \param[in] x1, y1, x2, y2 The TLegend coordinates
488/// \param[in] title The legend title. By default it is " "
489/// \param[in] option The TLegend option
490///
491/// The caller program owns the returned TLegend.
492///
493/// If the pad contains some TMultiGraph or THStack the individual
494/// graphs or histograms in them are added to the TLegend.
495///
496/// ### Automatic placement of the legend
497/// If `x1` is equal to `x2` and `y1` is equal to `y2` the legend will be automatically
498/// placed to avoid overlapping with the existing primitives already displayed.
499/// `x1` is considered as the width of the legend and `y1` the height. By default
500/// the legend is automatically placed with width = `x1`= `x2` = 0.3 and
501/// height = `y1`= `y2` = 0.21.
502
504 const char* title, Option_t *option)
505{
507 if (!lop) return 0;
508 TLegend *leg=0;
509 TIter next(lop);
510 TString mes;
511 TObject *o=0;
512 TString opt("");
513 while( (o=next()) ) {
516 ( !(o->InheritsFrom(TFrame::Class())) && !(o->InheritsFrom(TPave::Class())) )) {
517 if (!leg) leg = new TLegend(x1, y1, x2, y2, title);
518 if (o->InheritsFrom(TNamed::Class()) && strlen(((TNamed *)o)->GetTitle()))
519 mes = ((TNamed *)o)->GetTitle();
520 else if (strlen(o->GetName()))
521 mes = o->GetName();
522 else
523 mes = o->ClassName();
524 if (strlen(option)) {
525 opt = option;
526 } else {
527 if (o->InheritsFrom(TAttLine::Class())) opt += "l";
528 if (o->InheritsFrom(TAttMarker::Class())) opt += "p";
529 if (o->InheritsFrom(TAttFill::Class())) opt += "f";
530 }
531 leg->AddEntry(o,mes.Data(),opt.Data());
532 } else if ( o->InheritsFrom(TMultiGraph::Class() ) ) {
533 if (!leg) leg = new TLegend(x1, y1, x2, y2, title);
534 TList * grlist = ((TMultiGraph *)o)->GetListOfGraphs();
535 TIter nextgraph(grlist);
536 TGraph * gr;
537 TObject * obj;
538 while ((obj = nextgraph())) {
539 gr = (TGraph*) obj;
540 if (strlen(gr->GetTitle())) mes = gr->GetTitle();
541 else if (strlen(gr->GetName())) mes = gr->GetName();
542 else mes = gr->ClassName();
543 if (strlen(option)) opt = option;
544 else opt = "lpf";
545 leg->AddEntry( obj, mes.Data(), opt );
546 }
547 } else if ( o->InheritsFrom(THStack::Class() ) ) {
548 if (!leg) leg = new TLegend(x1, y1, x2, y2, title);
549 TList * hlist = ((THStack *)o)->GetHists();
550 TIter nexthist(hlist);
551 TH1 * hist;
552 TObject * obj;
553 while ((obj = nexthist())) {
554 hist = (TH1*) obj;
555 if (strlen(hist->GetTitle())) mes = hist->GetTitle();
556 else if (strlen(hist->GetName())) mes = hist->GetName();
557 else mes = hist->ClassName();
558 if (strlen(option)) opt = option;
559 else opt = "lpf";
560 leg->AddEntry( obj, mes.Data(), opt );
561 }
562 }
563 opt = "";
564 }
565 if (leg) {
566 TVirtualPad *gpadsave;
567 gpadsave = gPad;
568 this->cd();
569 leg->Draw();
570 gpadsave->cd();
571 } else {
572 Info("BuildLegend(void)","No object to build a TLegend.");
573 }
574 return leg;
575}
576
577////////////////////////////////////////////////////////////////////////////////
578/// Set Current pad.
579///
580/// When a canvas/pad is divided via TPad::Divide, one can directly
581/// set the current path to one of the subdivisions.
582/// See TPad::Divide for the convention to number sub-pads.
583///
584/// Returns the new current pad, or 0 in case of failure.
585///
586/// For example:
587/// ~~~ {.cpp}
588/// c1.Divide(2,3); // create 6 pads (2 divisions along x, 3 along y).
589/// ~~~
590/// To set the current pad to the bottom right pad, do
591/// ~~~ {.cpp}
592/// c1.cd(6);
593/// ~~~
594/// Note1: c1.cd() is equivalent to c1.cd(0) and sets the current pad
595/// to c1 itself.
596///
597/// Note2: after a statement like c1.cd(6), the global variable gPad
598/// points to the current pad. One can use gPad to set attributes
599/// of the current pad.
600///
601/// Note3: One can get a pointer to one of the sub-pads of pad with:
602/// TPad *subpad = (TPad*)pad->GetPad(subpadnumber);
603
605{
606 if (!subpadnumber) {
607 gPad = this;
608 if (!gPad->IsBatch() && GetPainter()) GetPainter()->SelectDrawable(fPixmapID);
609 if (!fPrimitives) fPrimitives = new TList;
610 return gPad;
611 }
612
613 TObject *obj;
614 if (!fPrimitives) fPrimitives = new TList;
615 TIter next(fPrimitives);
616 while ((obj = next())) {
617 if (obj->InheritsFrom(TPad::Class())) {
618 Int_t n = ((TPad*)obj)->GetNumber();
619 if (n == subpadnumber) {
620 return ((TPad*)obj)->cd();
621 }
622 }
623 }
624 return 0;
625}
626
627////////////////////////////////////////////////////////////////////////////////
628/// Delete all pad primitives.
629///
630/// If the bit kClearAfterCR has been set for this pad, the Clear function
631/// will execute only after having pressed a CarriageReturn
632/// Set the bit with `mypad->SetBit(TPad::kClearAfterCR)`
633
635{
636 if (!IsEditable()) return;
637
639
640 if (!fPadPaint) {
642 if (fPrimitives) fPrimitives->Clear(option);
643 if (fFrame) {
644 if (fFrame->TestBit(kNotDeleted)) delete fFrame;
645 fFrame = nullptr;
646 }
647 }
648 if (fCanvas) fCanvas->Cleared(this);
649
650 cd();
651
652 if (TestBit(kClearAfterCR)) {
653 // Intentional do not use the return value of getchar,
654 // we just want to get it and forget it
655 getchar();
656 }
657
658 if (!gPad->IsBatch() && GetPainter()) GetPainter()->ClearDrawable();
659 if (gVirtualPS && gPad == gPad->GetCanvas()) gVirtualPS->NewPage();
660
662 fCrosshairPos = 0;
664 if (fCollideGrid) {
665 delete [] fCollideGrid;
666 fCollideGrid = nullptr;
667 fCGnx = 0;
668 fCGny = 0;
669 }
671}
672
673////////////////////////////////////////////////////////////////////////////////
674/// Clipping routine: Cohen Sutherland algorithm.
675///
676/// - If Clip ==2 the segment is outside the boundary.
677/// - If Clip ==1 the segment has one point outside the boundary.
678/// - If Clip ==0 the segment is inside the boundary.
679///
680/// \param[inout] x[],y[] Segment coordinates (2 points)
681/// \param[in] xclipl,yclipb,xclipr,yclipt Clipping boundary
682
683Int_t TPad::Clip(Float_t *x, Float_t *y, Float_t xclipl, Float_t yclipb, Float_t xclipr, Float_t yclipt)
684{
685 const Float_t kP=10000;
686 Int_t clip = 0;
687
688 for (Int_t i=0;i<2;i++) {
689 if (TMath::Abs(xclipl-x[i]) <= TMath::Abs(xclipr-xclipl)/kP) x[i] = xclipl;
690 if (TMath::Abs(xclipr-x[i]) <= TMath::Abs(xclipr-xclipl)/kP) x[i] = xclipr;
691 if (TMath::Abs(yclipb-y[i]) <= TMath::Abs(yclipt-yclipb)/kP) y[i] = yclipb;
692 if (TMath::Abs(yclipt-y[i]) <= TMath::Abs(yclipt-yclipb)/kP) y[i] = yclipt;
693 }
694
695 // Compute the first endpoint codes.
696 Int_t code1 = ClippingCode(x[0],y[0],xclipl,yclipb,xclipr,yclipt);
697 Int_t code2 = ClippingCode(x[1],y[1],xclipl,yclipb,xclipr,yclipt);
698
699 Double_t xt=0, yt=0;
700 Int_t clipped = 0; //this variable could be used in a future version
701 while(code1 + code2) {
702 clipped = 1;
703
704 // The line lies entirely outside the clipping boundary
705 if (code1&code2) {
706 clip = 2;
707 return clip;
708 }
709
710 // The line is subdivided into several parts
711 Int_t ic = code1;
712 if (ic == 0) ic = code2;
713 if (ic & 0x1) {
714 yt = y[0] + (y[1]-y[0])*(xclipl-x[0])/(x[1]-x[0]);
715 xt = xclipl;
716 }
717 if (ic & 0x2) {
718 yt = y[0] + (y[1]-y[0])*(xclipr-x[0])/(x[1]-x[0]);
719 xt = xclipr;
720 }
721 if (ic & 0x4) {
722 xt = x[0] + (x[1]-x[0])*(yclipb-y[0])/(y[1]-y[0]);
723 yt = yclipb;
724 }
725 if (ic & 0x8) {
726 xt = x[0] + (x[1]-x[0])*(yclipt-y[0])/(y[1]-y[0]);
727 yt = yclipt;
728 }
729 if (ic == code1) {
730 x[0] = xt;
731 y[0] = yt;
732 code1 = ClippingCode(xt,yt,xclipl,yclipb,xclipr,yclipt);
733 } else {
734 x[1] = xt;
735 y[1] = yt;
736 code2 = ClippingCode(xt,yt,xclipl,yclipb,xclipr,yclipt);
737 }
738 }
739 clip = clipped;
740 return clip;
741}
742
743/// @copydoc TPad::Clip(Float_t*,Float_t*,Float_t,Float_t,Float_t,Float_t)
744
746{
747 const Double_t kP=10000;
748 Int_t clip = 0;
749
750 for (Int_t i=0;i<2;i++) {
751 if (TMath::Abs(xclipl-x[i]) <= TMath::Abs(xclipr-xclipl)/kP) x[i] = xclipl;
752 if (TMath::Abs(xclipr-x[i]) <= TMath::Abs(xclipr-xclipl)/kP) x[i] = xclipr;
753 if (TMath::Abs(yclipb-y[i]) <= TMath::Abs(yclipt-yclipb)/kP) y[i] = yclipb;
754 if (TMath::Abs(yclipt-y[i]) <= TMath::Abs(yclipt-yclipb)/kP) y[i] = yclipt;
755 }
756
757 // Compute the first endpoint codes.
758 Int_t code1 = 0;
759 if (x[0] < xclipl) code1 = code1 | 0x1;
760 if (x[0] > xclipr) code1 = code1 | 0x2;
761 if (y[0] < yclipb) code1 = code1 | 0x4;
762 if (y[0] > yclipt) code1 = code1 | 0x8;
763 Int_t code2 = 0;
764 if (x[1] < xclipl) code2 = code2 | 0x1;
765 if (x[1] > xclipr) code2 = code2 | 0x2;
766 if (y[1] < yclipb) code2 = code2 | 0x4;
767 if (y[1] > yclipt) code2 = code2 | 0x8;
768
769 Double_t xt=0, yt=0;
770 Int_t clipped = 0; //this variable could be used in a future version
771 while(code1 + code2) {
772 clipped = 1;
773
774 // The line lies entirely outside the clipping boundary
775 if (code1&code2) {
776 clip = 2;
777 return clip;
778 }
779
780 // The line is subdivided into several parts
781 Int_t ic = code1;
782 if (ic == 0) ic = code2;
783 if (ic & 0x1) {
784 yt = y[0] + (y[1]-y[0])*(xclipl-x[0])/(x[1]-x[0]);
785 xt = xclipl;
786 }
787 if (ic & 0x2) {
788 yt = y[0] + (y[1]-y[0])*(xclipr-x[0])/(x[1]-x[0]);
789 xt = xclipr;
790 }
791 if (ic & 0x4) {
792 xt = x[0] + (x[1]-x[0])*(yclipb-y[0])/(y[1]-y[0]);
793 yt = yclipb;
794 }
795 if (ic & 0x8) {
796 xt = x[0] + (x[1]-x[0])*(yclipt-y[0])/(y[1]-y[0]);
797 yt = yclipt;
798 }
799 if (ic == code1) {
800 x[0] = xt;
801 y[0] = yt;
802 code1 = ClippingCode(xt,yt,xclipl,yclipb,xclipr,yclipt);
803 } else {
804 x[1] = xt;
805 y[1] = yt;
806 code2 = ClippingCode(xt,yt,xclipl,yclipb,xclipr,yclipt);
807 }
808 }
809 clip = clipped;
810 return clip;
811}
812
813////////////////////////////////////////////////////////////////////////////////
814/// Compute the endpoint codes for TPad::Clip.
815
817{
818 Int_t code = 0;
819 if (x < xcl1) code = code | 0x1;
820 if (x > xcl2) code = code | 0x2;
821 if (y < ycl1) code = code | 0x4;
822 if (y > ycl2) code = code | 0x8;
823 return code;
824}
825
826////////////////////////////////////////////////////////////////////////////////
827/// Clip polygon using the Sutherland-Hodgman algorithm.
828///
829/// \param[in] n Number of points in the polygon to
830/// be clipped
831/// \param[in] x,y Polygon x[n], y[n] do be clipped vertices
832/// \param[in] xclipl,yclipb,xclipr,yclipt Clipping boundary
833/// \param[out] nn Number of points in xc and yc
834/// \param[out] xc,yc Clipped polygon vertices. The Int_t
835/// returned by this function is
836/// the number of points in the clipped
837/// polygon. These vectors must
838/// be allocated by the calling function.
839/// A size of 2*n for each is
840/// enough.
841///
842/// Sutherland and Hodgman's polygon-clipping algorithm uses a divide-and-conquer
843/// strategy: It solves a series of simple and identical problems that, when
844/// combined, solve the overall problem. The simple problem is to clip a polygon
845/// against a single infinite clip edge. Four clip edges, each defining one boundary
846/// of the clip rectangle, successively clip a polygon against a clip rectangle.
847///
848/// Steps of Sutherland-Hodgman's polygon-clipping algorithm:
849///
850/// * Polygons can be clipped against each edge of the window one at a time.
851/// Windows/edge intersections, if any, are easy to find since the X or Y coordinates
852/// are already known.
853/// * Vertices which are kept after clipping against one window edge are saved for
854/// clipping against the remaining edges.
855/// * Note that the number of vertices usually changes and will often increases.
856///
857/// The clip boundary determines a visible and invisible region. The edges from
858/// vertex i to vertex i+1 can be one of four types:
859///
860/// * Case 1 : Wholly inside visible region - save endpoint
861/// * Case 2 : Exit visible region - save the intersection
862/// * Case 3 : Wholly outside visible region - save nothing
863/// * Case 4 : Enter visible region - save intersection and endpoint
864
866{
867 Int_t nc, nc2;
868 Double_t x1, y1, x2, y2, slope; // Segment to be clipped
869
870 Double_t *xc2 = new Double_t[nn];
871 Double_t *yc2 = new Double_t[nn];
872
873 // Clip against the left boundary
874 x1 = x[n-1]; y1 = y[n-1];
875 nc2 = 0;
876 Int_t i;
877 for (i=0; i<n; i++) {
878 x2 = x[i]; y2 = y[i];
879 if (x1 == x2) {
880 slope = 0;
881 } else {
882 slope = (y2-y1)/(x2-x1);
883 }
884 if (x1 >= xclipl) {
885 if (x2 < xclipl) {
886 xc2[nc2] = xclipl; yc2[nc2++] = slope*(xclipl-x1)+y1;
887 } else {
888 xc2[nc2] = x2; yc2[nc2++] = y2;
889 }
890 } else {
891 if (x2 >= xclipl) {
892 xc2[nc2] = xclipl; yc2[nc2++] = slope*(xclipl-x1)+y1;
893 xc2[nc2] = x2; yc2[nc2++] = y2;
894 }
895 }
896 x1 = x2; y1 = y2;
897 }
898
899 // Clip against the top boundary
900 x1 = xc2[nc2-1]; y1 = yc2[nc2-1];
901 nc = 0;
902 for (i=0; i<nc2; i++) {
903 x2 = xc2[i]; y2 = yc2[i];
904 if (y1 == y2) {
905 slope = 0;
906 } else {
907 slope = (x2-x1)/(y2-y1);
908 }
909 if (y1 <= yclipt) {
910 if (y2 > yclipt) {
911 xc[nc] = x1+(yclipt-y1)*slope; yc[nc++] = yclipt;
912 } else {
913 xc[nc] = x2; yc[nc++] = y2;
914 }
915 } else {
916 if (y2 <= yclipt) {
917 xc[nc] = x1+(yclipt-y1)*slope; yc[nc++] = yclipt;
918 xc[nc] = x2; yc[nc++] = y2;
919 }
920 }
921 x1 = x2; y1 = y2;
922 }
923
924 if (nc>0) {
925
926 // Clip against the right boundary
927 x1 = xc[nc-1]; y1 = yc[nc-1];
928 nc2 = 0;
929 for (i=0; i<nc; i++) {
930 x2 = xc[i]; y2 = yc[i];
931 if (x1 == x2) {
932 slope = 0;
933 } else {
934 slope = (y2-y1)/(x2-x1);
935 }
936 if (x1 <= xclipr) {
937 if (x2 > xclipr) {
938 xc2[nc2] = xclipr; yc2[nc2++] = slope*(xclipr-x1)+y1;
939 } else {
940 xc2[nc2] = x2; yc2[nc2++] = y2;
941 }
942 } else {
943 if (x2 <= xclipr) {
944 xc2[nc2] = xclipr; yc2[nc2++] = slope*(xclipr-x1)+y1;
945 xc2[nc2] = x2; yc2[nc2++] = y2;
946 }
947 }
948 x1 = x2; y1 = y2;
949 }
950
951 // Clip against the bottom boundary
952 x1 = xc2[nc2-1]; y1 = yc2[nc2-1];
953 nc = 0;
954 for (i=0; i<nc2; i++) {
955 x2 = xc2[i]; y2 = yc2[i];
956 if (y1 == y2) {
957 slope = 0;
958 } else {
959 slope = (x2-x1)/(y2-y1);
960 }
961 if (y1 >= yclipb) {
962 if (y2 < yclipb) {
963 xc[nc] = x1+(yclipb-y1)*slope; yc[nc++] = yclipb;
964 } else {
965 xc[nc] = x2; yc[nc++] = y2;
966 }
967 } else {
968 if (y2 >= yclipb) {
969 xc[nc] = x1+(yclipb-y1)*slope; yc[nc++] = yclipb;
970 xc[nc] = x2; yc[nc++] = y2;
971 }
972 }
973 x1 = x2; y1 = y2;
974 }
975 }
976
977 delete [] xc2;
978 delete [] yc2;
979
980 if (nc < 3) nc =0;
981 return nc;
982}
983
984////////////////////////////////////////////////////////////////////////////////
985/// Delete all primitives in pad and pad itself.
986/// Pad cannot be used anymore after this call.
987/// Emits signal "Closed()".
988
990{
991 if (!TestBit(kNotDeleted)) return;
992 if (!fMother) return;
993 if (!fMother->TestBit(kNotDeleted)) return;
994
995 if (fPrimitives)
997 if (fView) {
998 if (fView->TestBit(kNotDeleted)) delete fView;
999 fView = nullptr;
1000 }
1001 if (fFrame) {
1002 if (fFrame->TestBit(kNotDeleted)) delete fFrame;
1003 fFrame = nullptr;
1004 }
1005
1006 // emit signal
1007 if (IsA() != TCanvas::Class())
1008 Closed();
1009
1010 if (fPixmapID != -1) {
1011 if (gPad) {
1012 if (!gPad->IsBatch() && GetPainter())
1014 }
1015 fPixmapID = -1;
1016
1017 if (!gROOT->GetListOfCanvases()) return;
1018 if (fMother == this) {
1019 gROOT->GetListOfCanvases()->Remove(this);
1020 return; // in case of TCanvas
1021 }
1022
1023 // remove from the mother's list of primitives
1024 if (fMother) {
1027
1028 if (gPad == this) fMother->cd();
1029 }
1030 if (fCanvas) {
1031 if (fCanvas->GetPadSave() == this)
1033 if (fCanvas->GetSelectedPad() == this)
1035 if (fCanvas->GetClickSelectedPad() == this)
1037 }
1038 }
1039
1040 fMother = nullptr;
1041 if (gROOT->GetSelectedPad() == this) gROOT->SetSelectedPad(nullptr);
1042}
1043
1044////////////////////////////////////////////////////////////////////////////////
1045/// Copy the pixmap of the pad to the canvas.
1046
1048{
1049 int px, py;
1050 XYtoAbsPixel(fX1, fY2, px, py);
1051
1052 if (fPixmapID != -1 && GetPainter())
1053 GetPainter()->CopyDrawable(fPixmapID, px, py);
1054
1055 if (this == gPad) HighLight(gPad->GetHighLightColor());
1056}
1057
1058////////////////////////////////////////////////////////////////////////////////
1059/// Copy the sub-pixmaps of the pad to the canvas.
1060
1062{
1063 TObject *obj;
1064 if (!fPrimitives) fPrimitives = new TList;
1065 TIter next(GetListOfPrimitives());
1066 while ((obj = next())) {
1067 if (obj->InheritsFrom(TPad::Class())) {
1068 ((TPad*)obj)->CopyPixmap();
1069 ((TPad*)obj)->CopyPixmaps();
1070 }
1071 }
1072}
1073
1074////////////////////////////////////////////////////////////////////////////////
1075/// Remove TExec name from the list of Execs.
1076
1077void TPad::DeleteExec(const char *name)
1078{
1079 if (!fExecs) fExecs = new TList;
1081 if (!ex) return;
1082 fExecs->Remove(ex);
1083 delete ex;
1084}
1085
1086////////////////////////////////////////////////////////////////////////////////
1087/// Compute distance from point px,py to a box.
1088///
1089/// Compute the closest distance of approach from point px,py to the
1090/// edges of this pad.
1091/// The distance is computed in pixels units.
1092
1094{
1095 Int_t pxl, pyl, pxt, pyt;
1096 Int_t px1 = gPad->XtoAbsPixel(fX1);
1097 Int_t py1 = gPad->YtoAbsPixel(fY1);
1098 Int_t px2 = gPad->XtoAbsPixel(fX2);
1099 Int_t py2 = gPad->YtoAbsPixel(fY2);
1100 if (px1 < px2) {pxl = px1; pxt = px2;}
1101 else {pxl = px2; pxt = px1;}
1102 if (py1 < py2) {pyl = py1; pyt = py2;}
1103 else {pyl = py2; pyt = py1;}
1104
1105 // Are we inside the box?
1106 // ======================
1107 if ( (px > pxl && px < pxt) && (py > pyl && py < pyt) ) {
1108 if (GetFillStyle()) return 0; //*-* if pad is filled
1109 }
1110
1111 // Are we on the edges?
1112 // ====================
1113 Int_t dxl = TMath::Abs(px - pxl);
1114 if (py < pyl) dxl += pyl - py;
1115 if (py > pyt) dxl += py - pyt;
1116 Int_t dxt = TMath::Abs(px - pxt);
1117 if (py < pyl) dxt += pyl - py;
1118 if (py > pyt) dxt += py - pyt;
1119 Int_t dyl = TMath::Abs(py - pyl);
1120 if (px < pxl) dyl += pxl - px;
1121 if (px > pxt) dyl += px - pxt;
1122 Int_t dyt = TMath::Abs(py - pyt);
1123 if (px < pxl) dyt += pxl - px;
1124 if (px > pxt) dyt += px - pxt;
1125
1126 Int_t distance = dxl;
1127 if (dxt < distance) distance = dxt;
1128 if (dyl < distance) distance = dyl;
1129 if (dyt < distance) distance = dyt;
1130
1131 return distance - Int_t(0.5*fLineWidth);
1132}
1133
1134////////////////////////////////////////////////////////////////////////////////
1135/// Automatic pad generation by division.
1136///
1137/// - The current canvas is divided in nx by ny equal divisions (pads).
1138/// - xmargin is the space along x between pads in percent of canvas.
1139/// - ymargin is the space along y between pads in percent of canvas.
1140/// - color is the color of the new pads. If 0, color is the canvas color.
1141///
1142/// Pads are automatically named `canvasname_n` where `n` is the division number
1143/// starting from top left pad.
1144///
1145/// Example if canvasname=c1 , nx=2, ny=3:
1146///
1147/// \image html gpad_pad3.png
1148///
1149/// Once a pad is divided into sub-pads, one can set the current pad
1150/// to a subpad with a given division number as illustrated above
1151/// with TPad::cd(subpad_number).
1152///
1153/// For example, to set the current pad to c1_4, one can do:
1154/// ~~~ {.cpp}
1155/// c1->cd(4)
1156/// ~~~
1157/// __Note1:__ c1.cd() is equivalent to c1.cd(0) and sets the current pad
1158/// to c1 itself.
1159///
1160/// __Note2:__ after a statement like c1.cd(6), the global variable gPad
1161/// points to the current pad. One can use gPad to set attributes
1162/// of the current pad.
1163///
1164/// __Note3:__ in case xmargin <=0 and ymargin <= 0, there is no space
1165/// between pads. The current pad margins are recomputed to
1166/// optimize the layout.
1167
1168void TPad::Divide(Int_t nx, Int_t ny, Float_t xmargin, Float_t ymargin, Int_t color)
1169{
1170 if (!IsEditable()) return;
1171
1172
1173 if (gThreadXAR) {
1174 void *arr[7];
1175 arr[1] = this; arr[2] = (void*)&nx;arr[3] = (void*)& ny;
1176 arr[4] = (void*)&xmargin; arr[5] = (void *)& ymargin; arr[6] = (void *)&color;
1177 if ((*gThreadXAR)("PDCD", 7, arr, 0)) return;
1178 }
1179
1180 TPad *padsav = (TPad*)gPad;
1181 cd();
1182 if (nx <= 0) nx = 1;
1183 if (ny <= 0) ny = 1;
1184 Int_t ix,iy;
1185 Double_t x1,y1,x2,y2;
1186 Double_t dx,dy;
1187 TPad *pad;
1188 Int_t nchname = strlen(GetName())+6;
1189 Int_t nchtitle = strlen(GetTitle())+6;
1190 char *name = new char [nchname];
1191 char *title = new char [nchtitle];
1192 Int_t n = 0;
1193 if (color == 0) color = GetFillColor();
1194 if (xmargin > 0 && ymargin > 0) {
1195 //general case
1196 dy = 1/Double_t(ny);
1197 dx = 1/Double_t(nx);
1198 for (iy=0;iy<ny;iy++) {
1199 y2 = 1 - iy*dy - ymargin;
1200 y1 = y2 - dy + 2*ymargin;
1201 if (y1 < 0) y1 = 0;
1202 if (y1 > y2) continue;
1203 for (ix=0;ix<nx;ix++) {
1204 x1 = ix*dx + xmargin;
1205 x2 = x1 +dx -2*xmargin;
1206 if (x1 > x2) continue;
1207 n++;
1208 snprintf(name,nchname,"%s_%d",GetName(),n);
1209 pad = new TPad(name,name,x1,y1,x2,y2,color);
1210 pad->SetNumber(n);
1211 pad->Draw();
1212 }
1213 }
1214 } else {
1215 // special case when xmargin <= 0 && ymargin <= 0
1216 Double_t xl = GetLeftMargin();
1217 Double_t xr = GetRightMargin();
1219 Double_t yt = GetTopMargin();
1220 xl /= (1-xl+xr)*nx;
1221 xr /= (1-xl+xr)*nx;
1222 yb /= (1-yb+yt)*ny;
1223 yt /= (1-yb+yt)*ny;
1224 SetLeftMargin(xl);
1225 SetRightMargin(xr);
1226 SetBottomMargin(yb);
1227 SetTopMargin(yt);
1228 dx = (1-xl-xr)/nx;
1229 dy = (1-yb-yt)/ny;
1230 Int_t number = 0;
1231 for (Int_t i=0;i<nx;i++) {
1232 x1 = i*dx+xl;
1233 x2 = x1 + dx;
1234 if (i == 0) x1 = 0;
1235 if (i == nx-1) x2 = 1-xr;
1236 for (Int_t j=0;j<ny;j++) {
1237 number = j*nx + i +1;
1238 y2 = 1 -j*dy -yt;
1239 y1 = y2 - dy;
1240 if (j == 0) y2 = 1-yt;
1241 if (j == ny-1) y1 = 0;
1242 snprintf(name,nchname,"%s_%d",GetName(),number);
1243 snprintf(title,nchtitle,"%s_%d",GetTitle(),number);
1244 pad = new TPad(name,title,x1,y1,x2,y2);
1245 pad->SetNumber(number);
1246 pad->SetBorderMode(0);
1247 if (i == 0) pad->SetLeftMargin(xl*nx);
1248 else pad->SetLeftMargin(0);
1249 pad->SetRightMargin(0);
1250 pad->SetTopMargin(0);
1251 if (j == ny-1) pad->SetBottomMargin(yb*ny);
1252 else pad->SetBottomMargin(0);
1253 pad->Draw();
1254 }
1255 }
1256 }
1257 delete [] name;
1258 delete [] title;
1259 Modified();
1260 if (padsav) padsav->cd();
1261}
1262
1263////////////////////////////////////////////////////////////////////////////////
1264/// "n" is the total number of sub-pads. The number of sub-pads along the X
1265/// and Y axis are computed according to the square root of n.
1266
1267void TPad::DivideSquare(Int_t n, Float_t xmargin, Float_t ymargin, Int_t color)
1268{
1269 Int_t w = 1, h = 1;
1270 if (!fCanvas) {
1271 Error("DivideSquare", "No canvas associated with this pad.");
1272 return;
1273 }
1277 if (w*h < n) w++;
1278 } else {
1281 if (w*h < n) h++;
1282 }
1283
1284 Divide( w, h, xmargin, ymargin, color);
1285}
1286
1287////////////////////////////////////////////////////////////////////////////////
1288/// Draw Pad in Current pad (re-parent pad if necessary).
1289
1291{
1292 // if no canvas opened yet create a default canvas
1293 if (!gPad) {
1294 gROOT->MakeDefCanvas();
1295 }
1296
1297 // pad cannot be in itself and it can only be in one other pad at a time
1298 if (!fPrimitives) fPrimitives = new TList;
1299 if (gPad != this) {
1302 TPad *oldMother = fMother;
1303 fCanvas = gPad->GetCanvas();
1304 //
1305 fMother = (TPad*)gPad;
1306 if (oldMother != fMother || fPixmapID == -1) ResizePad();
1307 }
1308
1309 Paint();
1310
1311 if (gPad->IsRetained() && gPad != this && fMother)
1312 if (fMother->GetListOfPrimitives()) fMother->GetListOfPrimitives()->Add(this, option);
1313}
1314
1315////////////////////////////////////////////////////////////////////////////////
1316/// Draw class inheritance tree of the class to which obj belongs.
1317///
1318/// If a class B inherits from a class A, description of B is drawn
1319/// on the right side of description of A.
1320///
1321/// Member functions overridden by B are shown in class A with a blue line
1322/// crossing-out the corresponding member function.
1323
1324void TPad::DrawClassObject(const TObject *classobj, Option_t *option)
1325{
1326 if (!classobj) return;
1327 char dname[256];
1328 const Int_t kMAXLEVELS = 10;
1329 TClass *clevel[kMAXLEVELS], *cl, *cll;
1330 TBaseClass *base, *cinherit;
1331 TText *ptext = 0;
1332 TString opt=option;
1333 Double_t x,y,dy,y1,v1,v2,dv;
1334 Int_t nd,nf,nc,nkd,nkf,i,j;
1335 TPaveText *pt;
1336 Int_t maxlev = 4;
1337 if (opt.Contains("2")) maxlev = 2;
1338 if (opt.Contains("3")) maxlev = 3;
1339 if (opt.Contains("5")) maxlev = 5;
1340 if (opt.Contains("6")) maxlev = 6;
1341 if (opt.Contains("7")) maxlev = 7;
1342
1343 // Clear and Set Pad range
1344 Double_t xpad = 20.5;
1345 Double_t ypad = 27.5;
1346 Clear();
1347 Range(0,0,xpad,ypad);
1348
1349 // Find number of levels
1350 Int_t nlevel = 0;
1351 TClass *obj = (TClass*)classobj;
1352 clevel[nlevel] = obj;
1353 TList *lbase = obj->GetListOfBases();
1354 while(lbase) {
1355 base = (TBaseClass*)lbase->First();
1356 if (!base) break;
1357 if ( base->GetClassPointer() == 0) break;
1358 nlevel++;
1359 clevel[nlevel] = base->GetClassPointer();
1360 lbase = clevel[nlevel]->GetListOfBases();
1361 if (nlevel >= maxlev-1) break;
1362 }
1363 Int_t maxelem = 0;
1364 Int_t ncdraw = 0;
1365 Int_t ilevel, nelem;
1366 for (ilevel=nlevel;ilevel>=0;ilevel--) {
1367 cl = clevel[ilevel];
1368 nelem = cl->GetNdata() + cl->GetNmethods();
1369 if (nelem > maxelem) maxelem = nelem;
1370 nc = (nelem/50) + 1;
1371 ncdraw += nc;
1372 }
1373
1374 Double_t tsizcm = 0.40;
1375 Double_t x1 = 0.25;
1376 Double_t x2 = 0;
1377 Double_t dx = 3.5;
1378 if (ncdraw > 4) {
1379 dx = dx - 0.42*Double_t(ncdraw-5);
1380 if (dx < 1.3) dx = 1.3;
1381 tsizcm = tsizcm - 0.03*Double_t(ncdraw-5);
1382 if (tsizcm < 0.27) tsizcm = 0.27;
1383 }
1384 Double_t tsiz = 1.2*tsizcm/ypad;
1385
1386 // Now loop on levels
1387 for (ilevel=nlevel;ilevel>=0;ilevel--) {
1388 cl = clevel[ilevel];
1389 nelem = cl->GetNdata() + cl->GetNmethods();
1390 if (nelem > maxelem) maxelem = nelem;
1391 nc = (nelem/50) + 1;
1392 dy = 0.45;
1393 if (ilevel < nlevel) x1 = x2 + 0.5;
1394 x2 = x1 + nc*dx;
1395 v2 = ypad - 0.5;
1396 lbase = cl->GetListOfBases();
1397 cinherit = 0;
1398 if (lbase) cinherit = (TBaseClass*)lbase->First();
1399
1400 do {
1401 nd = cl->GetNdata();
1402 nf = cl->GetNmethods() - 2; //do not show default constructor and destructor
1403 if (cl->GetListOfMethods()->FindObject("Dictionary")) {
1404 nf -= 6; // do not count the Dictionary/ClassDef functions
1405 }
1406 nkf= nf/nc +1;
1407 nkd= nd/nc +1;
1408 if (nd == 0) nkd=0;
1409 if (nf == 0) nkf=0;
1410 y1 = v2 - 0.7;
1411 v1 = y1 - Double_t(nkf+nkd+nc-1)*dy;
1412 dv = v2 - v1;
1413
1414 // Create a new PaveText
1415 pt = new TPaveText(x1,v1,x2,v2);
1417 pt->SetFillColor(19);
1418 pt->Draw();
1419 pt->SetTextColor(4);
1420 pt->SetTextFont(61);
1421 pt->SetTextAlign(12);
1422 pt->SetTextSize(tsiz);
1423 TBox *box = pt->AddBox(0,(y1+0.01-v1)/dv,0,(v2-0.01-v1)/dv);
1424 if (box) box->SetFillColor(17);
1425 pt->AddLine(0,(y1-v1)/dv,0,(y1-v1)/dv);
1426 TText *title = pt->AddText(0.5,(0.5*(y1+v2)-v1)/dv,(char*)cl->GetName());
1427 title->SetTextAlign(22);
1428 title->SetTextSize(0.6*(v2-y1)/ypad);
1429
1430 // Draw data Members
1431 i = 0;
1432 x = 0.03;
1433 y = y1 + 0.5*dy;
1434 TDataMember *d;
1435 TIter nextd(cl->GetListOfDataMembers());
1436 while ((d = (TDataMember *) nextd())) {
1437 if (i >= nkd) { i = 1; y = y1 - 0.5*dy; x += 1/Double_t(nc); }
1438 else { i++; y -= dy; }
1439
1440 // Take in account the room the array index will occupy
1441
1442 Int_t dim = d->GetArrayDim();
1443 Int_t indx = 0;
1444 snprintf(dname,256,"%s",d->GetName());
1445 Int_t ldname = 0;
1446 while (indx < dim ){
1447 ldname = strlen(dname);
1448 snprintf(&dname[ldname],256-ldname,"[%d]",d->GetMaxIndex(indx));
1449 indx++;
1450 }
1451 pt->AddText(x,(y-v1)/dv,dname);
1452 }
1453
1454 // Draw a separator line
1455 Double_t ysep;
1456 if (nd) {
1457 ysep = y1 - Double_t(nkd)*dy;
1458 pt->AddLine(0,(ysep-v1)/dv,0,(ysep-v1)/dv);
1459 ysep -= 0.5*dy;
1460 } else ysep = y1;
1461
1462 // Draw Member Functions
1463 Int_t fcount = 0;
1464 i = 0;
1465 x = 0.03;
1466 y = ysep + 0.5*dy;
1467 TMethod *m;
1468 TIter nextm(cl->GetListOfMethods());
1469 while ((m = (TMethod *) nextm())) {
1470 if (
1471 !strcmp( m->GetName(), "Dictionary" ) ||
1472 !strcmp( m->GetName(), "Class_Version" ) ||
1473 !strcmp( m->GetName(), "DeclFileName" ) ||
1474 !strcmp( m->GetName(), "DeclFileLine" ) ||
1475 !strcmp( m->GetName(), "ImplFileName" ) ||
1476 !strcmp( m->GetName(), "ImplFileLine" )
1477 ) continue;
1478 fcount++;
1479 if (fcount > nf) break;
1480 if (i >= nkf) { i = 1; y = ysep - 0.5*dy; x += 1/Double_t(nc); }
1481 else { i++; y -= dy; }
1482
1483 ptext = pt->AddText(x,(y-v1)/dv,m->GetName());
1484 // Check if method is overloaded in a derived class
1485 // If yes, Change the color of the text to blue
1486 for (j=ilevel-1;j>=0;j--) {
1487 if (cl == clevel[ilevel]) {
1488 if (clevel[j]->GetMethodAny((char*)m->GetName())) {
1489 ptext->SetTextColor(15);
1490 break;
1491 }
1492 }
1493 }
1494 }
1495
1496 // Draw second inheritance classes for this class
1497 cll = 0;
1498 if (cinherit) {
1499 cinherit = (TBaseClass*)lbase->After(cinherit);
1500 if (cinherit) {
1501 cl = cinherit->GetClassPointer();
1502 cll = cl;
1503 v2 = v1 -0.4;
1504 dy = 0.35;
1505 }
1506 }
1507 } while (cll);
1508 }
1509 Update();
1510}
1511
1512////////////////////////////////////////////////////////////////////////////////
1513/// Function called to draw a crosshair in the canvas
1514///
1515/// Example:
1516/// ~~~ {.cpp}
1517/// Root > TFile f("hsimple.root");
1518/// Root > hpxpy.Draw();
1519/// Root > c1.SetCrosshair();
1520/// ~~~
1521/// When moving the mouse in the canvas, a crosshair is drawn
1522///
1523/// - if the canvas fCrosshair = 1 , the crosshair spans the full canvas
1524/// - if the canvas fCrosshair > 1 , the crosshair spans only the pad
1525
1527{
1528 if (gPad->GetEvent() == kMouseEnter) return;
1529
1530 TPad *cpad = (TPad*)gPad;
1531 TCanvas *canvas = cpad->GetCanvas();
1532 canvas->FeedbackMode(kTRUE);
1533
1534 //erase old position and draw a line at current position
1535 Int_t pxmin,pxmax,pymin,pymax,pxold,pyold,px,py;
1536 pxold = fCrosshairPos%10000;
1537 pyold = fCrosshairPos/10000;
1538 px = cpad->GetEventX();
1539 py = cpad->GetEventY()+1;
1540 if (canvas->GetCrosshair() > 1) { //crosshair only in the current pad
1541 pxmin = cpad->XtoAbsPixel(fX1);
1542 pxmax = cpad->XtoAbsPixel(fX2);
1543 pymin = cpad->YtoAbsPixel(fY1);
1544 pymax = cpad->YtoAbsPixel(fY2);
1545 } else { //default; crosshair spans the full canvas
1546 pxmin = 0;
1547 pxmax = canvas->GetWw();
1548 pymin = 0;
1549 pymax = cpad->GetWh();
1550 }
1551#ifndef R__HAS_COCOA
1552 // Not needed, no XOR with Cocoa.
1553 if(pxold) gVirtualX->DrawLine(pxold,pymin,pxold,pymax);
1554 if(pyold) gVirtualX->DrawLine(pxmin,pyold,pxmax,pyold);
1555#endif // R__HAS_COCOA
1556 if (cpad->GetEvent() == kButton1Down ||
1557 cpad->GetEvent() == kButton1Up ||
1558 cpad->GetEvent() == kMouseLeave) {
1559 fCrosshairPos = 0;
1560 return;
1561 }
1562 gVirtualX->DrawLine(px,pymin,px,pymax);
1563 gVirtualX->DrawLine(pxmin,py,pxmax,py);
1564 fCrosshairPos = px + 10000*py;
1565}
1566
1567////////////////////////////////////////////////////////////////////////////////
1568/// Draw an empty pad frame with X and Y axis.
1569///
1570/// \return The pointer to the histogram used to draw the frame.
1571///
1572/// \param[in] xmin X axis lower limit
1573/// \param[in] xmax X axis upper limit
1574/// \param[in] ymin Y axis lower limit
1575/// \param[in] ymax Y axis upper limit
1576/// \param[in] title Pad title.If title is of the form "stringt;stringx;stringy"
1577/// the pad title is set to stringt, the x axis title to
1578/// stringx, the y axis title to stringy.
1579///
1580/// #### Example:
1581///
1582/// Begin_Macro(source)
1583/// {
1584/// auto c = new TCanvas("c","c",200,10,500,300);
1585///
1586/// const Int_t n = 50;
1587/// auto g = new TGraph();
1588/// for (Int_t i=0;i<n;i++) g->SetPoint(i,i*0.1,100*sin(i*0.1+0.2));
1589///
1590/// auto frame = c->DrawFrame(0, -110, 2, 110);
1591/// frame->GetXaxis()->SetTitle("X axis");
1592///
1593/// g->Draw("L*");
1594/// }
1595/// End_Macro
1596
1598{
1599 if (!IsEditable()) return 0;
1600 TPad *padsav = (TPad*)gPad;
1601 if (this != padsav) {
1602 Warning("DrawFrame","Must be called for the current pad only");
1603 return padsav->DrawFrame(xmin,ymin,xmax,ymax,title);
1604 }
1605
1606 cd();
1607
1608 TH1F *hframe = (TH1F*)FindObject("hframe");
1609 if (hframe) delete hframe;
1610 Int_t nbins = 1000;
1611 //if log scale in X, use variable bin size linear with log(x)
1612 //this gives a better precision when zooming on the axis
1613 if (fLogx && xmin > 0 && xmax > xmin) {
1614 Double_t xminl = TMath::Log(xmin);
1615 Double_t xmaxl = TMath::Log(xmax);
1616 Double_t dx = (xmaxl-xminl)/nbins;
1617 Double_t *xbins = new Double_t[nbins+1];
1618 xbins[0] = xmin;
1619 for (Int_t i=1;i<=nbins;i++) {
1620 xbins[i] = TMath::Exp(xminl+i*dx);
1621 }
1622 hframe = new TH1F("hframe",title,nbins,xbins);
1623 delete [] xbins;
1624 } else {
1625 hframe = new TH1F("hframe",title,nbins,xmin,xmax);
1626 }
1627 hframe->SetBit(TH1::kNoStats);
1628 hframe->SetBit(kCanDelete);
1629 hframe->SetMinimum(ymin);
1630 hframe->SetMaximum(ymax);
1631 hframe->GetYaxis()->SetLimits(ymin,ymax);
1632 hframe->SetDirectory(0);
1633 hframe->Draw(" ");
1634 Update();
1635 if (padsav) padsav->cd();
1636 return hframe;
1637}
1638
1639////////////////////////////////////////////////////////////////////////////////
1640/// Static function to Display Color Table in a pad.
1641
1643{
1644 Int_t i, j;
1645 Int_t color;
1646 Double_t xlow, ylow, xup, yup, hs, ws;
1647 Double_t x1, y1, x2, y2;
1648 x1 = y1 = 0;
1649 x2 = y2 = 20;
1650
1651 gPad->SetFillColor(0);
1652 gPad->Clear();
1653 gPad->Range(x1,y1,x2,y2);
1654
1655 TText *text = new TText(0,0,"");
1656 text->SetTextFont(61);
1657 text->SetTextSize(0.07);
1658 text->SetTextAlign(22);
1659
1660 TBox *box = new TBox();
1661
1662 // Draw color table boxes.
1663 hs = (y2-y1)/Double_t(5);
1664 ws = (x2-x1)/Double_t(10);
1665 for (i=0;i<10;i++) {
1666 xlow = x1 + ws*(Double_t(i)+0.1);
1667 xup = x1 + ws*(Double_t(i)+0.9);
1668 for (j=0;j<5;j++) {
1669 ylow = y1 + hs*(Double_t(j)+0.1);
1670 yup = y1 + hs*(Double_t(j)+0.9);
1671 color = 10*j + i;
1672 box->SetFillStyle(1001);
1673 box->SetFillColor(color);
1674 box->DrawBox(xlow, ylow, xup, yup);
1675 box->SetFillStyle(0);
1676 box->SetLineColor(1);
1677 box->DrawBox(xlow, ylow, xup, yup);
1678 if (color == 1) text->SetTextColor(0);
1679 else text->SetTextColor(1);
1680 text->DrawText(0.5*(xlow+xup), 0.5*(ylow+yup), Form("%d",color));
1681 }
1682 }
1683}
1684
1685////////////////////////////////////////////////////////////////////////////////
1686/// Execute action corresponding to one event.
1687///
1688/// This member function is called when a TPad object is clicked.
1689///
1690/// If the mouse is clicked in one of the 4 corners of the pad (pA,pB,pC,pD)
1691/// the pad is resized with the rubber rectangle.
1692///
1693/// If the mouse is clicked inside the pad, the pad is moved.
1694///
1695/// If the mouse is clicked on the 4 edges (pL,pR,pTop,pBot), the pad is scaled
1696/// parallel to this edge.
1697///
1698/// \image html gpad_pad4.png
1699///
1700/// Note that this function duplicates on purpose the functionality
1701/// already implemented in TBox::ExecuteEvent.
1702/// If somebody modifies this function, may be similar changes should also
1703/// be applied to TBox::ExecuteEvent.
1704
1706{
1707 const Int_t kMaxDiff = 5;
1708 const Int_t kMinSize = 20;
1709 static Int_t pxorg, pyorg;
1710 static Int_t px1, px2, py1, py2, pxl, pyl, pxt, pyt, pxold, pyold;
1711 static Int_t px1p, px2p, py1p, py2p, pxlp, pylp, pxtp, pytp;
1712 static Bool_t pA, pB, pC, pD, pTop, pL, pR, pBot, pINSIDE;
1713 Int_t wx, wy;
1714 Bool_t opaque = OpaqueMoving();
1715 Bool_t ropaque = OpaqueResizing();
1716 Bool_t fixedr = HasFixedAspectRatio();
1717
1718 if (!IsEditable() && event != kMouseEnter) return;
1719 TVirtualPad *parent = GetMother();
1720 if (!parent->IsEditable()) return;
1721
1723
1724 if (fXlowNDC < 0 && event != kButton1Down) return;
1725 if (fYlowNDC < 0 && event != kButton1Down) return;
1726
1727 // keep old mouse position
1728 if (event == kButton1Down) {
1729 pxorg = px;
1730 pyorg = py;
1731 }
1732
1733 Int_t newcode = gROOT->GetEditorMode();
1734 if (newcode)
1735 pA = pB = pC = pD = pTop = pL = pR = pBot = pINSIDE = kFALSE;
1736 switch (newcode) {
1737 case kPad:
1739 break;
1740 case kMarker:
1741 case kText:
1742 TCreatePrimitives::Text(event,px,py,newcode);
1743 break;
1744 case kLine:
1746 break;
1747 case kArrow:
1749 break;
1750 case kCurlyLine:
1752 break;
1753 case kCurlyArc:
1755 break;
1756 case kPolyLine:
1758 break;
1759 case kCutG:
1761 break;
1762 case kArc:
1764 break;
1765 case kEllipse:
1767 break;
1768 case kButton:
1769 case kPave:
1770 case kPaveLabel:
1771 case kPaveText:
1772 case kPavesText:
1773 case kDiamond:
1774 TCreatePrimitives::Pave(event,px,py,newcode);
1775 return;
1776 default:
1777 break;
1778 }
1779 if (newcode) return;
1780
1781 switch (event) {
1782
1783 case kMouseEnter:
1784 if (fTip)
1786 break;
1787
1788 case kArrowKeyPress:
1789 case kButton1Down:
1790
1793
1794 GetPainter()->SetLineColor(-1);
1795 TAttLine::Modify(); //Change line attributes only if necessary
1796 if (GetFillColor())
1798 else
1801
1802 // No break !!!
1803
1804 case kMouseMotion:
1805
1806 px1 = XtoAbsPixel(fX1);
1807 py1 = YtoAbsPixel(fY1);
1808 px2 = XtoAbsPixel(fX2);
1809 py2 = YtoAbsPixel(fY2);
1810
1811 if (px1 < px2) {
1812 pxl = px1;
1813 pxt = px2;
1814 } else {
1815 pxl = px2;
1816 pxt = px1;
1817 }
1818 if (py1 < py2) {
1819 pyl = py1;
1820 pyt = py2;
1821 } else {
1822 pyl = py2;
1823 pyt = py1;
1824 }
1825
1826 px1p = parent->XtoAbsPixel(parent->GetX1()) + parent->GetBorderSize();
1827 py1p = parent->YtoAbsPixel(parent->GetY1()) - parent->GetBorderSize();
1828 px2p = parent->XtoAbsPixel(parent->GetX2()) - parent->GetBorderSize();
1829 py2p = parent->YtoAbsPixel(parent->GetY2()) + parent->GetBorderSize();
1830
1831 if (px1p < px2p) {
1832 pxlp = px1p;
1833 pxtp = px2p;
1834 } else {
1835 pxlp = px2p;
1836 pxtp = px1p;
1837 }
1838 if (py1p < py2p) {
1839 pylp = py1p;
1840 pytp = py2p;
1841 } else {
1842 pylp = py2p;
1843 pytp = py1p;
1844 }
1845
1846 pA = pB = pC = pD = pTop = pL = pR = pBot = pINSIDE = kFALSE;
1847
1848 // case pA
1849 if (TMath::Abs(px - pxl) <= kMaxDiff && TMath::Abs(py - pyl) <= kMaxDiff) {
1850 pxold = pxl; pyold = pyl; pA = kTRUE;
1852 }
1853 // case pB
1854 if (TMath::Abs(px - pxt) <= kMaxDiff && TMath::Abs(py - pyl) <= kMaxDiff) {
1855 pxold = pxt; pyold = pyl; pB = kTRUE;
1857 }
1858 // case pC
1859 if (TMath::Abs(px - pxt) <= kMaxDiff && TMath::Abs(py - pyt) <= kMaxDiff) {
1860 pxold = pxt; pyold = pyt; pC = kTRUE;
1862 }
1863 // case pD
1864 if (TMath::Abs(px - pxl) <= kMaxDiff && TMath::Abs(py - pyt) <= kMaxDiff) {
1865 pxold = pxl; pyold = pyt; pD = kTRUE;
1867 }
1868
1869 if ((px > pxl+kMaxDiff && px < pxt-kMaxDiff) &&
1870 TMath::Abs(py - pyl) < kMaxDiff) { // top edge
1871 pxold = pxl; pyold = pyl; pTop = kTRUE;
1873 }
1874
1875 if ((px > pxl+kMaxDiff && px < pxt-kMaxDiff) &&
1876 TMath::Abs(py - pyt) < kMaxDiff) { // bottom edge
1877 pxold = pxt; pyold = pyt; pBot = kTRUE;
1879 }
1880
1881 if ((py > pyl+kMaxDiff && py < pyt-kMaxDiff) &&
1882 TMath::Abs(px - pxl) < kMaxDiff) { // left edge
1883 pxold = pxl; pyold = pyl; pL = kTRUE;
1885 }
1886
1887 if ((py > pyl+kMaxDiff && py < pyt-kMaxDiff) &&
1888 TMath::Abs(px - pxt) < kMaxDiff) { // right edge
1889 pxold = pxt; pyold = pyt; pR = kTRUE;
1891 }
1892
1893 if ((px > pxl+kMaxDiff && px < pxt-kMaxDiff) &&
1894 (py > pyl+kMaxDiff && py < pyt-kMaxDiff)) { // inside box
1895 pxold = px; pyold = py; pINSIDE = kTRUE;
1896 if (event == kButton1Down)
1898 else
1900 }
1901
1902 fResizing = kFALSE;
1903 if (pA || pB || pC || pD || pTop || pL || pR || pBot)
1904 fResizing = kTRUE;
1905
1906 if (!pA && !pB && !pC && !pD && !pTop && !pL && !pR && !pBot && !pINSIDE)
1908
1909 break;
1910
1911 case kArrowKeyRelease:
1912 case kButton1Motion:
1913
1914 if (TestBit(kCannotMove)) break;
1915 wx = wy = 0;
1916
1917 if (pA) {
1918 if (!ropaque) gVirtualX->DrawBox(pxold, pyt, pxt, pyold, TVirtualX::kHollow);
1919 if (px > pxt-kMinSize) { px = pxt-kMinSize; wx = px; }
1920 if (py > pyt-kMinSize) { py = pyt-kMinSize; wy = py; }
1921 if (px < pxlp) { px = pxlp; wx = px; }
1922 if (py < pylp) { py = pylp; wy = py; }
1923 if (fixedr) {
1924 Double_t dy = Double_t(TMath::Abs(pxt-px))/parent->UtoPixel(1.) /
1926 Int_t npy2 = pyt - TMath::Abs(parent->VtoAbsPixel(dy) -
1927 parent->VtoAbsPixel(0));
1928 if (npy2 < pylp) {
1929 px = pxold;
1930 py = pyold;
1931 } else
1932 py = npy2;
1933
1934 wx = wy = 0;
1935 }
1936 if (!ropaque) gVirtualX->DrawBox(px, pyt, pxt, py, TVirtualX::kHollow);
1937 }
1938 if (pB) {
1939 if (!ropaque) gVirtualX->DrawBox(pxl , pyt, pxold, pyold, TVirtualX::kHollow);
1940 if (px < pxl+kMinSize) { px = pxl+kMinSize; wx = px; }
1941 if (py > pyt-kMinSize) { py = pyt-kMinSize; wy = py; }
1942 if (px > pxtp) { px = pxtp; wx = px; }
1943 if (py < pylp) { py = pylp; wy = py; }
1944 if (fixedr) {
1945 Double_t dy = Double_t(TMath::Abs(pxl-px))/parent->UtoPixel(1.) /
1947 Int_t npy2 = pyt - TMath::Abs(parent->VtoAbsPixel(dy) -
1948 parent->VtoAbsPixel(0));
1949 if (npy2 < pylp) {
1950 px = pxold;
1951 py = pyold;
1952 } else
1953 py = npy2;
1954
1955 wx = wy = 0;
1956 }
1957 if (!ropaque) gVirtualX->DrawBox(pxl , pyt, px , py, TVirtualX::kHollow);
1958 }
1959 if (pC) {
1960 if (!ropaque) gVirtualX->DrawBox(pxl , pyl, pxold, pyold, TVirtualX::kHollow);
1961 if (px < pxl+kMinSize) { px = pxl+kMinSize; wx = px; }
1962 if (py < pyl+kMinSize) { py = pyl+kMinSize; wy = py; }
1963 if (px > pxtp) { px = pxtp; wx = px; }
1964 if (py > pytp) { py = pytp; wy = py; }
1965 if (fixedr) {
1966 Double_t dy = Double_t(TMath::Abs(pxl-px))/parent->UtoPixel(1.) /
1968 Int_t npy2 = pyl + TMath::Abs(parent->VtoAbsPixel(dy) -
1969 parent->VtoAbsPixel(0));
1970 if (npy2 > pytp) {
1971 px = pxold;
1972 py = pyold;
1973 } else
1974 py = npy2;
1975
1976 wx = wy = 0;
1977 }
1978 if (!ropaque) gVirtualX->DrawBox(pxl, pyl, px, py, TVirtualX::kHollow);
1979 }
1980 if (pD) {
1981 if (!ropaque) gVirtualX->DrawBox(pxold, pyold, pxt, pyl, TVirtualX::kHollow);
1982 if (px > pxt-kMinSize) { px = pxt-kMinSize; wx = px; }
1983 if (py < pyl+kMinSize) { py = pyl+kMinSize; wy = py; }
1984 if (px < pxlp) { px = pxlp; wx = px; }
1985 if (py > pytp) { py = pytp; wy = py; }
1986 if (fixedr) {
1987 Double_t dy = Double_t(TMath::Abs(pxt-px))/parent->UtoPixel(1.) /
1989 Int_t npy2 = pyl + TMath::Abs(parent->VtoAbsPixel(dy) -
1990 parent->VtoAbsPixel(0));
1991 if (npy2 > pytp) {
1992 px = pxold;
1993 py = pyold;
1994 } else
1995 py = npy2;
1996
1997 wx = wy = 0;
1998 }
1999 if (!ropaque) gVirtualX->DrawBox(px, py, pxt, pyl, TVirtualX::kHollow);
2000 }
2001 if (pTop) {
2002 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2003 py2 += py - pyold;
2004 if (py2 > py1-kMinSize) { py2 = py1-kMinSize; wy = py2; }
2005 if (py2 < py2p) { py2 = py2p; wy = py2; }
2006 if (fixedr) {
2007 Double_t dx = Double_t(TMath::Abs(py2-py1))/parent->VtoPixel(0) *
2009 Int_t npx2 = px1 + parent->UtoPixel(dx);
2010 if (npx2 > px2p)
2011 py2 -= py - pyold;
2012 else
2013 px2 = npx2;
2014 }
2015 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2016 }
2017 if (pBot) {
2018 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2019 py1 += py - pyold;
2020 if (py1 < py2+kMinSize) { py1 = py2+kMinSize; wy = py1; }
2021 if (py1 > py1p) { py1 = py1p; wy = py1; }
2022 if (fixedr) {
2023 Double_t dx = Double_t(TMath::Abs(py2-py1))/parent->VtoPixel(0) *
2025 Int_t npx2 = px1 + parent->UtoPixel(dx);
2026 if (npx2 > px2p)
2027 py1 -= py - pyold;
2028 else
2029 px2 = npx2;
2030 }
2031 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2032 }
2033 if (pL) {
2034 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2035 px1 += px - pxold;
2036 if (px1 > px2-kMinSize) { px1 = px2-kMinSize; wx = px1; }
2037 if (px1 < px1p) { px1 = px1p; wx = px1; }
2038 if (fixedr) {
2039 Double_t dy = Double_t(TMath::Abs(px2-px1))/parent->UtoPixel(1.) /
2041 Int_t npy2 = py1 - TMath::Abs(parent->VtoAbsPixel(dy) -
2042 parent->VtoAbsPixel(0));
2043 if (npy2 < py2p)
2044 px1 -= px - pxold;
2045 else
2046 py2 = npy2;
2047 }
2048 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2049 }
2050 if (pR) {
2051 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2052 px2 += px - pxold;
2053 if (px2 < px1+kMinSize) { px2 = px1+kMinSize; wx = px2; }
2054 if (px2 > px2p) { px2 = px2p; wx = px2; }
2055 if (fixedr) {
2056 Double_t dy = Double_t(TMath::Abs(px2-px1))/parent->UtoPixel(1.) /
2058 Int_t npy2 = py1 - TMath::Abs(parent->VtoAbsPixel(dy) -
2059 parent->VtoAbsPixel(0));
2060 if (npy2 < py2p)
2061 px2 -= px - pxold;
2062 else
2063 py2 = npy2;
2064 }
2065 if (!ropaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow);
2066 }
2067 if (pINSIDE) {
2068 if (!opaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow); // draw the old box
2069 Int_t dx = px - pxold;
2070 Int_t dy = py - pyold;
2071 px1 += dx; py1 += dy; px2 += dx; py2 += dy;
2072 if (px1 < px1p) { dx = px1p - px1; px1 += dx; px2 += dx; wx = px+dx; }
2073 if (px2 > px2p) { dx = px2 - px2p; px1 -= dx; px2 -= dx; wx = px-dx; }
2074 if (py1 > py1p) { dy = py1 - py1p; py1 -= dy; py2 -= dy; wy = py-dy; }
2075 if (py2 < py2p) { dy = py2p - py2; py1 += dy; py2 += dy; wy = py+dy; }
2076 if (!opaque) gVirtualX->DrawBox(px1, py1, px2, py2, TVirtualX::kHollow); // draw the new box
2077 }
2078
2079 if (wx || wy) {
2080 if (wx) px = wx;
2081 if (wy) py = wy;
2082 gVirtualX->Warp(px, py);
2083 }
2084
2085 pxold = px;
2086 pyold = py;
2087
2088 Double_t x1, y1, x2, y2;
2089 x1 = x2 = y1 = y2 = 0;
2090
2091 if ((!fResizing && opaque) || (fResizing && ropaque)) {
2092 if (pA) {
2093 x1 = AbsPixeltoX(pxold);
2094 y1 = AbsPixeltoY(pyt);
2095 x2 = AbsPixeltoX(pxt);
2096 y2 = AbsPixeltoY(pyold);
2097 }
2098 if (pB) {
2099 x1 = AbsPixeltoX(pxl);
2100 y1 = AbsPixeltoY(pyt);
2101 x2 = AbsPixeltoX(pxold);
2102 y2 = AbsPixeltoY(pyold);
2103 }
2104 if (pC) {
2105 x1 = AbsPixeltoX(pxl);
2106 y1 = AbsPixeltoY(pyold);
2107 x2 = AbsPixeltoX(pxold);
2108 y2 = AbsPixeltoY(pyl);
2109 }
2110 if (pD) {
2111 x1 = AbsPixeltoX(pxold);
2112 y1 = AbsPixeltoY(pyold);
2113 x2 = AbsPixeltoX(pxt);
2114 y2 = AbsPixeltoY(pyl);
2115 }
2116 if (pTop || pBot || pL || pR || pINSIDE) {
2117 x1 = AbsPixeltoX(px1);
2118 y1 = AbsPixeltoY(py1);
2119 x2 = AbsPixeltoX(px2);
2120 y2 = AbsPixeltoY(py2);
2121 }
2122
2123 if (px != pxorg || py != pyorg) {
2124
2125 // Get parent corners pixels coordinates
2126 Int_t parentpx1 = fMother->XtoAbsPixel(parent->GetX1());
2127 Int_t parentpx2 = fMother->XtoAbsPixel(parent->GetX2());
2128 Int_t parentpy1 = fMother->YtoAbsPixel(parent->GetY1());
2129 Int_t parentpy2 = fMother->YtoAbsPixel(parent->GetY2());
2130
2131 // Get pad new corners pixels coordinates
2132 Int_t apx1 = XtoAbsPixel(x1); if (apx1 < parentpx1) {apx1 = parentpx1; }
2133 Int_t apx2 = XtoAbsPixel(x2); if (apx2 > parentpx2) {apx2 = parentpx2; }
2134 Int_t apy1 = YtoAbsPixel(y1); if (apy1 > parentpy1) {apy1 = parentpy1; }
2135 Int_t apy2 = YtoAbsPixel(y2); if (apy2 < parentpy2) {apy2 = parentpy2; }
2136
2137 // Compute new pad positions in the NDC space of parent
2138 fXlowNDC = Double_t(apx1 - parentpx1)/Double_t(parentpx2 - parentpx1);
2139 fYlowNDC = Double_t(apy1 - parentpy1)/Double_t(parentpy2 - parentpy1);
2140 fWNDC = Double_t(apx2 - apx1)/Double_t(parentpx2 - parentpx1);
2141 fHNDC = Double_t(apy2 - apy1)/Double_t(parentpy2 - parentpy1);
2142 }
2143
2144 // Reset pad parameters and recompute conversion coefficients
2145 ResizePad();
2146
2147 if (pINSIDE) gPad->ShowGuidelines(this, event);
2148 if (pTop) gPad->ShowGuidelines(this, event, 't', true);
2149 if (pBot) gPad->ShowGuidelines(this, event, 'b', true);
2150 if (pL) gPad->ShowGuidelines(this, event, 'l', true);
2151 if (pR) gPad->ShowGuidelines(this, event, 'r', true);
2152 if (pA) gPad->ShowGuidelines(this, event, '1', true);
2153 if (pB) gPad->ShowGuidelines(this, event, '2', true);
2154 if (pC) gPad->ShowGuidelines(this, event, '3', true);
2155 if (pD) gPad->ShowGuidelines(this, event, '4', true);
2156
2157 Modified(kTRUE);
2158 }
2159
2160 break;
2161
2162 case kButton1Up:
2163
2164 if (gROOT->IsEscaped()) {
2165 gROOT->SetEscape(kFALSE);
2166 break;
2167 }
2168
2169 if (opaque||ropaque) {
2170 ShowGuidelines(this, event);
2171 } else {
2172 x1 = x2 = y1 = y2 = 0;
2173
2174 if (pA) {
2175 x1 = AbsPixeltoX(pxold);
2176 y1 = AbsPixeltoY(pyt);
2177 x2 = AbsPixeltoX(pxt);
2178 y2 = AbsPixeltoY(pyold);
2179 }
2180 if (pB) {
2181 x1 = AbsPixeltoX(pxl);
2182 y1 = AbsPixeltoY(pyt);
2183 x2 = AbsPixeltoX(pxold);
2184 y2 = AbsPixeltoY(pyold);
2185 }
2186 if (pC) {
2187 x1 = AbsPixeltoX(pxl);
2188 y1 = AbsPixeltoY(pyold);
2189 x2 = AbsPixeltoX(pxold);
2190 y2 = AbsPixeltoY(pyl);
2191 }
2192 if (pD) {
2193 x1 = AbsPixeltoX(pxold);
2194 y1 = AbsPixeltoY(pyold);
2195 x2 = AbsPixeltoX(pxt);
2196 y2 = AbsPixeltoY(pyl);
2197 }
2198 if (pTop || pBot || pL || pR || pINSIDE) {
2199 x1 = AbsPixeltoX(px1);
2200 y1 = AbsPixeltoY(py1);
2201 x2 = AbsPixeltoX(px2);
2202 y2 = AbsPixeltoY(py2);
2203 }
2204
2205 if (pA || pB || pC || pD || pTop || pL || pR || pBot)
2206 Modified(kTRUE);
2207
2208 gVirtualX->SetLineColor(-1);
2209 gVirtualX->SetLineWidth(-1);
2210
2211 if (px != pxorg || py != pyorg) {
2212
2213 // Get parent corners pixels coordinates
2214 Int_t parentpx1 = fMother->XtoAbsPixel(parent->GetX1());
2215 Int_t parentpx2 = fMother->XtoAbsPixel(parent->GetX2());
2216 Int_t parentpy1 = fMother->YtoAbsPixel(parent->GetY1());
2217 Int_t parentpy2 = fMother->YtoAbsPixel(parent->GetY2());
2218
2219 // Get pad new corners pixels coordinates
2220 Int_t apx1 = XtoAbsPixel(x1); if (apx1 < parentpx1) {apx1 = parentpx1; }
2221 Int_t apx2 = XtoAbsPixel(x2); if (apx2 > parentpx2) {apx2 = parentpx2; }
2222 Int_t apy1 = YtoAbsPixel(y1); if (apy1 > parentpy1) {apy1 = parentpy1; }
2223 Int_t apy2 = YtoAbsPixel(y2); if (apy2 < parentpy2) {apy2 = parentpy2; }
2224
2225 // Compute new pad positions in the NDC space of parent
2226 fXlowNDC = Double_t(apx1 - parentpx1)/Double_t(parentpx2 - parentpx1);
2227 fYlowNDC = Double_t(apy1 - parentpy1)/Double_t(parentpy2 - parentpy1);
2228 fWNDC = Double_t(apx2 - apx1)/Double_t(parentpx2 - parentpx1);
2229 fHNDC = Double_t(apy2 - apy1)/Double_t(parentpy2 - parentpy1);
2230 }
2231
2232 // Reset pad parameters and recompute conversion coefficients
2233 ResizePad();
2234
2235
2236 // emit signal
2237 RangeChanged();
2238 }
2239
2240 break;
2241
2242 case kButton1Locate:
2243
2244 ExecuteEvent(kButton1Down, px, py);
2245
2246 while (1) {
2247 px = py = 0;
2248 event = gVirtualX->RequestLocator(1, 1, px, py);
2249
2251
2252 if (event != -1) { // button is released
2253 ExecuteEvent(kButton1Up, px, py);
2254 return;
2255 }
2256 }
2257
2258 case kButton2Down:
2259
2260 Pop();
2261 break;
2262
2263 }
2264}
2265
2266////////////////////////////////////////////////////////////////////////////////
2267/// Execute action corresponding to one event for a TAxis object
2268/// (called by TAxis::ExecuteEvent.)
2269/// This member function is called when an axis is clicked with the locator
2270///
2271/// The axis range is set between the position where the mouse is pressed
2272/// and the position where it is released.
2273///
2274/// If the mouse position is outside the current axis range when it is released
2275/// the axis is unzoomed with the corresponding proportions.
2276///
2277/// Note that the mouse does not need to be in the pad or even canvas
2278/// when it is released.
2279
2281{
2282 if (!IsEditable()) return;
2283 if (!axis) return;
2285
2286 TView *view = GetView();
2287 static Int_t axisNumber;
2288 static Double_t ratio1, ratio2;
2289 static Int_t px1old, py1old, px2old, py2old;
2290 Int_t bin1, bin2, first, last;
2291 Double_t temp, xmin,xmax;
2292 Bool_t opaque = gPad->OpaqueMoving();
2293 static TBox *zoombox;
2294 Double_t zbx1=0,zbx2=0,zby1=0,zby2=0;
2295
2296 // The CONT4 option, used to paint TH2, is a special case; it uses a 3D
2297 // drawing technique to paint a 2D plot.
2298 TString opt = axis->GetParent()->GetDrawOption();
2299 opt.ToLower();
2300 Bool_t kCont4 = kFALSE;
2301 if (strstr(opt,"cont4")) {
2302 view = 0;
2303 kCont4 = kTRUE;
2304 }
2305
2306 switch (event) {
2307
2308 case kButton1Down:
2309 axisNumber = 1;
2310 if (!strcmp(axis->GetName(),"xaxis")) {
2311 axisNumber = 1;
2312 if (!IsVertical()) axisNumber = 2;
2313 }
2314 if (!strcmp(axis->GetName(),"yaxis")) {
2315 axisNumber = 2;
2316 if (!IsVertical()) axisNumber = 1;
2317 }
2318 if (!strcmp(axis->GetName(),"zaxis")) {
2319 axisNumber = 3;
2320 }
2321 if (view) {
2322 view->GetDistancetoAxis(axisNumber, px, py, ratio1);
2323 } else {
2324 if (axisNumber == 1) {
2325 ratio1 = (AbsPixeltoX(px) - GetUxmin())/(GetUxmax() - GetUxmin());
2326 px1old = XtoAbsPixel(GetUxmin()+ratio1*(GetUxmax() - GetUxmin()));
2327 py1old = YtoAbsPixel(GetUymin());
2328 px2old = px1old;
2329 py2old = YtoAbsPixel(GetUymax());
2330 } else if (axisNumber == 2) {
2331 ratio1 = (AbsPixeltoY(py) - GetUymin())/(GetUymax() - GetUymin());
2332 py1old = YtoAbsPixel(GetUymin()+ratio1*(GetUymax() - GetUymin()));
2333 px1old = XtoAbsPixel(GetUxmin());
2334 px2old = XtoAbsPixel(GetUxmax());
2335 py2old = py1old;
2336 } else {
2337 ratio1 = (AbsPixeltoY(py) - GetUymin())/(GetUymax() - GetUymin());
2338 py1old = YtoAbsPixel(GetUymin()+ratio1*(GetUymax() - GetUymin()));
2339 px1old = XtoAbsPixel(GetUxmax());
2340 px2old = XtoAbsPixel(GetX2());
2341 py2old = py1old;
2342 }
2343 if (!opaque) {
2344 gVirtualX->DrawBox(px1old, py1old, px2old, py2old, TVirtualX::kHollow);
2345 } else {
2346 if (axisNumber == 1) {
2347 zbx1 = AbsPixeltoX(px1old);
2348 zbx2 = AbsPixeltoX(px2old);
2349 zby1 = GetUymin();
2350 zby2 = GetUymax();
2351 } else if (axisNumber == 2) {
2352 zbx1 = GetUxmin();
2353 zbx2 = GetUxmax();
2354 zby1 = AbsPixeltoY(py1old);
2355 zby2 = AbsPixeltoY(py2old);
2356 }
2357 if (GetLogx()) {
2358 zbx1 = TMath::Power(10,zbx1);
2359 zbx2 = TMath::Power(10,zbx2);
2360 }
2361 if (GetLogy()) {
2362 zby1 = TMath::Power(10,zby1);
2363 zby2 = TMath::Power(10,zby2);
2364 }
2365 zoombox = new TBox(zbx1, zby1, zbx2, zby2);
2366 Int_t ci = TColor::GetColor("#7d7dff");
2367 TColor *zoomcolor = gROOT->GetColor(ci);
2368 if (!TCanvas::SupportAlpha() || !zoomcolor) zoombox->SetFillStyle(3002);
2369 else zoomcolor->SetAlpha(0.5);
2370 zoombox->SetFillColor(ci);
2371 zoombox->Draw();
2372 gPad->Modified();
2373 gPad->Update();
2374 }
2375 }
2376 if (!opaque) gVirtualX->SetLineColor(-1);
2377 // No break !!!
2378
2379 case kButton1Motion:
2380 if (view) {
2381 view->GetDistancetoAxis(axisNumber, px, py, ratio2);
2382 } else {
2383 if (!opaque) gVirtualX->DrawBox(px1old, py1old, px2old, py2old, TVirtualX::kHollow);
2384 if (axisNumber == 1) {
2385 ratio2 = (AbsPixeltoX(px) - GetUxmin())/(GetUxmax() - GetUxmin());
2386 px2old = XtoAbsPixel(GetUxmin()+ratio2*(GetUxmax() - GetUxmin()));
2387 } else {
2388 ratio2 = (AbsPixeltoY(py) - GetUymin())/(GetUymax() - GetUymin());
2389 py2old = YtoAbsPixel(GetUymin()+ratio2*(GetUymax() - GetUymin()));
2390 }
2391 if (!opaque) {
2392 gVirtualX->DrawBox(px1old, py1old, px2old, py2old, TVirtualX::kHollow);
2393 } else {
2394 if (axisNumber == 1) {
2395 zbx1 = AbsPixeltoX(px1old);
2396 zbx2 = AbsPixeltoX(px2old);
2397 zby1 = GetUymin();
2398 zby2 = GetUymax();
2399 } else if (axisNumber == 2) {
2400 zbx1 = GetUxmin();
2401 zbx2 = GetUxmax();
2402 zby1 = AbsPixeltoY(py1old);
2403 zby2 = AbsPixeltoY(py2old);
2404 }
2405 if (GetLogx()) {
2406 zbx1 = TMath::Power(10,zbx1);
2407 zbx2 = TMath::Power(10,zbx2);
2408 }
2409 if (GetLogy()) {
2410 zby1 = TMath::Power(10,zby1);
2411 zby2 = TMath::Power(10,zby2);
2412 }
2413 if (zoombox) {
2414 zoombox->SetX1(zbx1);
2415 zoombox->SetY1(zby1);
2416 zoombox->SetX2(zbx2);
2417 zoombox->SetY2(zby2);
2418 }
2419 gPad->Modified();
2420 gPad->Update();
2421 }
2422 }
2423 break;
2424
2425 case kWheelUp:
2426 bin1 = axis->GetFirst()+1;
2427 bin2 = axis->GetLast()-1;
2428 bin1 = TMath::Max(bin1, 1);
2429 bin2 = TMath::Min(bin2, axis->GetNbins());
2430 if (bin2>bin1) {
2431 axis->SetRange(bin1,bin2);
2432 gPad->Modified();
2433 gPad->Update();
2434 }
2435 break;
2436
2437 case kWheelDown:
2438 bin1 = axis->GetFirst()-1;
2439 bin2 = axis->GetLast()+1;
2440 bin1 = TMath::Max(bin1, 1);
2441 bin2 = TMath::Min(bin2, axis->GetNbins());
2442 if (bin2>bin1) {
2443 axis->SetRange(bin1,bin2);
2444 gPad->Modified();
2445 gPad->Update();
2446 }
2447 break;
2448
2449 case kButton1Up:
2450 if (gROOT->IsEscaped()) {
2451 gROOT->SetEscape(kFALSE);
2452 if (opaque && zoombox) {
2453 zoombox->Delete();
2454 zoombox = 0;
2455 }
2456 break;
2457 }
2458
2459 if (view) {
2460 view->GetDistancetoAxis(axisNumber, px, py, ratio2);
2461 if (ratio1 > ratio2) {
2462 temp = ratio1;
2463 ratio1 = ratio2;
2464 ratio2 = temp;
2465 }
2466 if (ratio2 - ratio1 > 0.05) {
2467 TH1 *hobj = (TH1*)axis->GetParent();
2468 if (axisNumber == 3 && hobj && hobj->GetDimension() != 3) {
2469 Float_t zmin = hobj->GetMinimum();
2470 Float_t zmax = hobj->GetMaximum();
2471 if(GetLogz()){
2472 if (zmin <= 0 && zmax > 0) zmin = TMath::Min((Double_t)1,
2473 (Double_t)0.001*zmax);
2474 zmin = TMath::Log10(zmin);
2475 zmax = TMath::Log10(zmax);
2476 }
2477 Float_t newmin = zmin + (zmax-zmin)*ratio1;
2478 Float_t newmax = zmin + (zmax-zmin)*ratio2;
2479 if (newmin < zmin) newmin = hobj->GetBinContent(hobj->GetMinimumBin());
2480 if (newmax > zmax) newmax = hobj->GetBinContent(hobj->GetMaximumBin());
2481 if (GetLogz()){
2482 newmin = TMath::Exp(2.302585092994*newmin);
2483 newmax = TMath::Exp(2.302585092994*newmax);
2484 }
2485 hobj->SetMinimum(newmin);
2486 hobj->SetMaximum(newmax);
2487 hobj->SetBit(TH1::kIsZoomed);
2488 } else {
2489 first = axis->GetFirst();
2490 last = axis->GetLast();
2491 bin1 = first + Int_t((last-first+1)*ratio1);
2492 bin2 = first + Int_t((last-first+1)*ratio2);
2493 bin1 = TMath::Max(bin1, 1);
2494 bin2 = TMath::Min(bin2, axis->GetNbins());
2495 axis->SetRange(bin1, bin2);
2496 }
2497 delete view;
2498 SetView(0);
2499 Modified(kTRUE);
2500 }
2501 } else {
2502 if (axisNumber == 1) {
2503 ratio2 = (AbsPixeltoX(px) - GetUxmin())/(GetUxmax() - GetUxmin());
2504 xmin = GetUxmin() +ratio1*(GetUxmax() - GetUxmin());
2505 xmax = GetUxmin() +ratio2*(GetUxmax() - GetUxmin());
2506 if (GetLogx() && !kCont4) {
2507 xmin = PadtoX(xmin);
2508 xmax = PadtoX(xmax);
2509 }
2510 } else if (axisNumber == 2) {
2511 ratio2 = (AbsPixeltoY(py) - GetUymin())/(GetUymax() - GetUymin());
2512 xmin = GetUymin() +ratio1*(GetUymax() - GetUymin());
2513 xmax = GetUymin() +ratio2*(GetUymax() - GetUymin());
2514 if (GetLogy() && !kCont4) {
2515 xmin = PadtoY(xmin);
2516 xmax = PadtoY(xmax);
2517 }
2518 } else {
2519 ratio2 = (AbsPixeltoY(py) - GetUymin())/(GetUymax() - GetUymin());
2520 xmin = ratio1;
2521 xmax = ratio2;
2522 }
2523 if (xmin > xmax) {
2524 temp = xmin;
2525 xmin = xmax;
2526 xmax = temp;
2527 temp = ratio1;
2528 ratio1 = ratio2;
2529 ratio2 = temp;
2530 }
2531
2532 // xmin and xmax need to be adjusted in case of CONT4.
2533 if (kCont4) {
2534 Double_t low = axis->GetBinLowEdge(axis->GetFirst());
2535 Double_t up = axis->GetBinUpEdge(axis->GetLast());
2536 Double_t xmi = GetUxmin();
2537 Double_t xma = GetUxmax();
2538 xmin = ((xmin-xmi)/(xma-xmi))*(up-low)+low;
2539 xmax = ((xmax-xmi)/(xma-xmi))*(up-low)+low;
2540 }
2541
2542 if (!strcmp(axis->GetName(),"xaxis")) axisNumber = 1;
2543 if (!strcmp(axis->GetName(),"yaxis")) axisNumber = 2;
2544 if (ratio2 - ratio1 > 0.05) {
2545 //update object owning this axis
2546 TH1 *hobj1 = (TH1*)axis->GetParent();
2547 bin1 = axis->FindFixBin(xmin);
2548 bin2 = axis->FindFixBin(xmax);
2549 bin1 = TMath::Max(bin1, 1);
2550 bin2 = TMath::Min(bin2, axis->GetNbins());
2551 if (axisNumber == 1) axis->SetRange(bin1,bin2);
2552 if (axisNumber == 2 && hobj1) {
2553 if (hobj1->GetDimension() == 1) {
2554 if (hobj1->GetNormFactor() != 0) {
2555 Double_t norm = hobj1->GetSumOfWeights()/hobj1->GetNormFactor();
2556 xmin *= norm;
2557 xmax *= norm;
2558 }
2559 hobj1->SetMinimum(xmin);
2560 hobj1->SetMaximum(xmax);
2561 hobj1->SetBit(TH1::kIsZoomed);
2562 } else {
2563 axis->SetRange(bin1,bin2);
2564 }
2565 }
2566 //update all histograms in the pad
2567 TIter next(GetListOfPrimitives());
2568 TObject *obj;
2569 while ((obj= next())) {
2570 if (!obj->InheritsFrom(TH1::Class())) continue;
2571 TH1 *hobj = (TH1*)obj;
2572 if (hobj == hobj1) continue;
2573 bin1 = hobj->GetXaxis()->FindFixBin(xmin);
2574 bin2 = hobj->GetXaxis()->FindFixBin(xmax);
2575 if (axisNumber == 1) {
2576 hobj->GetXaxis()->SetRange(bin1,bin2);
2577 } else if (axisNumber == 2) {
2578 if (hobj->GetDimension() == 1) {
2579 Double_t xxmin = xmin;
2580 Double_t xxmax = xmax;
2581 if (hobj->GetNormFactor() != 0) {
2582 Double_t norm = hobj->GetSumOfWeights()/hobj->GetNormFactor();
2583 xxmin *= norm;
2584 xxmax *= norm;
2585 }
2586 hobj->SetMinimum(xxmin);
2587 hobj->SetMaximum(xxmax);
2588 hobj->SetBit(TH1::kIsZoomed);
2589 } else {
2590 bin1 = hobj->GetYaxis()->FindFixBin(xmin);
2591 bin2 = hobj->GetYaxis()->FindFixBin(xmax);
2592 hobj->GetYaxis()->SetRange(bin1,bin2);
2593 }
2594 }
2595 }
2596 Modified(kTRUE);
2597 }
2598 }
2599 if (!opaque) {
2600 gVirtualX->SetLineColor(-1);
2601 } else {
2602 if (zoombox) {
2603 zoombox->Delete();
2604 zoombox = 0;
2605 }
2606 }
2607 break;
2608 }
2609}
2610
2611////////////////////////////////////////////////////////////////////////////////
2612/// Search if object named name is inside this pad or in pads inside this pad.
2613///
2614/// In case name is in several sub-pads the first one is returned.
2615
2616TObject *TPad::FindObject(const char *name) const
2617{
2618 if (!fPrimitives) return nullptr;
2620 if (found) return found;
2621 TObject *cur;
2622 TIter next(GetListOfPrimitives());
2623 while ((cur = next())) {
2624 if (cur->InheritsFrom(TPad::Class())) {
2625 found = ((TPad*)cur)->FindObject(name);
2626 if (found) return found;
2627 }
2628 }
2629 return nullptr;
2630}
2631
2632////////////////////////////////////////////////////////////////////////////////
2633/// Search if obj is in pad or in pads inside this pad.
2634///
2635/// In case obj is in several sub-pads the first one is returned.
2636
2638{
2639 if (!fPrimitives) return nullptr;
2640 TObject *found = fPrimitives->FindObject(obj);
2641 if (found) return found;
2642 TObject *cur;
2643 TIter next(GetListOfPrimitives());
2644 while ((cur = next())) {
2645 if (cur->InheritsFrom(TPad::Class())) {
2646 found = ((TPad*)cur)->FindObject(obj);
2647 if (found) return found;
2648 }
2649 }
2650 return nullptr;
2651}
2652
2653////////////////////////////////////////////////////////////////////////////////
2654/// Get canvas identifier.
2655
2657{
2658 return fCanvas ? fCanvas->GetCanvasID() : -1;
2659}
2660
2661////////////////////////////////////////////////////////////////////////////////
2662/// Get canvas implementation pointer if any
2663
2665{
2666 return fCanvas ? fCanvas->GetCanvasImp() : nullptr;
2667}
2668
2669////////////////////////////////////////////////////////////////////////////////
2670/// Get Event.
2671
2673{
2674 return fCanvas ? fCanvas->GetEvent() : 0;
2675}
2676
2677////////////////////////////////////////////////////////////////////////////////
2678/// Get X event.
2679
2681{
2682 return fCanvas ? fCanvas->GetEventX() : 0;
2683}
2684
2685////////////////////////////////////////////////////////////////////////////////
2686/// Get Y event.
2687
2689{
2690 return fCanvas ? fCanvas->GetEventY() : 0;
2691}
2692
2693////////////////////////////////////////////////////////////////////////////////
2694/// Get virtual canvas.
2695
2697{
2698 return fCanvas ? (TVirtualPad*) fCanvas : nullptr;
2699}
2700
2701////////////////////////////////////////////////////////////////////////////////
2702/// Get highlight color.
2703
2705{
2706 return fCanvas ? fCanvas->GetHighLightColor() : 0;
2707}
2708
2709////////////////////////////////////////////////////////////////////////////////
2710/// Static function (see also TPad::SetMaxPickDistance)
2711
2713{
2714 return fgMaxPickDistance;
2715}
2716
2717////////////////////////////////////////////////////////////////////////////////
2718/// Get selected.
2719
2721{
2722 if (fCanvas == this) return nullptr;
2723 return fCanvas ? fCanvas->GetSelected() : nullptr;
2724}
2725
2726////////////////////////////////////////////////////////////////////////////////
2727/// Get selected pad.
2728
2730{
2731 if (fCanvas == this) return nullptr;
2732 return fCanvas ? fCanvas->GetSelectedPad() : nullptr;
2733}
2734
2735////////////////////////////////////////////////////////////////////////////////
2736/// Get save pad.
2737
2739{
2740 if (fCanvas == this) return nullptr;
2741 return fCanvas ? fCanvas->GetPadSave() : nullptr;
2742}
2743
2744////////////////////////////////////////////////////////////////////////////////
2745/// Get Wh.
2746
2748{
2749 return fCanvas ? fCanvas->GetWh() : 0;
2750}
2751
2752////////////////////////////////////////////////////////////////////////////////
2753/// Get Ww.
2754
2756{
2757 return fCanvas ? fCanvas->GetWw() : 0;
2758}
2759
2760////////////////////////////////////////////////////////////////////////////////
2761/// Hide tool tip depending on the event type. Typically tool tips
2762/// are hidden when event is not a kMouseEnter and not a kMouseMotion
2763/// event.
2764
2766{
2767 if (event != kMouseEnter && event != kMouseMotion && fTip)
2768 gPad->CloseToolTip(fTip);
2769}
2770
2771////////////////////////////////////////////////////////////////////////////////
2772/// Is pad in batch mode ?
2773
2775{
2776 return fCanvas ? fCanvas->IsBatch() : kFALSE;
2777}
2778
2779////////////////////////////////////////////////////////////////////////////////
2780/// Is pad retained ?
2781
2783{
2784 return fCanvas ? fCanvas->IsRetained() : kFALSE;
2785}
2786
2787////////////////////////////////////////////////////////////////////////////////
2788/// Is pad moving in opaque mode ?
2789
2791{
2792 return fCanvas ? fCanvas->OpaqueMoving() : kFALSE;
2793}
2794
2795////////////////////////////////////////////////////////////////////////////////
2796/// Is pad resizing in opaque mode ?
2797
2799{
2800 return fCanvas ? fCanvas->OpaqueResizing() : kFALSE;
2801}
2802
2803////////////////////////////////////////////////////////////////////////////////
2804/// Set pad in batch mode.
2805
2807{
2808 if (fCanvas) fCanvas->SetBatch(batch);
2809}
2810
2811////////////////////////////////////////////////////////////////////////////////
2812/// Set canvas size.
2813
2815{
2816 if (fCanvas) fCanvas->SetCanvasSize(ww,wh);
2817}
2818
2819////////////////////////////////////////////////////////////////////////////////
2820/// Set cursor type.
2821
2823{
2824 if (fCanvas) fCanvas->SetCursor(cursor);
2825}
2826
2827////////////////////////////////////////////////////////////////////////////////
2828/// Set double buffer mode ON or OFF.
2829
2831{
2832 if (fCanvas) fCanvas->SetDoubleBuffer(mode);
2833}
2834
2835////////////////////////////////////////////////////////////////////////////////
2836/// Set selected.
2837
2839{
2840 if (fCanvas) fCanvas->SetSelected(obj);
2841}
2842
2843////////////////////////////////////////////////////////////////////////////////
2844/// Update pad.
2845
2847{
2848 if (fCanvas) fCanvas->Update();
2849}
2850
2851////////////////////////////////////////////////////////////////////////////////
2852/// Get frame.
2853
2855{
2856 if (!fPrimitives) fPrimitives = new TList;
2858 if (!frame) frame = (TFrame*)GetListOfPrimitives()->FindObject("TFrame");
2859 fFrame = frame;
2860 if (!fFrame) {
2861 if (!frame) fFrame = new TFrame(0,0,1,1);
2862 Int_t framecolor = GetFrameFillColor();
2863 if (!framecolor) framecolor = GetFillColor();
2864 fFrame->SetFillColor(framecolor);
2871 }
2872 return fFrame;
2873}
2874
2875////////////////////////////////////////////////////////////////////////////////
2876/// Get primitive.
2877
2879{
2880 if (!fPrimitives) return nullptr;
2881 TIter next(fPrimitives);
2882 TObject *found, *obj;
2883 while ((obj=next())) {
2884 if (!strcmp(name, obj->GetName())) return obj;
2885 if (obj->InheritsFrom(TPad::Class())) continue;
2886 found = obj->FindObject(name);
2887 if (found) return found;
2888 }
2889 return nullptr;
2890}
2891
2892////////////////////////////////////////////////////////////////////////////////
2893/// Get a pointer to subpadnumber of this pad.
2894
2895TVirtualPad *TPad::GetPad(Int_t subpadnumber) const
2896{
2897 if (!subpadnumber) {
2898 return (TVirtualPad*)this;
2899 }
2900
2901 TObject *obj;
2902 if (!fPrimitives) return nullptr;
2903 TIter next(GetListOfPrimitives());
2904 while ((obj = next())) {
2905 if (obj->InheritsFrom(TVirtualPad::Class())) {
2906 TVirtualPad *pad = (TVirtualPad*)obj;
2907 if (pad->GetNumber() == subpadnumber) return pad;
2908 }
2909 }
2910 return nullptr;
2911}
2912
2913////////////////////////////////////////////////////////////////////////////////
2914/// Return lower and upper bounds of the pad in NDC coordinates.
2915
2916void TPad::GetPadPar(Double_t &xlow, Double_t &ylow, Double_t &xup, Double_t &yup)
2917{
2918 xlow = fXlowNDC;
2919 ylow = fYlowNDC;
2920 xup = fXlowNDC+fWNDC;
2921 yup = fYlowNDC+fHNDC;
2922}
2923
2924////////////////////////////////////////////////////////////////////////////////
2925/// Return pad world coordinates range.
2926
2928{
2929 x1 = fX1;
2930 y1 = fY1;
2931 x2 = fX2;
2932 y2 = fY2;
2933}
2934
2935////////////////////////////////////////////////////////////////////////////////
2936/// Return pad axis coordinates range.
2937
2939{
2940 xmin = fUxmin;
2941 ymin = fUymin;
2942 xmax = fUxmax;
2943 ymax = fUymax;
2944}
2945
2946////////////////////////////////////////////////////////////////////////////////
2947/// Highlight pad.
2948/// do not highlight when printing on Postscript
2949
2951{
2952 if (gVirtualPS && gVirtualPS->TestBit(kPrintingPS)) return;
2953
2954 if (color <= 0) return;
2955
2957
2958 // We do not want to have active(executable) buttons, etc highlighted
2959 // in this manner, unless we want to edit'em
2961 //When doing a DrawClone from the GUI you would do
2962 // - select an empty pad -
2963 // - right click on object -
2964 // - select DrawClone on menu -
2965 //
2966 // Without the SetSelectedPad(); in the HighLight function, the
2967 // above instruction lead to the clone to be drawn in the
2968 // same canvas as the original object. This is because the
2969 // 'right clicking' (via TCanvas::HandleInput) changes gPad
2970 // momentarily such that when DrawClone is called, it is
2971 // not the right value (for DrawClone). Should be FIXED.
2972 gROOT->SetSelectedPad(this);
2973 if (GetBorderMode()>0) {
2974 if (set) PaintBorder(-color, kFALSE);
2976 }
2977 }
2978
2980}
2981
2982////////////////////////////////////////////////////////////////////////////////
2983/// List all primitives in pad.
2984
2985void TPad::ls(Option_t *option) const
2986{
2988 std::cout <<IsA()->GetName()<<" fXlowNDC=" <<fXlowNDC<<" fYlowNDC="<<fYlowNDC<<" fWNDC="<<GetWNDC()<<" fHNDC="<<GetHNDC()
2989 <<" Name= "<<GetName()<<" Title= "<<GetTitle()<<" Option="<<option<<std::endl;
2991 if (!fPrimitives) return;
2992 fPrimitives->ls(option);
2994}
2995
2996////////////////////////////////////////////////////////////////////////////////
2997/// Increment (i==1) or set (i>1) the number of autocolor in the pad.
2998
3000{
3001 if (opt.Index("pfc")>=0 || opt.Index("plc")>=0 || opt.Index("pmc")>=0) {
3002 if (i==1) fNumPaletteColor++;
3003 else fNumPaletteColor = i;
3004 return fNumPaletteColor;
3005 } else {
3006 return 0;
3007 }
3008}
3009
3010////////////////////////////////////////////////////////////////////////////////
3011/// Get the next autocolor in the pad.
3012
3014{
3015 Int_t i = 0;
3016 Int_t ncolors = gStyle->GetNumberOfColors();
3017 if (fNumPaletteColor>1) {
3018 i = fNextPaletteColor*(ncolors/(fNumPaletteColor-1));
3019 if (i>=ncolors) i = ncolors-1;
3020 }
3023 return gStyle->GetColorPalette(i);
3024}
3025
3026////////////////////////////////////////////////////////////////////////////////
3027/// Initialise the grid used to find empty space when adding a box (Legend) in a pad
3028
3030{
3031 Int_t const cellSize = 10; // Sive of an individual grid cell in pixels.
3032
3033 if (fCGnx == 0 && fCGny == 0) {
3034 fCGnx = (Int_t)(gPad->GetWw())/cellSize;
3035 fCGny = (Int_t)(gPad->GetWh())/cellSize;
3036 } else {
3037 Int_t CGnx = (Int_t)(gPad->GetWw())/cellSize;
3038 Int_t CGny = (Int_t)(gPad->GetWh())/cellSize;
3039 if (fCGnx != CGnx || fCGny != CGny) {
3040 fCGnx = CGnx;
3041 fCGny = CGny;
3042 delete [] fCollideGrid;
3043 fCollideGrid = nullptr;
3044 }
3045 }
3046
3047 // Initialise the collide grid
3048 if (!fCollideGrid) {
3050 for (int i = 0; i<fCGnx; i++) {
3051 for (int j = 0; j<fCGny; j++) {
3052 fCollideGrid[i + j*fCGnx] = kTRUE;
3053 }
3054 }
3055 }
3056
3057 // Fill the collide grid
3059 if (!l) return;
3060 Int_t np = l->GetSize();
3061 TObject *o;
3062
3063 for (int i=0; i<np; i++) {
3064 o = (TObject *) l->At(i);
3065 if (o!=oi) {
3066 if (o->InheritsFrom(TFrame::Class())) { FillCollideGridTFrame(o); continue;}
3067 if (o->InheritsFrom(TBox::Class())) { FillCollideGridTBox(o); continue;}
3068 if (o->InheritsFrom(TH1::Class())) { FillCollideGridTH1(o); continue;}
3069 if (o->InheritsFrom(TGraph::Class())) { FillCollideGridTGraph(o); continue;}
3070 if (o->InheritsFrom(TMultiGraph::Class())) {
3071 TList * grlist = ((TMultiGraph *)o)->GetListOfGraphs();
3072 TIter nextgraph(grlist);
3073 TObject * og;
3074 while ((og = nextgraph())) FillCollideGridTGraph(og);
3075 }
3076 if (o->InheritsFrom(THStack::Class())) {
3077 TList * hlist = ((THStack *)o)->GetHists();
3078 TIter nexthist(hlist);
3079 TObject * oh;
3080 while ((oh = nexthist())) {
3082 }
3083 }
3084 }
3085 }
3086}
3087
3088////////////////////////////////////////////////////////////////////////////////
3089/// Check if a box of size w and h collide some primitives in the pad at
3090/// position i,j
3091
3093{
3094 for (int r=i; r<w+i; r++) {
3095 for (int c=j; c<h+j; c++) {
3096 if (!fCollideGrid[r + c*fCGnx]) return kTRUE;
3097 }
3098 }
3099 return kFALSE;
3100}
3101
3102////////////////////////////////////////////////////////////////////////////////
3103/// Place a box in NDC space
3104///
3105/// \return `true` if the box could be placed, `false` if not.
3106///
3107/// \param[in] o pointer to the box to be placed
3108/// \param[in] w box width to be placed
3109/// \param[in] h box height to be placed
3110/// \param[out] xl x position of the bottom left corner of the placed box
3111/// \param[out] yb y position of the bottom left corner of the placed box
3112
3114{
3115 FillCollideGrid(o);
3116
3117 Int_t iw = (int)(fCGnx*w);
3118 Int_t ih = (int)(fCGny*h);
3119
3120 Int_t nxmax = fCGnx-iw-1;
3121 Int_t nymax = fCGny-ih-1;
3122
3123 for (Int_t i = 0; i<nxmax; i++) {
3124 for (Int_t j = 0; j<=nymax; j++) {
3125 if (Collide(i,j,iw,ih)) {
3126 continue;
3127 } else {
3128 xl = (Double_t)(i)/(Double_t)(fCGnx);
3129 yb = (Double_t)(j)/(Double_t)(fCGny);
3130 return kTRUE;
3131 }
3132 }
3133 }
3134 return kFALSE;
3135}
3136
3137#define NotFree(i, j) fCollideGrid[TMath::Max(TMath::Min(i+j*fCGnx,fCGnx*fCGny),0)] = kFALSE;
3138
3139////////////////////////////////////////////////////////////////////////////////
3140/// Mark as "not free" the cells along a line.
3141
3143{
3144 NotFree(x1, y1);
3145 NotFree(x2, y2);
3146 Int_t i, j, xt, yt;
3147
3148 // horizontal lines
3149 if (y1==y2) {
3150 for (i=x1+1; i<x2; i++) NotFree(i,y1);
3151 return;
3152 }
3153
3154 // vertical lines
3155 if (x1==x2) {
3156 for (i=y1+1; i<y2; i++) NotFree(x1,i);
3157 return;
3158 }
3159
3160 // other lines
3161 if (TMath::Abs(x2-x1)>TMath::Abs(y2-y1)) {
3162 if (x1>x2) {
3163 xt = x1; x1 = x2; x2 = xt;
3164 yt = y1; y1 = y2; y2 = yt;
3165 }
3166 for (i=x1+1; i<x2; i++) {
3167 j = (Int_t)((Double_t)(y2-y1)*(Double_t)((i-x1)/(Double_t)(x2-x1))+y1);
3168 NotFree(i,j);
3169 NotFree(i,(j+1));
3170 }
3171 } else {
3172 if (y1>y2) {
3173 yt = y1; y1 = y2; y2 = yt;
3174 xt = x1; x1 = x2; x2 = xt;
3175 }
3176 for (j=y1+1; j<y2; j++) {
3177 i = (Int_t)((Double_t)(x2-x1)*(Double_t)((j-y1)/(Double_t)(y2-y1))+x1);
3178 NotFree(i,j);
3179 NotFree((i+1),j);
3180 }
3181 }
3182}
3183
3184////////////////////////////////////////////////////////////////////////////////
3186{
3187 TBox *b = (TBox *)o;
3188 if (fCGnx==0||fCGny==0) return;
3189 Double_t xs = (fX2-fX1)/fCGnx;
3190 Double_t ys = (fY2-fY1)/fCGny;
3191
3192 Int_t x1 = (Int_t)((b->GetX1()-fX1)/xs);
3193 Int_t x2 = (Int_t)((b->GetX2()-fX1)/xs);
3194 Int_t y1 = (Int_t)((b->GetY1()-fY1)/ys);
3195 Int_t y2 = (Int_t)((b->GetY2()-fY1)/ys);
3196 for (int i = x1; i<=x2; i++) {
3197 for (int j = y1; j<=y2; j++) NotFree(i, j);
3198 }
3199}
3200
3201////////////////////////////////////////////////////////////////////////////////
3203{
3204 TFrame *f = (TFrame *)o;
3205 if (fCGnx==0||fCGny==0) return;
3206 Double_t xs = (fX2-fX1)/fCGnx;
3207 Double_t ys = (fY2-fY1)/fCGny;
3208
3209 Int_t x1 = (Int_t)((f->GetX1()-fX1)/xs);
3210 Int_t x2 = (Int_t)((f->GetX2()-fX1)/xs);
3211 Int_t y1 = (Int_t)((f->GetY1()-fY1)/ys);
3212 Int_t y2 = (Int_t)((f->GetY2()-fY1)/ys);
3213 Int_t i;
3214
3215 for (i = x1; i<=x2; i++) {
3216 NotFree(i, y1);
3217 NotFree(i, (y1-1));
3218 NotFree(i, (y1-2));
3219 }
3220 for (i = y1; i<=y2; i++) {
3221 NotFree(x1, i);
3222 NotFree((x1-1), i);
3223 NotFree((x1-2), i);
3224 }
3225}
3226
3227////////////////////////////////////////////////////////////////////////////////
3229{
3230 TGraph *g = (TGraph *)o;
3231 if (fCGnx==0||fCGny==0) return;
3232 Double_t xs = (fX2-fX1)/fCGnx;
3233 Double_t ys = (fY2-fY1)/fCGny;
3234
3235 Int_t n = g->GetN();
3236 Int_t s = TMath::Max(n/10,1);
3237 Double_t x1, x2, y1, y2;
3238 for (Int_t i=s; i<n; i=i+s) {
3239 g->GetPoint(TMath::Max(0,i-s),x1,y1);
3240 g->GetPoint(i ,x2,y2);
3241 if (fLogx) {
3242 if (x1 > 0) x1 = TMath::Log10(x1);
3243 else x1 = fUxmin;
3244 if (x2 > 0) x2 = TMath::Log10(x2);
3245 else x2 = fUxmin;
3246 }
3247 if (fLogy) {
3248 if (y1 > 0) y1 = TMath::Log10(y1);
3249 else y1 = fUymin;
3250 if (y2 > 0) y2 = TMath::Log10(y2);
3251 else y2 = fUymin;
3252 }
3253 LineNotFree((int)((x1-fX1)/xs), (int)((x2-fX1)/xs),
3254 (int)((y1-fY1)/ys), (int)((y2-fY1)/ys));
3255 }
3256}
3257
3258////////////////////////////////////////////////////////////////////////////////
3260{
3261 TH1 *h = (TH1 *)o;
3262 if (fCGnx==0||fCGny==0) return;
3263 if (o->InheritsFrom(TH2::Class())) return;
3264 if (o->InheritsFrom(TH3::Class())) return;
3265
3266 TString name = h->GetName();
3267 if (name.Index("hframe") >= 0) return;
3268
3269 Double_t xs = (fX2-fX1)/fCGnx;
3270 Double_t ys = (fY2-fY1)/fCGny;
3271
3272 bool haserrors = false;
3273 TString drawOption = h->GetDrawOption();
3274 drawOption.ToLower();
3275 drawOption.ReplaceAll("same","");
3276
3277 if (drawOption.Index("hist") < 0) {
3278 if (drawOption.Index("e") >= 0) haserrors = true;
3279 }
3280
3281 Int_t nx = h->GetNbinsX();
3282 Int_t x1, y1, y2;
3283 Int_t i, j;
3284 Double_t x1l, y1l, y2l;
3285
3286 for (i = 1; i<nx; i++) {
3287 if (haserrors) {
3288 x1l = h->GetBinCenter(i);
3289 if (fLogx) {
3290 if (x1l > 0) x1l = TMath::Log10(x1l);
3291 else x1l = fUxmin;
3292 }
3293 x1 = (Int_t)((x1l-fX1)/xs);
3294 y1l = h->GetBinContent(i)-h->GetBinErrorLow(i);
3295 if (fLogy) {
3296 if (y1l > 0) y1l = TMath::Log10(y1l);
3297 else y1l = fUymin;
3298 }
3299 y1 = (Int_t)((y1l-fY1)/ys);
3300 y2l = h->GetBinContent(i)+h->GetBinErrorUp(i);
3301 if (fLogy) {
3302 if (y2l > 0) y2l = TMath::Log10(y2l);
3303 else y2l = fUymin;
3304 }
3305 y2 = (Int_t)((y2l-fY1)/ys);
3306 for (j=y1; j<=y2; j++) {
3307 NotFree(x1, j);
3308 }
3309 }
3310 x1l = h->GetBinLowEdge(i);
3311 if (fLogx) {
3312 if (x1l > 0) x1l = TMath::Log10(x1l);
3313 else x1l = fUxmin;
3314 }
3315 x1 = (Int_t)((x1l-fX1)/xs);
3316 y1l = h->GetBinContent(i);
3317 if (fLogy) {
3318 if (y1l > 0) y1l = TMath::Log10(y1l);
3319 else y1l = fUymin;
3320 }
3321 y1 = (Int_t)((y1l-fY1)/ys);
3322 NotFree(x1, y1);
3323 x1l = h->GetBinLowEdge(i)+h->GetBinWidth(i);
3324 if (fLogx) {
3325 if (x1l > 0) x1l = TMath::Log10(x1l);
3326 else x1l = fUxmin;
3327 }
3328 x1 = (int)((x1l-fX1)/xs);
3329 NotFree(x1, y1);
3330 }
3331
3332 // Extra objects in the list of function
3333 TPaveStats *ps = (TPaveStats*)h->GetListOfFunctions()->FindObject("stats");
3335}
3336
3337////////////////////////////////////////////////////////////////////////////////
3338/// This method draws the collide grid on top of the canvas. This is used for
3339/// debugging only. At some point it will be removed.
3340
3342{
3343 if (fCGnx==0||fCGny==0) return;
3344 auto box = new TBox();
3345 box->SetFillColorAlpha(kRed,0.5);
3346
3347 Double_t xs = (fX2-fX1)/fCGnx;
3348 Double_t ys = (fY2-fY1)/fCGny;
3349
3350 Double_t X1L, X2L, Y1L, Y2L;
3351 Double_t t = 0.15;
3352 Double_t Y1, Y2;
3353 Double_t X1 = fX1;
3354 Double_t X2 = X1+xs;
3355
3356 for (int i = 0; i<fCGnx; i++) {
3357 Y1 = fY1;
3358 Y2 = Y1+ys;
3359 for (int j = 0; j<fCGny; j++) {
3360 if (gPad->GetLogx()) {
3361 X1L = TMath::Power(10,X1);
3362 X2L = TMath::Power(10,X2);
3363 } else {
3364 X1L = X1;
3365 X2L = X2;
3366 }
3367 if (gPad->GetLogy()) {
3368 Y1L = TMath::Power(10,Y1);
3369 Y2L = TMath::Power(10,Y2);
3370 } else {
3371 Y1L = Y1;
3372 Y2L = Y2;
3373 }
3374 if (!fCollideGrid[i + j*fCGnx]) {
3375 box->SetFillColorAlpha(kBlack,t);
3376 box->DrawBox(X1L, Y1L, X2L, Y2L);
3377 } else {
3378 box->SetFillColorAlpha(kRed,t);
3379 box->DrawBox(X1L, Y1L, X2L, Y2L);
3380 }
3381 Y1 = Y2;
3382 Y2 = Y1+ys;
3383 if (t==0.15) t = 0.1;
3384 else t = 0.15;
3385 }
3386 X1 = X2;
3387 X2 = X1+xs;
3388 }
3389}
3390
3391
3392////////////////////////////////////////////////////////////////////////////////
3393/// Convert x from pad to X.
3394
3396{
3397 if (fLogx && x < 50) return Double_t(TMath::Exp(2.302585092994*x));
3398 return x;
3399}
3400
3401////////////////////////////////////////////////////////////////////////////////
3402/// Convert y from pad to Y.
3403
3405{
3406 if (fLogy && y < 50) return Double_t(TMath::Exp(2.302585092994*y));
3407 return y;
3408}
3409
3410////////////////////////////////////////////////////////////////////////////////
3411/// Convert x from X to pad.
3412
3414{
3415 if (fLogx) {
3416 if (x > 0) x = TMath::Log10(x);
3417 else x = fUxmin;
3418 }
3419 return x;
3420}
3421
3422////////////////////////////////////////////////////////////////////////////////
3423/// Convert y from Y to pad.
3424
3426{
3427 if (fLogy) {
3428 if (y > 0) y = TMath::Log10(y);
3429 else y = fUymin;
3430 }
3431 return y;
3432}
3433
3434////////////////////////////////////////////////////////////////////////////////
3435/// Paint all primitives in pad.
3436
3437void TPad::Paint(Option_t * /*option*/)
3438{
3439 if (!fPrimitives) fPrimitives = new TList;
3441 fViewer3D->PadPaint(this);
3443 if (GetGLDevice()!=-1 && gVirtualPS) {
3444 TPad *padsav = (TPad*)gPad;
3445 gPad = this;
3446 if (gGLManager) gGLManager->PrintViewer(GetViewer3D());
3447 gPad = padsav;
3448 }
3449 return;
3450 }
3451
3453
3454 TPad *padsav = (TPad*)gPad;
3455
3456 fPadPaint = 1;
3457 cd();
3458
3460 PaintDate();
3461
3463 TObject *obj;
3464
3465 Bool_t began3DScene = kFALSE;
3466 while (lnk) {
3467 obj = lnk->GetObject();
3468
3469 // Create a pad 3D viewer if none exists and we encounter a 3D shape
3470 if (!fViewer3D && obj->InheritsFrom(TAtt3D::Class())) {
3471 GetViewer3D("pad");
3472 }
3473
3474 // Open a 3D scene if required
3475 if (fViewer3D && !fViewer3D->BuildingScene()) {
3477 began3DScene = kTRUE;
3478 }
3479
3480 obj->Paint(lnk->GetOption());
3481 lnk = (TObjOptLink*)lnk->Next();
3482 }
3483
3484 if (padsav) padsav->cd();
3485 fPadPaint = 0;
3487
3488 // Close the 3D scene if we opened it. This must be done after modified
3489 // flag is cleared, as some viewers will invoke another paint by marking pad modified again
3490 if (began3DScene) {
3492 }
3493}
3494
3495////////////////////////////////////////////////////////////////////////////////
3496/// Paint the pad border.
3497/// Draw first a box as a normal filled box
3498
3500{
3501 if (color >= 0) {
3502 TAttLine::Modify(); //Change line attributes only if necessary
3503 TAttFill::Modify(); //Change fill area attributes only if necessary
3504
3505 //With Cocoa we have a transparency. But we also have
3506 //pixmaps, and if you just paint a new content over the old one
3507 //with alpha < 1., you'll be able to see the old content.
3508 if (!gROOT->IsBatch() && gVirtualX->InheritsFrom("TGCocoa") && GetPainter())
3510
3512 }
3513 if (color < 0) color = -color;
3514 // then paint 3d frame (depending on bordermode)
3515 if (IsTransparent()) return;
3516 // Paint a 3D frame around the pad.
3517
3518 if (fBorderMode == 0) return;
3519 Int_t bordersize = fBorderSize;
3520 if (bordersize <= 0) bordersize = 2;
3521
3522 const Double_t realBsX = bordersize / (GetAbsWNDC() * GetWw()) * (fX2 - fX1);
3523 const Double_t realBsY = bordersize / (GetAbsHNDC() * GetWh()) * (fY2 - fY1);
3524
3525 Short_t px1,py1,px2,py2;
3526 Double_t xl, xt, yl, yt;
3527
3528 // GetDarkColor() and GetLightColor() use GetFillColor()
3529 Color_t oldcolor = GetFillColor();
3530 SetFillColor(color);
3532 Color_t light = 0, dark = 0;
3533 if (color != 0) {
3534 light = TColor::GetColorBright(color);
3535 dark = TColor::GetColorDark(color);
3536 }
3537
3538 // Compute real left bottom & top right of the box in pixels
3539 px1 = XtoPixel(fX1); py1 = YtoPixel(fY1);
3540 px2 = XtoPixel(fX2); py2 = YtoPixel(fY2);
3541 if (px1 < px2) {xl = fX1; xt = fX2; }
3542 else {xl = fX2; xt = fX1;}
3543 if (py1 > py2) {yl = fY1; yt = fY2;}
3544 else {yl = fY2; yt = fY1;}
3545
3546 Double_t frameXs[7] = {}, frameYs[7] = {};
3547
3548 if (!IsBatch() && GetPainter()) {
3549 // Draw top&left part of the box
3550 frameXs[0] = xl; frameYs[0] = yl;
3551 frameXs[1] = xl + realBsX; frameYs[1] = yl + realBsY;
3552 frameXs[2] = frameXs[1]; frameYs[2] = yt - realBsY;
3553 frameXs[3] = xt - realBsX; frameYs[3] = frameYs[2];
3554 frameXs[4] = xt; frameYs[4] = yt;
3555 frameXs[5] = xl; frameYs[5] = yt;
3556 frameXs[6] = xl; frameYs[6] = yl;
3557
3558 if (fBorderMode == -1) GetPainter()->SetFillColor(dark);
3559 else GetPainter()->SetFillColor(light);
3560 GetPainter()->DrawFillArea(7, frameXs, frameYs);
3561
3562 // Draw bottom&right part of the box
3563 frameXs[0] = xl; frameYs[0] = yl;
3564 frameXs[1] = xl + realBsX; frameYs[1] = yl + realBsY;
3565 frameXs[2] = xt - realBsX; frameYs[2] = frameYs[1];
3566 frameXs[3] = frameXs[2]; frameYs[3] = yt - realBsY;
3567 frameXs[4] = xt; frameYs[4] = yt;
3568 frameXs[5] = xt; frameYs[5] = yl;
3569 frameXs[6] = xl; frameYs[6] = yl;
3570
3571 if (fBorderMode == -1) GetPainter()->SetFillColor(light);
3572 else GetPainter()->SetFillColor(dark);
3573 GetPainter()->DrawFillArea(7, frameXs, frameYs);
3574
3575 // If this pad is a button, highlight it
3576 if (InheritsFrom(TButton::Class()) && fBorderMode == -1) {
3577 if (TestBit(kFraming)) { // bit set in TButton::SetFraming
3578 if (GetFillColor() != 2) GetPainter()->SetLineColor(2);
3579 else GetPainter()->SetLineColor(4);
3580 GetPainter()->DrawBox(xl + realBsX, yl + realBsY, xt - realBsX, yt - realBsY, TVirtualPadPainter::kHollow);
3581 }
3582 }
3583 GetPainter()->SetFillColor(-1);
3584 SetFillColor(oldcolor);
3585 }
3586
3587 if (!tops) return;
3588
3589 PaintBorderPS(xl, yl, xt, yt, fBorderMode, bordersize, dark, light);
3590}
3591
3592////////////////////////////////////////////////////////////////////////////////
3593/// Paint a frame border with Postscript.
3594
3596{
3597 if (!gVirtualPS) return;
3598 gVirtualPS->DrawFrame(xl, yl, xt, yt, bmode,bsize,dark,light);
3599}
3600
3601////////////////////////////////////////////////////////////////////////////////
3602/// Paint the current date and time if the option date is on.
3603
3605{
3606 if (fCanvas == this && gStyle->GetOptDate()) {
3607 TDatime dt;
3608 const char *dates;
3609 char iso[16];
3610 if (gStyle->GetOptDate() < 10) {
3611 //by default use format like "Wed Sep 25 17:10:35 2002"
3612 dates = dt.AsString();
3613 } else if (gStyle->GetOptDate() < 20) {
3614 //use ISO format like 2002-09-25
3615 strlcpy(iso,dt.AsSQLString(),16);
3616 dates = iso;
3617 } else {
3618 //use ISO format like 2002-09-25 17:10:35
3619 dates = dt.AsSQLString();
3620 }
3621 TText tdate(gStyle->GetDateX(),gStyle->GetDateY(),dates);
3627 tdate.SetNDC();
3628 tdate.Paint();
3629 }
3630}
3631
3632////////////////////////////////////////////////////////////////////////////////
3633/// Paint histogram/graph frame.
3634
3636{
3637 if (!fPrimitives) fPrimitives = new TList;
3638 TList *glist = GetListOfPrimitives();
3639 TFrame *frame = GetFrame();
3640 frame->SetX1(xmin);
3641 frame->SetX2(xmax);
3642 frame->SetY1(ymin);
3643 frame->SetY2(ymax);
3644 if (!glist->FindObject(fFrame)) {
3645 glist->AddFirst(frame);
3647 }
3648 frame->Paint();
3649}
3650
3651////////////////////////////////////////////////////////////////////////////////
3652/// Traverse pad hierarchy and (re)paint only modified pads.
3653
3655{
3657 if (IsModified()) {
3658 fViewer3D->PadPaint(this);
3660 }
3661 TList *pList = GetListOfPrimitives();
3662 TObjOptLink *lnk = 0;
3663 if (pList) lnk = (TObjOptLink*)pList->FirstLink();
3664 TObject *obj;
3665 while (lnk) {
3666 obj = lnk->GetObject();
3667 if (obj->InheritsFrom(TPad::Class()))
3668 ((TPad*)obj)->PaintModified();
3669 lnk = (TObjOptLink*)lnk->Next();
3670 }
3671 return;
3672 }
3673
3675
3676 TPad *padsav = (TPad*)gPad;
3677 TVirtualPS *saveps = gVirtualPS;
3678 if (gVirtualPS) {
3680 }
3681 fPadPaint = 1;
3682 cd();
3683 if (IsModified() || IsTransparent()) {
3684 if ((fFillStyle < 3026) && (fFillStyle > 3000)) {
3685 if (!gPad->IsBatch() && GetPainter()) GetPainter()->ClearDrawable();
3686 }
3688 }
3689
3690 PaintDate();
3691
3692 TList *pList = GetListOfPrimitives();
3693 TObjOptLink *lnk = 0;
3694 if (pList) lnk = (TObjOptLink*)pList->FirstLink();
3695 TObject *obj;
3696
3697 Bool_t began3DScene = kFALSE;
3698
3699 while (lnk) {
3700 obj = lnk->GetObject();
3701 if (obj->InheritsFrom(TPad::Class())) {
3702 ((TPad*)obj)->PaintModified();
3703 } else if (IsModified() || IsTransparent()) {
3704
3705 // Create a pad 3D viewer if none exists and we encounter a
3706 // 3D shape
3707 if (!fViewer3D && obj->InheritsFrom(TAtt3D::Class())) {
3708 GetViewer3D("pad");
3709 }
3710
3711 // Open a 3D scene if required
3712 if (fViewer3D && !fViewer3D->BuildingScene()) {
3714 began3DScene = kTRUE;
3715 }
3716
3717 obj->Paint(lnk->GetOption());
3718 }
3719 lnk = (TObjOptLink*)lnk->Next();
3720 }
3721
3722 if (padsav) padsav->cd();
3723 fPadPaint = 0;
3725
3726 // This must be done after modified flag is cleared, as some
3727 // viewers will invoke another paint by marking pad modified again
3728 if (began3DScene) {
3730 }
3731
3732 gVirtualPS = saveps;
3733}
3734
3735////////////////////////////////////////////////////////////////////////////////
3736/// Paint box in CurrentPad World coordinates.
3737///
3738/// - if option[0] = 's' the box is forced to be paint with style=0
3739/// - if option[0] = 'l' the box contour is drawn
3740
3742{
3743 if (!gPad->IsBatch() && GetPainter()) {
3744 Int_t style0 = GetPainter()->GetFillStyle();
3745 Int_t style = style0;
3746 if (option[0] == 's') {
3748 style = 0;
3749 }
3750 if (style) {
3751 if (style > 3000 && style < 4000) {
3752 if (style < 3026) {
3753 // draw stipples with fFillColor foreground
3755 }
3756
3757 if (style >= 3100 && style < 4000) {
3758 Double_t xb[4], yb[4];
3759 xb[0] = x1; xb[1] = x1; xb[2] = x2; xb[3] = x2;
3760 yb[0] = y1; yb[1] = y2; yb[2] = y2; yb[3] = y1;
3761 PaintFillAreaHatches(4, xb, yb, style);
3762 return;
3763 }
3764 //special case for TAttFillCanvas
3765 if (GetPainter()->GetFillColor() == 10) {
3768 GetPainter()->SetFillColor(10);
3769 }
3770 } else if (style >= 4000 && style <= 4100) {
3771 // For style >=4000 we make the window transparent.
3772 // From 4000 to 4100 the window is 100% transparent to 100% opaque
3773
3774 //ignore this style option when this is the canvas itself
3775 if (this == fMother) {
3776 //It's clear, that virtual X checks a style (4000) and will render a hollow rect!
3777 const Style_t oldFillStyle = GetPainter()->GetFillStyle();
3778 if (gVirtualX->InheritsFrom("TGCocoa"))
3779 GetPainter()->SetFillStyle(1000);
3781 if (gVirtualX->InheritsFrom("TGCocoa"))
3782 GetPainter()->SetFillStyle(oldFillStyle);
3783 } else {
3784 //draw background by blitting all bottom pads
3785 int px, py;
3786 XYtoAbsPixel(fX1, fY2, px, py);
3787
3788 if (fMother) {
3790 CopyBackgroundPixmaps(fMother, this, px, py);
3791 }
3792
3793 GetPainter()->SetOpacity(style - 4000);
3794 }
3795 } else if (style >= 1000 && style <= 1999) {
3797 } else {
3799 }
3800 if (option[0] == 'l') GetPainter()->DrawBox(x1, y1, x2, y2, TVirtualPadPainter::kHollow);
3801 } else {
3803 if (option[0] == 's') GetPainter()->SetFillStyle(style0);
3804 }
3805 }
3806
3807 if (gVirtualPS) {
3808 Int_t style0 = gVirtualPS->GetFillStyle();
3809 if (option[0] == 's') {
3811 } else {
3812 if (style0 >= 3100 && style0 < 4000) {
3813 Double_t xb[4], yb[4];
3814 xb[0] = x1; xb[1] = x1; xb[2] = x2; xb[3] = x2;
3815 yb[0] = y1; yb[1] = y2; yb[2] = y2; yb[3] = y1;
3816 PaintFillAreaHatches(4, xb, yb, style0);
3817 return;
3818 }
3819 }
3820 gVirtualPS->DrawBox(x1, y1, x2, y2);
3821 if (option[0] == 'l') {
3823 gVirtualPS->DrawBox(x1, y1, x2, y2);
3824 }
3825 if (option[0] == 's' || option[0] == 'l') gVirtualPS->SetFillStyle(style0);
3826 }
3827
3828 Modified();
3829}
3830
3831////////////////////////////////////////////////////////////////////////////////
3832/// Copy pixmaps of pads laying below pad "stop" into pad "stop". This
3833/// gives the effect of pad "stop" being transparent.
3834
3836{
3837 if (!start) return;
3838 TObject *obj;
3839 if (!fPrimitives) fPrimitives = new TList;
3840 TIter next(start->GetListOfPrimitives());
3841 while ((obj = next())) {
3842 if (obj->InheritsFrom(TPad::Class())) {
3843 if (obj == stop) break;
3844 ((TPad*)obj)->CopyBackgroundPixmap(x, y);
3845 ((TPad*)obj)->CopyBackgroundPixmaps((TPad*)obj, stop, x, y);
3846 }
3847 }
3848}
3849
3850////////////////////////////////////////////////////////////////////////////////
3851/// Copy pixmap of this pad as background of the current pad.
3852
3854{
3855 int px, py;
3856 XYtoAbsPixel(fX1, fY2, px, py);
3857 if (GetPainter()) GetPainter()->CopyDrawable(GetPixmapID(), px-x, py-y);
3858}
3859
3860////////////////////////////////////////////////////////////////////////////////
3861
3863{
3864 Warning("TPad::PaintFillArea", "Float_t signature is obsolete. Use Double_t signature.");
3865}
3866
3867////////////////////////////////////////////////////////////////////////////////
3868/// Paint fill area in CurrentPad World coordinates.
3869
3871{
3872 if (nn <3) return;
3873 Int_t n=0;
3877 } else {
3878 xmin = fX1; ymin = fY1; xmax = fX2; ymax = fY2;
3879 }
3880
3881 Int_t nc = 2*nn+1;
3882 std::vector<Double_t> x(nc, 0.);
3883 std::vector<Double_t> y(nc, 0.);
3884
3885 n = ClipPolygon(nn, xx, yy, nc, &x.front(), &y.front(),xmin,ymin,xmax,ymax);
3886 if (!n)
3887 return;
3888
3889 // Paint the fill area with hatches
3890 Int_t fillstyle = GetPainter()?GetPainter()->GetFillStyle():1;
3891 if (gPad->IsBatch() && GetPainter() && gVirtualPS) fillstyle = gVirtualPS->GetFillStyle();
3892 if (fillstyle >= 3100 && fillstyle < 4000) {
3893 PaintFillAreaHatches(nn, &x.front(), &y.front(), fillstyle);
3894 return;
3895 }
3896
3897 if (!gPad->IsBatch() && GetPainter())
3898 // invoke the graphics subsystem
3899 GetPainter()->DrawFillArea(n, &x.front(), &y.front());
3900
3901 if (gVirtualPS)
3902 gVirtualPS->DrawPS(-n, &x.front(), &y.front());
3903
3904 Modified();
3905}
3906
3907////////////////////////////////////////////////////////////////////////////////
3908/// Paint fill area in CurrentPad NDC coordinates.
3909
3911{
3912 auto xw = new Double_t[n];
3913 auto yw = new Double_t[n];
3914 for (int i=0; i<n; i++) {
3915 xw[i] = fX1 + x[i]*(fX2 - fX1);
3916 yw[i] = fY1 + y[i]*(fY2 - fY1);
3917 }
3918 PaintFillArea(n, xw, yw, option);
3919 delete [] xw;
3920 delete [] yw;
3921}
3922
3923////////////////////////////////////////////////////////////////////////////////
3924/// This function paints hatched fill area according to the FillStyle value
3925/// The convention for the Hatch is the following:
3926///
3927/// `FillStyle = 3ijk`
3928///
3929/// - i (1-9) : specify the space between each hatch
3930/// 1 = minimum 9 = maximum
3931/// the final spacing is i*GetHatchesSpacing(). The hatches spacing
3932/// is set by SetHatchesSpacing()
3933/// - j (0-9) : specify angle between 0 and 90 degrees
3934/// * 0 = 0
3935/// * 1 = 10
3936/// * 2 = 20
3937/// * 3 = 30
3938/// * 4 = 45
3939/// * 5 = Not drawn
3940/// * 6 = 60
3941/// * 7 = 70
3942/// * 8 = 80
3943/// * 9 = 90
3944/// - k (0-9) : specify angle between 90 and 180 degrees
3945/// * 0 = 180
3946/// * 1 = 170
3947/// * 2 = 160
3948/// * 3 = 150
3949/// * 4 = 135
3950/// * 5 = Not drawn
3951/// * 6 = 120
3952/// * 7 = 110
3953/// * 8 = 100
3954/// * 9 = 90
3955
3957{
3958 static Double_t ang1[10] = { 0., 10., 20., 30., 45.,5., 60., 70., 80., 89.99};
3959 static Double_t ang2[10] = {180.,170.,160.,150.,135.,5.,120.,110.,100., 89.99};
3960
3961 Int_t fasi = FillStyle%1000;
3962 Int_t idSPA = (Int_t)(fasi/100);
3963 Int_t iAng2 = (Int_t)((fasi-100*idSPA)/10);
3964 Int_t iAng1 = fasi%10;
3965 Double_t dy = 0.003*(Double_t)(idSPA)*gStyle->GetHatchesSpacing();
3967 Short_t lws = 0;
3968 Int_t lss = 0;
3969 Int_t lcs = 0;
3970
3971 // Save the current line attributes
3972 if (!gPad->IsBatch() && GetPainter()) {
3973 lws = GetPainter()->GetLineWidth();
3974 lss = GetPainter()->GetLineStyle();
3975 lcs = GetPainter()->GetLineColor();
3976 } else {
3977 if (gVirtualPS) {
3978 lws = gVirtualPS->GetLineWidth();
3979 lss = gVirtualPS->GetLineStyle();
3980 lcs = gVirtualPS->GetLineColor();
3981 }
3982 }
3983
3984 // Change the current line attributes to draw the hatches
3985 if (!gPad->IsBatch() && GetPainter()) {
3989 }
3990 if (gVirtualPS) {
3994 }
3995
3996 // Draw the hatches
3997 if (ang1[iAng1] != 5.) PaintHatches(dy, ang1[iAng1], nn, xx, yy);
3998 if (ang2[iAng2] != 5.) PaintHatches(dy, ang2[iAng2], nn, xx, yy);
3999
4000 // Restore the line attributes
4001 if (!gPad->IsBatch() && GetPainter()) {
4002 GetPainter()->SetLineStyle(lss);
4003 GetPainter()->SetLineWidth(lws);
4004 GetPainter()->SetLineColor(lcs);
4005 }
4006 if (gVirtualPS) {
4010 }
4011}
4012
4013////////////////////////////////////////////////////////////////////////////////
4014/// This routine draw hatches inclined with the
4015/// angle "angle" and spaced of "dy" in normalized device
4016/// coordinates in the surface defined by n,xx,yy.
4017
4019 Int_t nn, Double_t *xx, Double_t *yy)
4020{
4021 Int_t i, i1, i2, nbi, m, inv;
4022 Double_t ratiox, ratioy, ymin, ymax, yrot, ycur;
4023 const Double_t angr = TMath::Pi()*(180.-angle)/180.;
4024 const Double_t epsil = 0.0001;
4025 const Int_t maxnbi = 100;
4026 Double_t xli[maxnbi], xlh[2], ylh[2], xt1, xt2, yt1, yt2;
4027 Double_t ll, x, y, x1, x2, y1, y2, a, b, xi, xip, xin, yi, yip;
4028
4029 Double_t rwxmin = gPad->GetX1();
4030 Double_t rwxmax = gPad->GetX2();
4031 Double_t rwymin = gPad->GetY1();
4032 Double_t rwymax = gPad->GetY2();
4033 ratiox = 1./(rwxmax-rwxmin);
4034 ratioy = 1./(rwymax-rwymin);
4035
4036 Double_t sina = TMath::Sin(angr), sinb;
4037 Double_t cosa = TMath::Cos(angr), cosb;
4038 if (TMath::Abs(cosa) <= epsil) cosa=0.;
4039 if (TMath::Abs(sina) <= epsil) sina=0.;
4040 sinb = -sina;
4041 cosb = cosa;
4042
4043 // Values needed to compute the hatches in TRUE normalized space (NDC)
4044 Int_t iw = (Int_t)gPad->GetWw();
4045 Int_t ih = (Int_t)gPad->GetWh();
4046 Double_t x1p,y1p,x2p,y2p;
4047 gPad->GetPadPar(x1p,y1p,x2p,y2p);
4048 iw = (Int_t)(iw*x2p)-(Int_t)(iw*x1p);
4049 ih = (Int_t)(ih*y2p)-(Int_t)(ih*y1p);
4050 Double_t wndc = TMath::Min(1.,(Double_t)iw/(Double_t)ih);
4051 Double_t hndc = TMath::Min(1.,(Double_t)ih/(Double_t)iw);
4052
4053 // Search ymin and ymax
4054 ymin = 1.;
4055 ymax = 0.;
4056 for (i=1; i<=nn; i++) {
4057 x = wndc*ratiox*(xx[i-1]-rwxmin);
4058 y = hndc*ratioy*(yy[i-1]-rwymin);
4059 yrot = sina*x+cosa*y;
4060 if (yrot > ymax) ymax = yrot;
4061 if (yrot < ymin) ymin = yrot;
4062 }
4063 ymax = (Double_t)((Int_t)(ymax/dy))*dy;
4064
4065 for (ycur=ymax; ycur>=ymin; ycur=ycur-dy) {
4066 nbi = 0;
4067 for (i=2; i<=nn+1; i++) {
4068 i2 = i;
4069 i1 = i-1;
4070 if (i == nn+1) i2=1;
4071 x1 = wndc*ratiox*(xx[i1-1]-rwxmin);
4072 y1 = hndc*ratioy*(yy[i1-1]-rwymin);
4073 x2 = wndc*ratiox*(xx[i2-1]-rwxmin);
4074 y2 = hndc*ratioy*(yy[i2-1]-rwymin);
4075 xt1 = cosa*x1-sina*y1;
4076 yt1 = sina*x1+cosa*y1;
4077 xt2 = cosa*x2-sina*y2;
4078 yt2 = sina*x2+cosa*y2;
4079
4080 // Line segment parallel to oy
4081 if (xt1 == xt2) {
4082 if (yt1 < yt2) {
4083 yi = yt1;
4084 yip = yt2;
4085 } else {
4086 yi = yt2;
4087 yip = yt1;
4088 }
4089 if ((yi <= ycur) && (ycur < yip)) {
4090 nbi++;
4091 if (nbi >= maxnbi) return;
4092 xli[nbi-1] = xt1;
4093 }
4094 continue;
4095 }
4096
4097 // Line segment parallel to ox
4098 if (yt1 == yt2) {
4099 if (yt1 == ycur) {
4100 nbi++;
4101 if (nbi >= maxnbi) return;
4102 xli[nbi-1] = xt1;
4103 nbi++;
4104 if (nbi >= maxnbi) return;
4105 xli[nbi-1] = xt2;
4106 }
4107 continue;
4108 }
4109
4110 // Other line segment
4111 a = (yt1-yt2)/(xt1-xt2);
4112 b = (yt2*xt1-xt2*yt1)/(xt1-xt2);
4113 if (xt1 < xt2) {
4114 xi = xt1;
4115 xip = xt2;
4116 } else {
4117 xi = xt2;
4118 xip = xt1;
4119 }
4120 xin = (ycur-b)/a;
4121 if ((xi <= xin) && (xin < xip) &&
4122 (TMath::Min(yt1,yt2) <= ycur) &&
4123 (ycur < TMath::Max(yt1,yt2))) {
4124 nbi++;
4125 if (nbi >= maxnbi) return;
4126 xli[nbi-1] = xin;
4127 }
4128 }
4129
4130 // Sorting of the x coordinates intersections
4131 inv = 0;
4132 m = nbi-1;
4133L30:
4134 for (i=1; i<=m; i++) {
4135 if (xli[i] < xli[i-1]) {
4136 inv++;
4137 ll = xli[i-1];
4138 xli[i-1] = xli[i];
4139 xli[i] = ll;
4140 }
4141 }
4142 m--;
4143 if (inv == 0) goto L50;
4144 inv = 0;
4145 goto L30;
4146
4147 // Draw the hatches
4148L50:
4149 if (nbi%2 != 0) continue;
4150
4151 for (i=1; i<=nbi; i=i+2) {
4152 // Rotate back the hatches
4153 xlh[0] = cosb*xli[i-1]-sinb*ycur;
4154 ylh[0] = sinb*xli[i-1]+cosb*ycur;
4155 xlh[1] = cosb*xli[i] -sinb*ycur;
4156 ylh[1] = sinb*xli[i] +cosb*ycur;
4157 // Convert hatches' positions from true NDC to WC
4158 xlh[0] = (xlh[0]/wndc)*(rwxmax-rwxmin)+rwxmin;
4159 ylh[0] = (ylh[0]/hndc)*(rwymax-rwymin)+rwymin;
4160 xlh[1] = (xlh[1]/wndc)*(rwxmax-rwxmin)+rwxmin;
4161 ylh[1] = (ylh[1]/hndc)*(rwymax-rwymin)+rwymin;
4162 gPad->PaintLine(xlh[0], ylh[0], xlh[1], ylh[1]);
4163 }
4164 }
4165}
4166
4167////////////////////////////////////////////////////////////////////////////////
4168/// Paint line in CurrentPad World coordinates.
4169
4171{
4172 Double_t x[2], y[2];
4173 x[0] = x1; x[1] = x2; y[0] = y1; y[1] = y2;
4174
4175 //If line is totally clipped, return
4177 if (Clip(x,y,fUxmin,fUymin,fUxmax,fUymax) == 2) return;
4178 } else {
4179 if (Clip(x,y,fX1,fY1,fX2,fY2) == 2) return;
4180 }
4181
4182 if (!gPad->IsBatch() && GetPainter())
4183 GetPainter()->DrawLine(x[0], y[0], x[1], y[1]);
4184
4185 if (gVirtualPS) {
4186 gVirtualPS->DrawPS(2, x, y);
4187 }
4188
4189 Modified();
4190}
4191
4192////////////////////////////////////////////////////////////////////////////////
4193/// Paint line in normalized coordinates.
4194
4196{
4197 static Double_t xw[2], yw[2];
4198 if (!gPad->IsBatch() && GetPainter())
4199 GetPainter()->DrawLineNDC(u1, v1, u2, v2);
4200
4201 if (gVirtualPS) {
4202 xw[0] = fX1 + u1*(fX2 - fX1);
4203 xw[1] = fX1 + u2*(fX2 - fX1);
4204 yw[0] = fY1 + v1*(fY2 - fY1);
4205 yw[1] = fY1 + v2*(fY2 - fY1);
4206 gVirtualPS->DrawPS(2, xw, yw);
4207 }
4208
4209 Modified();
4210}
4211
4212////////////////////////////////////////////////////////////////////////////////
4213/// Paint 3-D line in the CurrentPad.
4214
4216{
4217 if (!fView) return;
4218
4219 // convert from 3-D to 2-D pad coordinate system
4220 Double_t xpad[6];
4221 Double_t temp[3];
4222 Int_t i;
4223 for (i=0;i<3;i++) temp[i] = p1[i];
4224 fView->WCtoNDC(temp, &xpad[0]);
4225 for (i=0;i<3;i++) temp[i] = p2[i];
4226 fView->WCtoNDC(temp, &xpad[3]);
4227 PaintLine(xpad[0],xpad[1],xpad[3],xpad[4]);
4228}
4229
4230////////////////////////////////////////////////////////////////////////////////
4231/// Paint 3-D line in the CurrentPad.
4232
4234{
4235 //take into account perspective view
4236 if (!fView) return;
4237 // convert from 3-D to 2-D pad coordinate system
4238 Double_t xpad[6];
4239 Double_t temp[3];
4240 Int_t i;
4241 for (i=0;i<3;i++) temp[i] = p1[i];
4242 fView->WCtoNDC(temp, &xpad[0]);
4243 for (i=0;i<3;i++) temp[i] = p2[i];
4244 fView->WCtoNDC(temp, &xpad[3]);
4245 PaintLine(xpad[0],xpad[1],xpad[3],xpad[4]);
4246}
4247
4248////////////////////////////////////////////////////////////////////////////////
4249/// Paint polyline in CurrentPad World coordinates.
4250
4252{
4253 if (n < 2) return;
4254
4258 } else {
4259 xmin = fX1; ymin = fY1; xmax = fX2; ymax = fY2;
4260 }
4261 Int_t i, i1=-1,np=1;
4262 for (i=0; i<n-1; i++) {
4263 Double_t x1=x[i];
4264 Double_t y1=y[i];
4265 Double_t x2=x[i+1];
4266 Double_t y2=y[i+1];
4267 Int_t iclip = Clip(&x[i],&y[i],xmin,ymin,xmax,ymax);
4268 if (iclip == 2) {
4269 i1 = -1;
4270 continue;
4271 }
4272 np++;
4273 if (i1 < 0) i1 = i;
4274 if (iclip == 0 && i < n-2) continue;
4275 if (!gPad->IsBatch() && GetPainter())
4276 GetPainter()->DrawPolyLine(np, &x[i1], &y[i1]);
4277 if (gVirtualPS) {
4278 gVirtualPS->DrawPS(np, &x[i1], &y[i1]);
4279 }
4280 if (iclip) {
4281 x[i] = x1;
4282 y[i] = y1;
4283 x[i+1] = x2;
4284 y[i+1] = y2;
4285 }
4286 i1 = -1;
4287 np = 1;
4288 }
4289
4290 Modified();
4291}
4292
4293////////////////////////////////////////////////////////////////////////////////
4294/// Paint polyline in CurrentPad World coordinates.
4295///
4296/// If option[0] == 'C' no clipping
4297
4299{
4300 if (n < 2) return;
4301
4303 Bool_t mustClip = kTRUE;
4306 } else {
4307 xmin = fX1; ymin = fY1; xmax = fX2; ymax = fY2;
4308 if (option && (option[0] == 'C')) mustClip = kFALSE;
4309 }
4310
4311 Int_t i, i1=-1, np=1, iclip=0;
4312
4313 for (i=0; i < n-1; i++) {
4314 Double_t x1=x[i];
4315 Double_t y1=y[i];
4316 Double_t x2=x[i+1];
4317 Double_t y2=y[i+1];
4318 if (mustClip) {
4319 iclip = Clip(&x[i],&y[i],xmin,ymin,xmax,ymax);
4320 if (iclip == 2) {
4321 i1 = -1;
4322 continue;
4323 }
4324 }
4325 np++;
4326 if (i1 < 0) i1 = i;
4327 if (iclip == 0 && i < n-2) continue;
4328 if (!gPad->IsBatch() && GetPainter())
4329 GetPainter()->DrawPolyLine(np, &x[i1], &y[i1]);
4330 if (gVirtualPS) {
4331 gVirtualPS->DrawPS(np, &x[i1], &y[i1]);
4332 }
4333 if (iclip) {
4334 x[i] = x1;
4335 y[i] = y1;
4336 x[i+1] = x2;
4337 y[i+1] = y2;
4338 }
4339 i1 = -1;
4340 np = 1;
4341 }
4342
4343 Modified();
4344}
4345
4346////////////////////////////////////////////////////////////////////////////////
4347/// Paint polyline in CurrentPad NDC coordinates.
4348
4350{
4351 if (n <=0) return;
4352
4353 if (!gPad->IsBatch() && GetPainter())
4355
4356 if (gVirtualPS) {
4357 Double_t *xw = new Double_t[n];
4358 Double_t *yw = new Double_t[n];
4359 for (Int_t i=0; i<n; i++) {
4360 xw[i] = fX1 + x[i]*(fX2 - fX1);
4361 yw[i] = fY1 + y[i]*(fY2 - fY1);
4362 }
4363 gVirtualPS->DrawPS(n, xw, yw);
4364 delete [] xw;
4365 delete [] yw;
4366 }
4367 Modified();
4368}
4369
4370////////////////////////////////////////////////////////////////////////////////
4371/// Paint 3-D polyline in the CurrentPad.
4372
4374{
4375 if (!fView) return;
4376
4377 // Loop on each individual line
4378 for (Int_t i = 1; i < n; i++)
4379 PaintLine3D(&p[3*i-3], &p[3*i]);
4380
4381 Modified();
4382}
4383
4384////////////////////////////////////////////////////////////////////////////////
4385/// Paint polymarker in CurrentPad World coordinates.
4386
4388{
4389 Int_t n = TMath::Abs(nn);
4391 if (nn > 0 || TestBit(TGraph::kClipFrame)) {
4393 } else {
4394 xmin = fX1; ymin = fY1; xmax = fX2; ymax = fY2;
4395 }
4396 Int_t i,i1=-1,np=0;
4397 for (i=0; i<n; i++) {
4398 if (x[i] >= xmin && x[i] <= xmax && y[i] >= ymin && y[i] <= ymax) {
4399 np++;
4400 if (i1 < 0) i1 = i;
4401 if (i < n-1) continue;
4402 }
4403 if (np == 0) continue;
4404 if (!gPad->IsBatch() && GetPainter())
4405 GetPainter()->DrawPolyMarker(np, &x[i1], &y[i1]);
4406 if (gVirtualPS) {
4407 gVirtualPS->DrawPolyMarker(np, &x[i1], &y[i1]);
4408 }
4409 i1 = -1;
4410 np = 0;
4411 }
4412 Modified();
4413}
4414
4415////////////////////////////////////////////////////////////////////////////////
4416/// Paint polymarker in CurrentPad World coordinates.
4417
4419{
4420 Int_t n = TMath::Abs(nn);
4422 if (nn > 0 || TestBit(TGraph::kClipFrame)) {
4424 } else {
4425 xmin = fX1; ymin = fY1; xmax = fX2; ymax = fY2;
4426 }
4427 Int_t i,i1=-1,np=0;
4428 for (i=0; i<n; i++) {
4429 if (x[i] >= xmin && x[i] <= xmax && y[i] >= ymin && y[i] <= ymax) {
4430 np++;
4431 if (i1 < 0) i1 = i;
4432 if (i < n-1) continue;
4433 }
4434 if (np == 0) continue;
4435 if (!gPad->IsBatch() && GetPainter())
4436 GetPainter()->DrawPolyMarker(np, &x[i1], &y[i1]);
4437 if (gVirtualPS) {
4438 gVirtualPS->DrawPolyMarker(np, &x[i1], &y[i1]);
4439 }
4440 i1 = -1;
4441 np = 0;
4442 }
4443 Modified();
4444}
4445
4446////////////////////////////////////////////////////////////////////////////////
4447/// Paint text in CurrentPad World coordinates.
4448
4450{
4451 Modified();
4452
4453 if (!gPad->IsBatch() && GetPainter())
4455
4456 if (gVirtualPS) gVirtualPS->Text(x, y, text);
4457}
4458
4459////////////////////////////////////////////////////////////////////////////////
4460/// Paint text in CurrentPad World coordinates.
4461
4462void TPad::PaintText(Double_t x, Double_t y, const wchar_t *text)
4463{
4464 Modified();
4465
4466 if (!gPad->IsBatch() && GetPainter())
4468
4469 if (gVirtualPS) gVirtualPS->Text(x, y, text);
4470}
4471
4472////////////////////////////////////////////////////////////////////////////////
4473/// Paint text in CurrentPad NDC coordinates.
4474
4476{
4477 Modified();
4478
4479 if (!gPad->IsBatch() && GetPainter())
4481
4482 if (gVirtualPS) {
4483 Double_t x = fX1 + u*(fX2 - fX1);
4484 Double_t y = fY1 + v*(fY2 - fY1);
4485 gVirtualPS->Text(x, y, text);
4486 }
4487}
4488
4489////////////////////////////////////////////////////////////////////////////////
4490/// Paint text in CurrentPad NDC coordinates.
4491
4493{
4494 Modified();
4495
4496 if (!gPad->IsBatch() && GetPainter())
4498
4499 if (gVirtualPS) {
4500 Double_t x = fX1 + u*(fX2 - fX1);
4501 Double_t y = fY1 + v*(fY2 - fY1);
4502 gVirtualPS->Text(x, y, text);
4503 }
4504}
4505
4506////////////////////////////////////////////////////////////////////////////////
4507/// Search for an object at pixel position px,py.
4508///
4509/// Check if point is in this pad.
4510///
4511/// If yes, check if it is in one of the sub-pads
4512///
4513/// If found in the pad, compute closest distance of approach
4514/// to each primitive.
4515///
4516/// If one distance of approach is found to be within the limit Distancemaximum
4517/// the corresponding primitive is selected and the routine returns.
4518
4520{
4521 //the two following statements are necessary under NT (multithreaded)
4522 //when a TCanvas object is being created and a thread calling TPad::Pick
4523 //before the TPad constructor has completed in the other thread
4524 if (gPad == 0) return 0; //Andy Haas
4525 if (GetListOfPrimitives() == 0) return 0; //Andy Haas
4526
4527 Int_t dist;
4528 // Search if point is in pad itself
4529 Double_t x = AbsPixeltoX(px);
4530 Double_t y = AbsPixeltoY(py);
4531 if (this != gPad->GetCanvas()) {
4532 if (!((x >= fX1 && x <= fX2) && (y >= fY1 && y <= fY2))) return 0;
4533 }
4534
4535 // search for a primitive in this pad or its sub-pads
4536 static TObjOptLink dummyLink(0,""); //place holder for when no link available
4537 TPad *padsav = (TPad*)gPad;
4538 gPad = this; // since no drawing will be done, don't use cd() for efficiency reasons
4539 TPad *pick = 0;
4540 TPad *picked = this;
4541 pickobj = 0;
4543 dummyLink.SetObject(this);
4544 pickobj = &dummyLink;
4545 }
4546
4547 // Loop backwards over the list of primitives. The first non-pad primitive
4548 // found is the selected one. However, we have to keep going down the
4549 // list to see if there is maybe a pad overlaying the primitive. In that
4550 // case look into the pad for a possible primitive. Once a pad has been
4551 // found we can terminate the loop.
4552 Bool_t gotPrim = kFALSE; // true if found a non pad primitive
4554
4555 //We can have 3d stuff in pad. If canvas prefers to draw
4556 //such stuff with OpenGL, the selection of 3d objects is
4557 //a gl viewer business so, in first cycle we do not
4558 //call DistancetoPrimitive for TAtt3D descendants.
4559 //In case of gl we first try to select 2d object first.
4560
4561 while (lnk) {
4562 TObject *obj = lnk->GetObject();
4563
4564 //If canvas prefers GL, all 3d objects must be drawn/selected by
4565 //gl viewer
4566 if (obj->InheritsFrom(TAtt3D::Class()) && fEmbeddedGL) {
4567 lnk = lnk->Prev();
4568 continue;
4569 }
4570
4571 fPadPointer = obj;
4572 if (obj->InheritsFrom(TPad::Class())) {
4573 pick = ((TPad*)obj)->Pick(px, py, pickobj);
4574 if (pick) {
4575 picked = pick;
4576 break;
4577 }
4578 } else if (!gROOT->GetEditorMode()) {
4579 if (!gotPrim) {
4580 if (!obj->TestBit(kCannotPick)) {
4581 dist = obj->DistancetoPrimitive(px, py);
4582 if (dist < fgMaxPickDistance) {
4583 pickobj = lnk;
4584 gotPrim = kTRUE;
4585 if (dist == 0) break;
4586 }
4587 }
4588 }
4589 }
4590
4591 lnk = lnk->Prev();
4592 }
4593
4594 //if no primitive found, check if we have a TView
4595 //if yes, return the view except if you are in the lower or upper X range
4596 //of the pad.
4597 //In case canvas prefers gl, fView existence
4598 //automatically means viewer3d existence. (?)
4599
4600 if (fView && !gotPrim) {
4601 Double_t dx = 0.05*(fUxmax-fUxmin);
4602 if ((x > fUxmin + dx) && (x < fUxmax-dx)) {
4603
4604 if (fEmbeddedGL) {
4605 //No 2d stuff was selected, but we have gl-viewer. Let it select an object in
4606 //scene (or select itself). In any case it'll internally call
4607 //gPad->SetSelected(ptr) as, for example, hist painter does.
4608 py -= Int_t((1 - GetHNDC() - GetYlowNDC()) * GetWh());
4609 px -= Int_t(GetXlowNDC() * GetWw());
4611 }
4612 else
4613 dummyLink.SetObject(fView);
4614 }
4615 }
4616
4617 if (picked->InheritsFrom(TButton::Class())) {
4618 TButton *button = (TButton*)picked;
4619 if (!button->IsEditable()) pickobj = 0;
4620 }
4621
4622 if (TestBit(kCannotPick)) {
4623
4624 if (picked == this) {