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