Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMarker3DBox.cxx
Go to the documentation of this file.
1// @(#)root/g3d:$Id$
2// Author: Rene Brun , Olivier Couet 31/10/97
3
4/*************************************************************************
5 * Copyright (C) 1995-2004, 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 <iostream>
13#include "TROOT.h"
14#include "TBuffer.h"
15#include "TView.h"
16#include "TMarker3DBox.h"
17#include "TVirtualPad.h"
18#include "TH1.h"
19#include "TBuffer3D.h"
20#include "TBuffer3DTypes.h"
21#include "TVirtualViewer3D.h"
22#include "TGeometry.h"
23#include "TMath.h"
24
25#include <cassert>
26
28
29/** \class TMarker3DBox
30\ingroup g3d
31A special 3-D marker designed for event display.
32
33It has the following parameters:
34 - fX: X coordinate of the center of the box
35 - fY: Y coordinate of the center of the box
36 - fZ: Z coordinate of the center of the box
37 - fDx: half length in X
38 - fDy: half length in Y
39 - fDz: half length in Z
40 - fTheta: Angle of box z axis with respect to main Z axis
41 - fPhi: Angle of box x axis with respect to main Xaxis
42 - fRefObject: A reference to an object
43*/
44
45////////////////////////////////////////////////////////////////////////////////
46/// Marker3DBox default constructor
47
49{
50 fRefObject = nullptr;
51 fDx = 1;
52 fDy = 1;
53 fDz = 1;
54 fX = 0;
55 fY = 0;
56 fZ = 0;
57 fTheta = 0;
58 fPhi = 0;
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Marker3DBox normal constructor
64
66 Float_t dx, Float_t dy, Float_t dz,
67 Float_t theta, Float_t phi)
68 :TAttLine(1,1,1), TAttFill(1,0)
69{
70 fDx = dx;
71 fDy = dy;
72 fDz = dz;
73 fX = x;
74 fY = y;
75 fZ = z;
76 fTheta = theta;
77 fPhi = phi;
78 fRefObject = nullptr;
80}
81
82////////////////////////////////////////////////////////////////////////////////
83/// copy constructor
84
86 TObject(m3d),
87 TAttLine(m3d),
88 TAttFill(m3d),
89 TAtt3D(m3d),
90 fX(m3d.fX),
91 fY(m3d.fY),
92 fZ(m3d.fZ),
93 fDx(m3d.fDx),
94 fDy(m3d.fDy),
95 fDz(m3d.fDz),
96 fTheta(m3d.fTheta),
97 fPhi(m3d.fPhi),
98 fRefObject(m3d.fRefObject)
99{
100}
101
102////////////////////////////////////////////////////////////////////////////////
103/// assignment operator
104
106{
107 if(this!=&m3d) {
109 TAttLine::operator=(m3d);
110 TAttFill::operator=(m3d);
111 TAtt3D::operator=(m3d);
112 fX=m3d.fX;
113 fY=m3d.fY;
114 fZ=m3d.fZ;
115 fDx=m3d.fDx;
116 fDy=m3d.fDy;
117 fDz=m3d.fDz;
118 fTheta=m3d.fTheta;
119 fPhi=m3d.fPhi;
121 }
122 return *this;
123}
124
125////////////////////////////////////////////////////////////////////////////////
126/// Marker3DBox shape default destructor
127
129{
130}
131
132
133////////////////////////////////////////////////////////////////////////////////
134/// Compute distance from point px,py to a Marker3DBox
135///
136/// Compute the closest distance of approach from point px,py to each corner
137/// point of the Marker3DBox.
138
140{
141 const Int_t numPoints = 8;
142 Int_t dist = 9999;
143 Double_t points[3*numPoints];
144
145 TView *view = gPad->GetView();
146 if (!view) return dist;
147 const Int_t seg1[12] = {0,1,2,3,4,5,6,7,0,1,2,3};
148 const Int_t seg2[12] = {1,2,3,0,5,6,7,4,4,5,6,7};
149
151
152 Int_t i, i1, i2, dsegment;
154 Double_t xndc[3];
155 for (i = 0; i < 12; i++) {
156 i1 = 3*seg1[i];
157 view->WCtoNDC(&points[i1], xndc);
158 x1 = xndc[0];
159 y1 = xndc[1];
160
161 i2 = 3*seg2[i];
162 view->WCtoNDC(&points[i2], xndc);
163 x2 = xndc[0];
164 y2 = xndc[1];
165 dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
166 if (dsegment < dist) dist = dsegment;
167 }
168 if (dist < 5) {
169 gPad->SetCursor(kCross);
170 if (fRefObject) {gPad->SetSelected(fRefObject); return 0;}
171 }
172 return dist;
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Execute action corresponding to one event
177///
178/// This member function must be implemented to realize the action
179/// corresponding to the mouse click on the object in the window
180
182{
183 if (!gPad) return;
184 if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
185}
186
187////////////////////////////////////////////////////////////////////////////////
188/// Paint marker 3D box.
189
190void TMarker3DBox::Paint(Option_t * /* option */ )
191{
192 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
193
194 buffer.ClearSectionsValid();
195
196 // Section kCore
197
198 // If we are just a temporary object then no 'real object' to
199 // pass to viewer
200 if (TestBit(kTemporary)) {
201 buffer.fID = nullptr;
202 } else {
203 buffer.fID = this;
204 }
205 buffer.fColor = GetLineColor();
206 buffer.fTransparency = 0;
207 buffer.fLocalFrame = kFALSE;
209
210 // We fill kCore and kRawSizes on first pass and try with viewer
211 TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
212 if (!viewer3D) return;
213 Int_t reqSections = viewer3D->AddObject(buffer);
214 if (reqSections == TBuffer3D::kNone) {
215 return;
216 }
217
218 if (reqSections & TBuffer3D::kRawSizes) {
219 Int_t nbPnts = 8;
220 Int_t nbSegs = 12;
221 Int_t nbPols = 6;
222 if (!buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*nbPols)) {
223 return;
224 }
226 }
227
228 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
229 // Points
230 SetPoints(buffer.fPnts);
231
232 // Transform points
233 if (gGeometry && !buffer.fLocalFrame) {
234 Double_t dlocal[3];
235 Double_t dmaster[3];
236 for (UInt_t j=0; j<buffer.NbPnts(); j++) {
237 dlocal[0] = buffer.fPnts[3*j];
238 dlocal[1] = buffer.fPnts[3*j+1];
239 dlocal[2] = buffer.fPnts[3*j+2];
240 gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
241 buffer.fPnts[3*j] = dmaster[0];
242 buffer.fPnts[3*j+1] = dmaster[1];
243 buffer.fPnts[3*j+2] = dmaster[2];
244 }
245 }
246
247 // Basic colors: 0, 1, ... 8
248 Int_t c = (((GetLineColor()) %8) -1) * 4;
249 if (c < 0) c = 0;
250
251 // Segments
252 buffer.fSegs[ 0] = c ; buffer.fSegs[ 1] = 0 ; buffer.fSegs[ 2] = 1;
253 buffer.fSegs[ 3] = c+1 ; buffer.fSegs[ 4] = 1 ; buffer.fSegs[ 5] = 2;
254 buffer.fSegs[ 6] = c+1 ; buffer.fSegs[ 7] = 2 ; buffer.fSegs[ 8] = 3;
255 buffer.fSegs[ 9] = c ; buffer.fSegs[10] = 3 ; buffer.fSegs[11] = 0;
256 buffer.fSegs[12] = c+2 ; buffer.fSegs[13] = 4 ; buffer.fSegs[14] = 5;
257 buffer.fSegs[15] = c+2 ; buffer.fSegs[16] = 5 ; buffer.fSegs[17] = 6;
258 buffer.fSegs[18] = c+3 ; buffer.fSegs[19] = 6 ; buffer.fSegs[20] = 7;
259 buffer.fSegs[21] = c+3 ; buffer.fSegs[22] = 7 ; buffer.fSegs[23] = 4;
260 buffer.fSegs[24] = c ; buffer.fSegs[25] = 0 ; buffer.fSegs[26] = 4;
261 buffer.fSegs[27] = c+2 ; buffer.fSegs[28] = 1 ; buffer.fSegs[29] = 5;
262 buffer.fSegs[30] = c+1 ; buffer.fSegs[31] = 2 ; buffer.fSegs[32] = 6;
263 buffer.fSegs[33] = c+3 ; buffer.fSegs[34] = 3 ; buffer.fSegs[35] = 7;
264
265 // Polygons
266 buffer.fPols[ 0] = c ; buffer.fPols[ 1] = 4 ; buffer.fPols[ 2] = 0;
267 buffer.fPols[ 3] = 9 ; buffer.fPols[ 4] = 4 ; buffer.fPols[ 5] = 8;
268 buffer.fPols[ 6] = c+1 ; buffer.fPols[ 7] = 4 ; buffer.fPols[ 8] = 1;
269 buffer.fPols[ 9] = 10 ; buffer.fPols[10] = 5 ; buffer.fPols[11] = 9;
270 buffer.fPols[12] = c ; buffer.fPols[13] = 4 ; buffer.fPols[14] = 2;
271 buffer.fPols[15] = 11 ; buffer.fPols[16] = 6 ; buffer.fPols[17] = 10;
272 buffer.fPols[18] = c+1 ; buffer.fPols[19] = 4 ; buffer.fPols[20] = 3;
273 buffer.fPols[21] = 8 ; buffer.fPols[22] = 7 ; buffer.fPols[23] = 11;
274 buffer.fPols[24] = c+2 ; buffer.fPols[25] = 4 ; buffer.fPols[26] = 0;
275 buffer.fPols[27] = 3 ; buffer.fPols[28] = 2 ; buffer.fPols[29] = 1;
276 buffer.fPols[30] = c+3 ; buffer.fPols[31] = 4 ; buffer.fPols[32] = 4;
277 buffer.fPols[33] = 5 ; buffer.fPols[34] = 6 ; buffer.fPols[35] = 7;
278
280
283 }
284
285 viewer3D->AddObject(buffer);
286}
287
288////////////////////////////////////////////////////////////////////////////////
289/// Paint 3-d histogram h with marker3dboxes
290
292{
293 Int_t bin,ix,iy,iz;
294 Double_t xmin,xmax,ymin,ymax,zmin,zmax,wmin,wmax,w;
295 TAxis *xaxis = h->GetXaxis();
296 TAxis *yaxis = h->GetYaxis();
297 TAxis *zaxis = h->GetZaxis();
298
299 wmin = h->GetMinimum();
300 wmax = h->GetMaximum();
301
302 //Create or modify 3-d view object
303 TView *view = gPad->GetView();
304 if (!view) {
305 gPad->Range(-1,-1,1,1);
306 view = TView::CreateView(1,nullptr,nullptr);
307 if (!view) return;
308 }
309 view->SetRange(xaxis->GetBinLowEdge(xaxis->GetFirst()),
310 yaxis->GetBinLowEdge(yaxis->GetFirst()),
311 zaxis->GetBinLowEdge(zaxis->GetFirst()),
312 xaxis->GetBinUpEdge(xaxis->GetLast()),
313 yaxis->GetBinUpEdge(yaxis->GetLast()),
314 zaxis->GetBinUpEdge(zaxis->GetLast()));
315
316 view->PadRange(gPad->GetFrameFillColor());
317
318 //Draw TMarker3DBox with size proportional to cell content
319 TMarker3DBox m3;
321 m3.SetRefObject(h);
322 m3.SetDirection(0,0);
323 m3.SetLineColor(h->GetMarkerColor());
324 Double_t scale;
325 for (ix=xaxis->GetFirst();ix<=xaxis->GetLast();ix++) {
326 xmin = h->GetXaxis()->GetBinLowEdge(ix);
327 xmax = xmin + h->GetXaxis()->GetBinWidth(ix);
328 for (iy=yaxis->GetFirst();iy<=yaxis->GetLast();iy++) {
329 ymin = h->GetYaxis()->GetBinLowEdge(iy);
330 ymax = ymin + h->GetYaxis()->GetBinWidth(iy);
331 for (iz=zaxis->GetFirst();iz<=zaxis->GetLast();iz++) {
332 zmin = h->GetZaxis()->GetBinLowEdge(iz);
333 zmax = zmin + h->GetZaxis()->GetBinWidth(iz);
334 bin = h->GetBin(ix,iy,iz);
335 w = h->GetBinContent(bin);
336 if (w < wmin) continue;
337 if (w > wmax) w = wmax;
338 scale = (TMath::Power((w-wmin)/(wmax-wmin),1./3.))/2.;
339 if (scale == 0) continue;
340 m3.SetPosition(0.5*(xmin+xmax),0.5*(ymin+ymax),0.5*(zmin+zmax));
341 m3.SetSize(scale*(xmax-xmin),scale*(ymax-ymin),scale*(zmax-zmin));
342 m3.Paint(option);
343 }
344 }
345 }
346}
347
348////////////////////////////////////////////////////////////////////////////////
349/// Save primitive as a C++ statement(s) on output stream out
350
351void TMarker3DBox::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
352{
353 out<<" "<<std::endl;
354 if (gROOT->ClassSaved(TMarker3DBox::Class())) {
355 out<<" ";
356 } else {
357 out<<" TMarker3DBox *";
358 }
359 out<<"marker3DBox = new TMarker3DBox("<<fX<<","
360 <<fY<<","
361 <<fZ<<","
362 <<fDx<<","
363 <<fDy<<","
364 <<fDz<<","
365 <<fTheta<<","
366 <<fPhi<<");"<<std::endl;
367
368 SaveLineAttributes(out,"marker3DBox",1,1,1);
369 SaveFillAttributes(out,"marker3DBox",1,0);
370
371 out<<" marker3DBox->Draw();"<<std::endl;
372}
373
374////////////////////////////////////////////////////////////////////////////////
375/// Set direction.
376
378{
379 fTheta = theta;
380 fPhi = phi;
381}
382
383////////////////////////////////////////////////////////////////////////////////
384/// Set size.
385
387{
388 fDx = dx;
389 fDy = dy;
390 fDz = dz;
391}
392
393////////////////////////////////////////////////////////////////////////////////
394/// Set position.
395
397{
398 fX = x;
399 fY = y;
400 fZ = z;
401}
402
403////////////////////////////////////////////////////////////////////////////////
404/// Set points.
405
407{
408 if (points) {
409 points[ 0] = -fDx ; points[ 1] = -fDy ; points[ 2] = -fDz;
410 points[ 3] = -fDx ; points[ 4] = fDy ; points[ 5] = -fDz;
411 points[ 6] = fDx ; points[ 7] = fDy ; points[ 8] = -fDz;
412 points[ 9] = fDx ; points[10] = -fDy ; points[11] = -fDz;
413 points[12] = -fDx ; points[13] = -fDy ; points[14] = fDz;
414 points[15] = -fDx ; points[16] = fDy ; points[17] = fDz;
415 points[18] = fDx ; points[19] = fDy ; points[20] = fDz;
416 points[21] = fDx ; points[22] = -fDy ; points[23] = fDz;
417
418 Double_t x, y, z;
419 const Double_t kPI = TMath::Pi();
420 Double_t theta = fTheta*kPI/180;
421 Double_t phi = fPhi*kPI/180;
422 Double_t sinth = TMath::Sin(theta);
423 Double_t costh = TMath::Cos(theta);
424 Double_t sinfi = TMath::Sin(phi);
425 Double_t cosfi = TMath::Cos(phi);
426
427 // Matrix to convert from fruit frame to master frame
428 Double_t m[9];
429 m[0] = costh * cosfi; m[1] = -sinfi; m[2] = sinth*cosfi;
430 m[3] = costh * sinfi; m[4] = cosfi; m[5] = sinth*sinfi;
431 m[6] = -sinth; m[7] = 0; m[8] = costh;
432 for (Int_t i = 0; i < 8; i++) {
433 x = points[3*i];
434 y = points[3*i+1];
435 z = points[3*i+2];
436
437 points[3*i] = fX + m[0] * x + m[1] * y + m[2] * z;
438 points[3*i+1] = fY + m[3] * x + m[4] * y + m[5] * z;
439 points[3*i+2] = fZ + m[6] * x + m[7] * y + m[8] * z;
440 }
441 }
442}
443
444////////////////////////////////////////////////////////////////////////////////
445/// Stream an object of class TMarker3DBox.
446
448{
449 if (R__b.IsReading()) {
450 UInt_t R__s, R__c;
451 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
452 if (R__v > 1) {
453 R__b.ReadClassBuffer(TMarker3DBox::Class(), this, R__v, R__s, R__c);
454 return;
455 }
456 //====process old versions before automatic schema evolution
457 TObject::Streamer(R__b);
458 TAttLine::Streamer(R__b);
459 TAttFill::Streamer(R__b);
460 TAtt3D::Streamer(R__b);
461 R__b >> fX;
462 R__b >> fY;
463 R__b >> fZ;
464 R__b >> fDx;
465 R__b >> fDy;
466 R__b >> fDz;
467 R__b >> fTheta;
468 R__b >> fPhi;
469 R__b >> fRefObject;
470 R__b.CheckByteCount(R__s, R__c, TMarker3DBox::IsA());
471 //====end of old versions
472
473 } else {
475 }
476}
@ kCross
Definition GuiTypes.h:374
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
short Version_t
Definition RtypesCore.h:65
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
constexpr Double_t kPI
Definition TEllipse.cxx:24
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t wmin
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t wmax
Option_t Option_t TPoint TPoint const char y1
R__EXTERN TGeometry * gGeometry
Definition TGeometry.h:158
float xmin
float ymin
float xmax
float ymax
#define gROOT
Definition TROOT.h:406
#define gPad
Use this attribute class when an object should have 3D capabilities.
Definition TAtt3D.h:19
virtual void Streamer(TBuffer &)
Fill Area Attributes class.
Definition TAttFill.h:19
virtual void Streamer(TBuffer &)
virtual void Modify()
Change current fill area attributes if necessary.
Definition TAttFill.cxx:216
virtual void SaveFillAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
Save fill attributes as C++ statement(s) on output stream out.
Definition TAttFill.cxx:239
Line Attributes class.
Definition TAttLine.h:18
virtual void Streamer(TBuffer &)
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void Modify()
Change current line attributes if necessary.
Definition TAttLine.cxx:247
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition TAttLine.cxx:211
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition TAttLine.cxx:275
Class to manage histogram axis.
Definition TAxis.h:32
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition TAxis.cxx:513
Int_t GetLast() const
Return last bin on the axis i.e.
Definition TAxis.cxx:464
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:523
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition TAxis.cxx:453
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:115
UInt_t NbPnts() const
Definition TBuffer3D.h:80
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
void ClearSectionsValid()
Clear any sections marked valid.
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:114
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Int_t fColor
Definition TBuffer3D.h:88
Short_t fTransparency
Definition TBuffer3D.h:89
TObject * fID
Definition TBuffer3D.h:87
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:113
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Version_t ReadVersion(UInt_t *start=nullptr, UInt_t *bcnt=nullptr, const TClass *cl=nullptr)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual void Local2Master(Double_t *local, Double_t *master)
Convert one point from local system to master reference system.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
A special 3-D marker designed for event display.
virtual void SetPoints(Double_t *buff) const
Set points.
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
virtual void SetSize(Float_t dx, Float_t dy, Float_t dz)
Set size.
TMarker3DBox & operator=(const TMarker3DBox &)
assignment operator
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a Marker3DBox.
static void PaintH3(TH1 *h, Option_t *option)
Paint 3-d histogram h with marker3dboxes.
~TMarker3DBox() override
Marker3DBox shape default destructor.
TClass * IsA() const override
void Paint(Option_t *option) override
Paint marker 3D box.
static TClass * Class()
virtual void SetDirection(Float_t theta, Float_t phi)
Set direction.
virtual void SetPosition(Float_t x, Float_t y, Float_t z)
Set position.
TMarker3DBox()
Marker3DBox default constructor.
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
virtual void SetRefObject(TObject *obj=nullptr)
void Streamer(TBuffer &) override
Stream an object of class TMarker3DBox.
Float_t fTheta
TObject * fRefObject
Mother of all ROOT objects.
Definition TObject.h:41
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition TObject.h:296
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
virtual void Streamer(TBuffer &)
Stream an object of class TObject.
Definition TObject.cxx:906
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:798
See TView3D.
Definition TView.h:25
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
static TView * CreateView(Int_t system=1, const Double_t *rmin=nullptr, const Double_t *rmax=nullptr)
Create a concrete default 3-d view via the plug-in manager.
Definition TView.cxx:27
virtual void PadRange(Int_t rback)=0
virtual void SetRange(const Double_t *min, const Double_t *max)=0
Abstract 3D shapes viewer.
virtual Int_t AddObject(const TBuffer3D &buffer, Bool_t *addChildren=nullptr)=0
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:725
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:598
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:592
TMarker m
Definition textangle.C:8