Logo ROOT   6.16/01
Reference Guide
TView3D.cxx
Go to the documentation of this file.
1// @(#)root/g3d:$Id$
2// Author: Rene Brun, Nenad Buncic, Evgueni Tcherniaev, Olivier Couet 18/08/95
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 "RConfigure.h"
13
14#include "TVirtualPad.h"
15#include "TView3D.h"
16#include "TAxis3D.h"
17#include "TPolyLine3D.h"
18#include "TVirtualX.h"
19#include "TROOT.h"
20#include "TClass.h"
21#include "TList.h"
22#include "TPluginManager.h"
23#include "TMath.h"
24
25// Remove when TView3Der3DPad fix in ExecuteRotateView() is removed
26#include "TVirtualViewer3D.h"
27
29
30//const Int_t kPerspective = BIT(14);
31
32const Int_t kCARTESIAN = 1;
33const Int_t kPOLAR = 2;
34const Double_t kRad = 3.14159265358979323846/180.0;
35
36/** \class TView3D
37\ingroup g3d
38The 3D view class.
39
40This package was originally written by Evgueni Tcherniaev from IHEP/Protvino.
41
42The original Fortran implementation was adapted to HIGZ/PAW by Olivier Couet and
43Evgueni Tcherniaev.
44
45This View class is a subset of the original system. It has been converted to a
46C++ class by Rene Brun.
47
48TView3D creates a 3-D view in the current pad. In this 3D view Lego and Surface
49plots can be drawn and also 3D polyline and markers. Most of the time a TView3D
50is created automatically when a 3D object needs to be painted in a pad (for
51instance a Lego or a Surface plot).
52
53In some case a TView3D should be explicitly. For instance to paint a 3D simple
54scene composed of simple objects like polylines and polymarkers.
55The following macro gives an example:
56
57Begin_Macro(source)
58{
59 cV3D = new TCanvas("cV3D","PolyLine3D & PolyMarker3D Window",200,10,500,500);
60
61 // Creating a view
62 TView3D *view = (TView3D*) TView::CreateView(1);
63 view->SetRange(5,5,5,25,25,25);
64
65 // Create a first PolyLine3D
66 TPolyLine3D *pl3d1 = new TPolyLine3D(6);
67 pl3d1->SetPoint(0, 10, 20, 10);
68 pl3d1->SetPoint(1, 15, 15, 15);
69 pl3d1->SetPoint(2, 20, 20, 20);
70 pl3d1->SetPoint(3, 20, 10, 20);
71 pl3d1->SetPoint(4, 10, 10, 20);
72 pl3d1->SetPoint(5, 10, 10, 10);
73
74 // Create a first PolyMarker3D
75 TPolyMarker3D *pm3d1 = new TPolyMarker3D(9);
76 pm3d1->SetPoint( 0, 10, 10, 10);
77 pm3d1->SetPoint( 1, 20, 20, 20);
78 pm3d1->SetPoint( 2, 10, 20, 20);
79 pm3d1->SetPoint( 3, 10, 10, 20);
80 pm3d1->SetPoint( 4, 20, 20, 10);
81 pm3d1->SetPoint( 5, 20, 10, 10);
82 pm3d1->SetPoint( 6, 20, 10, 20);
83 pm3d1->SetPoint( 7, 10, 20, 10);
84 pm3d1->SetPoint( 8, 15, 15, 15);
85 pm3d1->SetMarkerSize(2);
86 pm3d1->SetMarkerColor(4);
87 pm3d1->SetMarkerStyle(2);
88
89 // Draw
90 pl3d1->Draw();
91 pm3d1->Draw();
92}
93End_Macro
94
95
96Several coordinate systems are available:
97
98 - Cartesian
99 - Polar
100 - Cylindrical
101 - Spherical
102 - PseudoRapidity/Phi
103*/
104
105////////////////////////////////////////////////////////////////////////////////
106/// Default constructor
107
109{
110 fSystem = 0;
111 fOutline = 0;
115
116 fPsi = 0;
117 Int_t i;
118 for (i = 0; i < 3; i++) {
119 fRmin[i] = 0;
120 fRmax[i] = 1;
121 fX1[i] = fX2[i] = fY1[i] = fY2[i] = fZ1[i] = fZ2[i] = 0;
122 }
123
124 if (gPad) {
125 fLongitude = -90 - gPad->GetPhi();
126 fLatitude = 90 - gPad->GetTheta();
127 } else {
128 fLongitude = 0;
129 fLatitude = 0;
130 }
131 Int_t irep = 1;
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// TView3D constructor
137///
138/// Creates a 3-D view in the current pad
139/// rmin[3], rmax[3] are the limits of the object depending on
140/// the selected coordinate system
141///
142/// Before drawing a 3-D object in a pad, a 3-D view must be created.
143/// Note that a view is automatically created when drawing legos or surfaces.
144///
145/// The coordinate system is selected via system:
146/// - system = 1 Cartesian
147/// - system = 2 Polar
148/// - system = 3 Cylindrical
149/// - system = 4 Spherical
150/// - system = 5 PseudoRapidity/Phi
151
152TView3D::TView3D(Int_t system, const Double_t *rmin, const Double_t *rmax) : TView()
153{
154 Int_t irep;
155
157
158 fSystem = system;
159 fOutline = 0;
163
164 if (system == kCARTESIAN || system == kPOLAR || system == 11) fPsi = 0;
165 else fPsi = 90;
166
167 // By default pad range in 3-D view is (-1,-1,1,1), so ...
168 if (gPad) gPad->Range(-1, -1, 1, 1);
170
171 Int_t i;
172 for (i = 0; i < 3; i++) {
173 if (rmin) fRmin[i] = rmin[i];
174 else fRmin[i] = 0;
175 if (rmax) fRmax[i] = rmax[i];
176 else fRmax[i] = 1;
177 fX1[i] = fX2[i] = fY1[i] = fY2[i] = fZ1[i] = fZ2[i] = 0;
178 }
179
180 if (gPad) {
181 fLongitude = -90 - gPad->GetPhi();
182 fLatitude = 90 - gPad->GetTheta();
183 } else {
184 fLongitude = 0;
185 fLatitude = 0;
186 }
188
189 if (gPad) gPad->SetView(this);
190 if (system == 11) SetPerspective();
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// Copy constructor.
195
197 :TView(tv),
198 fLatitude(tv.fLatitude),
199 fLongitude(tv.fLongitude),
200 fPsi(tv.fPsi),
201 fDview(tv.fDview),
202 fDproj(tv.fDproj),
203 fUpix(tv.fUpix),
204 fVpix(tv.fVpix),
205 fSystem(tv.fSystem),
206 fOutline(tv.fOutline),
207 fDefaultOutline(tv.fDefaultOutline),
208 fAutoRange(tv.fAutoRange),
209 fChanged(tv.fChanged)
210{
211 for (Int_t i=0; i<16; i++) {
212 fTN[i]=tv.fTN[i];
213 fTB[i]=tv.fTB[i];
214 fTnorm[i]=tv.fTnorm[i];
215 fTback[i]=tv.fTback[i];
216 }
217 for(Int_t i=0; i<3; i++) {
218 fRmax[i]=tv.fRmax[i];
219 fRmin[i]=tv.fRmin[i];
220 fX1[i]=tv.fX1[i];
221 fX2[i]=tv.fX2[i];
222 fY1[i]=tv.fY1[i];
223 fY2[i]=tv.fY2[i];
224 fZ1[i]=tv.fZ1[i];
225 fZ2[i]=tv.fZ2[i];
226 }
227 for(Int_t i=0; i<4; i++)
228 fUVcoord[i]=tv.fUVcoord[i];
229}
230
231////////////////////////////////////////////////////////////////////////////////
232/// Assignment operator.
233
235{
236 if (this!=&tv) {
240 fPsi=tv.fPsi;
241 fDview=tv.fDview;
242 fDproj=tv.fDproj;
243 fUpix=tv.fUpix;
244 fVpix=tv.fVpix;
245 fSystem=tv.fSystem;
250 for(Int_t i=0; i<16; i++) {
251 fTN[i]=tv.fTN[i];
252 fTB[i]=tv.fTB[i];
253 fTnorm[i]=tv.fTnorm[i];
254 fTback[i]=tv.fTback[i];
255 }
256 for(Int_t i=0; i<3; i++) {
257 fRmax[i]=tv.fRmax[i];
258 fRmin[i]=tv.fRmin[i];
259 fX1[i]=tv.fX1[i];
260 fX2[i]=tv.fX2[i];
261 fY1[i]=tv.fY1[i];
262 fY2[i]=tv.fY2[i];
263 fZ1[i]=tv.fZ1[i];
264 fZ2[i]=tv.fZ2[i];
265 }
266 for(Int_t i=0; i<4; i++)
267 fUVcoord[i]=tv.fUVcoord[i];
268 }
269 return *this;
270}
271
272////////////////////////////////////////////////////////////////////////////////
273/// TView3D default destructor.
274
276{
277 if (fOutline) fOutline->Delete();
278 delete fOutline;
279 fOutline = 0;
280}
281
282////////////////////////////////////////////////////////////////////////////////
283/// Define axis vertices.
284///
285/// Input:
286/// - ANG - angle between X and Y axis (not used anymore)
287///
288/// Output:
289/// - AV(3,8) - axis vertices
290/// - IX1 - 1st point of X-axis (x-min)
291/// - IX2 - 2nd point of X-axis (x-max)
292/// - IY1 - 1st point of Y-axis (y-min)
293/// - IY2 - 2nd point of Y-axis (y-max)
294/// - IZ1 - 1st point of Z-axis (z-min)
295/// - IZ2 - 2nd point of Z-axis (z-max)
296
297void TView3D::AxisVertex(Double_t, Double_t *av, Int_t &ix1, Int_t &ix2, Int_t &iy1, Int_t &iy2, Int_t &iz1, Int_t &iz2)
298{
299 Double_t p[8][3] = {
300 { fRmin[0], fRmin[1], fRmin[2] },
301 { fRmax[0], fRmin[1], fRmin[2] },
302 { fRmax[0], fRmax[1], fRmin[2] },
303 { fRmin[0], fRmax[1], fRmin[2] },
304 { fRmin[0], fRmin[1], fRmax[2] },
305 { fRmax[0], fRmin[1], fRmax[2] },
306 { fRmax[0], fRmax[1], fRmax[2] },
307 { fRmin[0], fRmax[1], fRmax[2] }
308 };
309 Int_t inodes[4][8] = {
310 { 2,3,4,1, 6,7,8,5 }, // x+, y+
311 { 3,4,1,2, 7,8,5,6 }, // x-, y+
312 { 1,2,3,4, 5,6,7,8 }, // x+, y-
313 { 4,1,2,3, 8,5,6,7 } // x-, y-
314 };
315 Int_t ixyminmax[16][4] = { // 8
316 { 3,2, 1,2 }, // x+, y+, z+, z-up 5 / \ 7
317 { 2,1, 3,2 }, // x-, y+, z+, z-up |\6/|
318 { 1,2, 2,3 }, // x+, y-, z+, z-up | | | Top view
319 { 2,3, 2,1 }, // x-, y-, z+, z-up 1 \|/ 3
320 // 2 6
321 { 4,1, 4,3 }, // x+, y+, z-, z-up 5 /|\ 7
322 { 3,4, 4,1 }, // x-, y+, z-, z-up | | |
323 { 4,3, 1,4 }, // x+, y-, z-, z-up Bottom view |/2\|
324 { 1,4, 3,4 }, // x-, y-, z-, z-up 1 \ / 3
325 // 2 4
326 { 8,5, 8,7 }, // x+, y+, z+, z-down 3 /|\ 1
327 { 7,8, 8,5 }, // x-, y+, z+, z-down | | |
328 { 8,7, 5,8 }, // x+, y-, z+, z-down |/6\| Bottom view
329 { 5,8, 7,8 }, // x-, y-, z+, z-down 7 \ / 5
330 // 8 4
331 { 7,6, 5,6 }, // x+, y+, z-, z-down 3 / \ 1
332 { 6,5, 7,6 }, // x-, y+, z-, z-down |\2/|
333 { 5,6, 6,7 }, // x+, y-, z-, z-down Top view | | |
334 { 6,7, 6,5 } // x-, y-, z-, z-down 7 \|/ 5
335 }; // 6
336
337 // Set vertices
338 Int_t icase = 0;
339 if (fTnorm[ 8] <= 0) icase += 1; // z projection of (1,0,0)
340 if (fTnorm[ 9] <= 0) icase += 2; // z projection of (0,1,0)
341 for (Int_t i=0; i<8; ++i) {
342 Int_t k = inodes[icase][i] - 1;
343 av[i*3+0] = p[k][0];
344 av[i*3+1] = p[k][1];
345 av[i*3+2] = p[k][2];
346 }
347
348 // Set indices for min and max
349 if (fTnorm[10] < 0) icase += 4; // z projection of (0,0,1)
350 if (fTnorm[ 6] < 0) icase += 8; // y projection of (0,0,1)
351 ix1 = ixyminmax[icase][0];
352 ix2 = ixyminmax[icase][1];
353 iy1 = ixyminmax[icase][2];
354 iy2 = ixyminmax[icase][3];
355 iz1 = (icase < 8) ? 1 : 3;
356 iz2 = (icase < 8) ? 5 : 7;
357}
358
359////////////////////////////////////////////////////////////////////////////////
360/// Define perspective view.
361///
362/// Compute transformation matrix from world coordinates
363/// to normalized coordinates (-1 to +1)
364///
365/// Input :
366/// - theta, phi - spherical angles giving the direction of projection
367/// - psi - screen rotation angle
368/// - cov[3] - center of view
369/// - dview - distance from COV to COP (center of projection)
370/// - umin, umax, vmin, vmax - view window in projection plane
371/// - dproj - distance from COP to projection plane
372/// - bcut, fcut - backward/forward range w.r.t projection plane (fcut<=0)
373///
374/// Output :
375/// - nper[16] - normalizing transformation
376///
377/// compute tr+rot to get COV in origin, view vector parallel to -Z axis, up
378/// vector parallel to Y.
379///
380/// ~~~ {.cpp}
381/// ^Yv UP ^ proj. plane
382/// | | /|
383/// | | / |
384/// | dproj / x--- center of window (COW)
385/// COV |----------|--x--|------------> Zv
386/// / | VRP'z
387/// / ---> | /
388/// / VPN |/
389/// Xv
390/// ~~~
391///
392/// 1. translate COP to origin of MARS : Tper = T(-copx, -copy, -copz)
393/// 2. rotate VPN : R = Rz(-psi)*Rx(-theta)*Rz(-phi) (inverse Euler)
394/// 3. left-handed screen reference to right-handed one of MARS : Trl
395///
396/// T12 = Tper*R*Trl
397
398
400{
401 Double_t t12[16];
402 Double_t cov[3];
403 Int_t i;
404 for (i=0; i<3; i++) cov[i] = 0.5*(fRmax[i]+fRmin[i]);
405
412
413 t12[0] = c1*c3 - s1*c2*s3;
414 t12[4] = c1*s3 + s1*c2*c3;
415 t12[8] = s1*s2;
416 t12[3] = 0;
417
418 t12[1] = -s1*c3 - c1*c2*s3;
419 t12[5] = -s1*s3 + c1*c2*c3;
420 t12[9] = c1*s2;
421 t12[7] = 0;
422
423 t12[2] = s2*s3;
424 t12[6] = -s2*c3;
425 t12[10] = c2; // contains Trl
426 t12[11] = 0;
427
428 // translate with -COP (before rotation):
429 t12[12] = -(cov[0]*t12[0]+cov[1]*t12[4]+cov[2]*t12[8]);
430 t12[13] = -(cov[0]*t12[1]+cov[1]*t12[5]+cov[2]*t12[9]);
431 t12[14] = -(cov[0]*t12[2]+cov[1]*t12[6]+cov[2]*t12[10]);
432 t12[15] = 1;
433
434 // translate with (0, 0, -dview) after rotation
435
436 t12[14] -= fDview;
437
438 // reflection on Z :
439 t12[2] *= -1;
440 t12[6] *= -1;
441 t12[10] *= -1;
442 t12[14] *= -1;
443
444 // Now we shear the center of window from (0.5*(umin+umax), 0.5*(vmin+vmax), dproj)
445 // to (0, 0, dproj)
446
447 Double_t a2 = -fUVcoord[0]/fDproj; // shear coef. on x
448 Double_t b2 = -fUVcoord[1]/fDproj; // shear coef. on y
449
450 // | 1 0 0 0 |
451 // SHz(a2,b2) = | 0 1 0 0 |
452 // | a2 b2 1 0 |
453 // | 0 0 0 1 |
454
455 fTnorm[0] = t12[0] + a2*t12[2];
456 fTnorm[1] = t12[1] + b2*t12[2];
457 fTnorm[2] = t12[2];
458 fTnorm[3] = 0;
459
460 fTnorm[4] = t12[4] + a2*t12[6];
461 fTnorm[5] = t12[5] + b2*t12[6];
462 fTnorm[6] = t12[6];
463 fTnorm[7] = 0;
464
465 fTnorm[8] = t12[8] + a2*t12[10];
466 fTnorm[9] = t12[9] + b2*t12[10];
467 fTnorm[10] = t12[10];
468 fTnorm[11] = 0;
469
470 fTnorm[12] = t12[12] + a2*t12[14];
471 fTnorm[13] = t12[13] + b2*t12[14];
472 fTnorm[14] = t12[14];
473 fTnorm[15] = 1;
474
475 // Scale so that the view volume becomes the canonical one
476 //
477 // Sper = (2/(umax-umin), 2/(vmax-vmin), 1/dproj
478 //
479 Double_t sz = 1./fDproj;
480 Double_t sx = 1./fUVcoord[2];
481 Double_t sy = 1./fUVcoord[3];
482
483 fTnorm[0] *= sx;
484 fTnorm[4] *= sx;
485 fTnorm[8] *= sx;
486 fTnorm[1] *= sy;
487 fTnorm[5] *= sy;
488 fTnorm[9] *= sy;
489 fTnorm[2] *= sz;
490 fTnorm[6] *= sz;
491 fTnorm[10] *= sz;
492 fTnorm[12] *= sx;
493 fTnorm[13] *= sy;
494 fTnorm[14] *= sz;
495}
496
497////////////////////////////////////////////////////////////////////////////////
498/// Define view direction (in spherical coordinates)
499///
500/// Compute transformation matrix from world coordinates
501/// to normalized coordinates (-1 to +1)
502///
503/// Input:
504/// - S(3) - scale factors
505/// - C(3) - centre of scope
506/// - COSPHI - longitude COS
507/// - SINPHI - longitude SIN
508/// - COSTHE - latitude COS (angle between +Z and view direction.)
509/// - SINTHE - latitude SIN
510/// - COSPSI - screen plane rotation angle COS
511/// - SINPSI - screen plane rotation angle SIN
512
514 Double_t cosphi, Double_t sinphi,
515 Double_t costhe, Double_t sinthe,
516 Double_t cospsi, Double_t sinpsi,
517 Double_t *tnorm, Double_t *tback)
518{
519 if (IsPerspective()) {
521 return;
522 }
523 Int_t i, k;
524 Double_t tran[16] /* was [4][4] */, rota[16] /* was [4][4] */;
525 Double_t c1, c2, c3, s1, s2, s3, scalex, scaley, scalez;
526
527 // Parameter adjustments
528 tback -= 5;
529 tnorm -= 5;
530
531 scalex = s[0];
532 scaley = s[1];
533 scalez = s[2];
534
535 //*-*- S E T T R A N S L A T I O N M A T R I X
536 tran[0] = 1 / scalex;
537 tran[1] = 0;
538 tran[2] = 0;
539 tran[3] = -c[0] / scalex;
540
541 tran[4] = 0;
542 tran[5] = 1 / scaley;
543 tran[6] = 0;
544 tran[7] = -c[1] / scaley;
545
546 tran[8] = 0;
547 tran[9] = 0;
548 tran[10] = 1 / scalez;
549 tran[11] = -c[2] / scalez;
550
551 tran[12] = 0;
552 tran[13] = 0;
553 tran[14] = 0;
554 tran[15] = 1;
555
556 //*-*- S E T R O T A T I O N M A T R I X
557 // ( C(PSI) S(PSI) 0) (1 0 0 ) ( C(90+PHI) S(90+PHI) 0)
558 // (-S(PSI) C(PSI) 0) * (0 C(THETA) S(THETA)) * (-S(90+PHI) C(90+PHI) 0)
559 // ( 0 0 1) (0 -S(THETA) C(THETA)) ( 0 0 1)
560 c1 = cospsi;
561 s1 = sinpsi;
562 c2 = costhe;
563 s2 = sinthe;
564 c3 = -sinphi;
565 s3 = cosphi;
566
567 rota[0] = c1*c3 - s1*c2*s3;
568 rota[1] = c1*s3 + s1*c2*c3;
569 rota[2] = s1*s2;
570 rota[3] = 0;
571
572 rota[4] = -s1*c3 - c1* c2*s3;
573 rota[5] = -s1*s3 + c1* c2*c3;
574 rota[6] = c1*s2;
575 rota[7] = 0;
576
577 rota[8] = s2*s3;
578 rota[9] = -s2*c3;
579 rota[10] = c2;
580 rota[11] = 0;
581
582 rota[12] = 0;
583 rota[13] = 0;
584 rota[14] = 0;
585 rota[15] = 1;
586
587 //*-*- F I N D T R A N S F O R M A T I O N M A T R I X
588 for (i = 1; i <= 3; ++i) {
589 for (k = 1; k <= 4; ++k) {
590 tnorm[k + (i << 2)] = rota[(i << 2) - 4]*tran[k - 1] + rota[(i
591 << 2) - 3]*tran[k + 3] + rota[(i << 2) - 2]*tran[k +7]
592 + rota[(i << 2) - 1]*tran[k + 11];
593 }
594 }
595
596 //*-*- S E T B A C K T R A N S L A T I O N M A T R I X
597 tran[0] = scalex;
598 tran[3] = c[0];
599
600 tran[5] = scaley;
601 tran[7] = c[1];
602
603 tran[10] = scalez;
604 tran[11] = c[2];
605
606 //*-*- F I N D B A C K T R A N S F O R M A T I O N
607 for (i = 1; i <= 3; ++i) {
608 for (k = 1; k <= 4; ++k) {
609 tback[k + (i << 2)] = tran[(i << 2) - 4]*rota[(k << 2) - 4] +
610 tran[(i << 2) - 3]*rota[(k << 2) - 3] + tran[(i << 2) -2]
611 *rota[(k << 2) - 2] + tran[(i << 2) - 1]*rota[(k <<2) - 1];
612 }
613 }
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Draw the outline of a cube while rotating a 3-d object in the pad.
618
620{
621 TPolyLine3D::DrawOutlineCube(outline,rmin,rmax);
622}
623
624////////////////////////////////////////////////////////////////////////////////
625/// Execute action corresponding to one event.
626
628{
629 ExecuteRotateView(event,px,py);
630}
631
632////////////////////////////////////////////////////////////////////////////////
633/// Execute action corresponding to one event.
634///
635/// This member function is called when a object is clicked with the locator
636///
637/// If Left button clicked in the object area, while the button is kept down
638/// the cube representing the surrounding frame for the corresponding
639/// new latitude and longitude position is drawn.
640
642{
643 static Int_t system, framewasdrawn;
644 static Double_t xrange, yrange, xmin, ymin, longitude1, latitude1, longitude2, latitude2;
645 static Double_t newlatitude, newlongitude, oldlatitude, oldlongitude;
646 Double_t dlatitude, dlongitude, x, y;
647 Int_t irep = 0;
648 Double_t psideg;
649 Bool_t opaque = gPad->OpaqueMoving();
650
651 // All coordinates transformation are from absolute to relative
652 if (!gPad->IsEditable()) return;
653 gPad->AbsCoordinates(kTRUE);
654
655 switch (event) {
656
657 case kKeyPress :
658 fChanged = kTRUE;
659 MoveViewCommand(Char_t(px), py);
660 break;
661
662 case kMouseMotion:
663 gPad->SetCursor(kRotate);
664 break;
665
666 case kButton1Down:
667
668 // remember position of the cube
669 xmin = gPad->GetX1();
670 ymin = gPad->GetY1();
671 xrange = gPad->GetX2() - xmin;
672 yrange = gPad->GetY2() - ymin;
673 x = gPad->PixeltoX(px);
674 y = gPad->PixeltoY(py);
675 system = GetSystem();
676 framewasdrawn = 0;
677 if (system == kCARTESIAN || system == kPOLAR || system == 11 || IsPerspective()) {
678 longitude1 = 180*(x-xmin)/xrange;
679 latitude1 = 90*(y-ymin)/yrange;
680 } else {
681 latitude1 = 90*(x-xmin)/xrange;
682 longitude1 = 180*(y-ymin)/yrange;
683 }
684 newlongitude = oldlongitude = -90 - gPad->GetPhi();
685 newlatitude = oldlatitude = 90 - gPad->GetTheta();
686 psideg = GetPsi();
687
688 // if outline isn't set, make it look like a cube
689 if(!fOutline)
691 break;
692
693 case kButton1Motion:
694 {
695 // draw the surrounding frame for the current mouse position
696 // first: Erase old frame
697 fChanged = kTRUE;
698 if (framewasdrawn && !opaque) fOutline->Paint();
699 framewasdrawn = 1;
700 x = gPad->PixeltoX(px);
701 y = gPad->PixeltoY(py);
702 if (system == kCARTESIAN || system == kPOLAR || system == 11 || IsPerspective()) {
703 longitude2 = 180*(x-xmin)/xrange;
704 latitude2 = 90*(y-ymin)/yrange;
705 } else {
706 latitude2 = 90*(x-xmin)/xrange;
707 longitude2 = 180*(y-ymin)/yrange;
708 }
709 dlongitude = longitude2 - longitude1;
710 dlatitude = latitude2 - latitude1;
711 newlatitude = oldlatitude + dlatitude;
712 newlongitude = oldlongitude - dlongitude;
713 psideg = GetPsi();
714 ResetView(newlongitude, newlatitude, psideg, irep);
715 if (!opaque) {
716 fOutline->Paint();
717 } else {
718 psideg = GetPsi();
719 SetView(newlongitude, newlatitude, psideg, irep);
720 gPad->SetPhi(-90-newlongitude);
721 gPad->SetTheta(90-newlatitude);
722 gPad->Modified(kTRUE);
723 }
724
725 break;
726 }
727 case kButton1Up:
728 if (gROOT->IsEscaped()) {
729 gROOT->SetEscape(kFALSE);
730 if (opaque) {
731 psideg = GetPsi();
732 SetView(oldlongitude, oldlatitude, psideg, irep);
733 gPad->SetPhi(-90-oldlongitude);
734 gPad->SetTheta(90-oldlatitude);
735 gPad->Modified(kTRUE);
736 }
737 break;
738 }
739
740 // Temporary fix for 2D drawing problems on pad. fOutline contains
741 // a TPolyLine3D object for the rotation box. This will be painted
742 // through a newly created TView3Der3DPad instance, which is left
743 // behind on pad. This remaining creates 2D drawing problems.
744 //
745 // This is a TEMPORARY fix - will be removed when proper multiple viewers
746 // on pad problems are resolved.
747 if (gPad) {
748 TVirtualViewer3D *viewer = gPad->GetViewer3D();
749 if (viewer && !strcmp(viewer->IsA()->GetName(),"TView3Der3DPad")) {
750 gPad->ReleaseViewer3D();
751 delete viewer;
752 }
753 }
754 // End fix
755
756 // Recompute new view matrix and redraw
757 psideg = GetPsi();
758 SetView(newlongitude, newlatitude, psideg, irep);
759 gPad->SetPhi(-90-newlongitude);
760 gPad->SetTheta(90-newlatitude);
761 gPad->Modified(kTRUE);
762
763 // Set line color, style and width
764 gVirtualX->SetLineColor(-1);
765 gVirtualX->SetLineStyle(-1);
766 gVirtualX->SetLineWidth(-1);
767 break;
768 }
769
770 // set back to default transformation mode
771 gPad->AbsCoordinates(kFALSE);
772}
773
774////////////////////////////////////////////////////////////////////////////////
775/// Find Z component of NORMAL in normalized coordinates.
776///
777/// Input:
778/// - X - X-component of NORMAL
779/// - Y - Y-component of NORMAL
780/// - Z - Z-component of NORMAL
781///
782/// Output:
783/// - ZN - Z-component of NORMAL in normalized coordinates
784
786{
787 zn = x*(fTN[1] * fTN[6] - fTN[2] * fTN[5]) + y*(fTN[2] * fTN[4] -
788 fTN[0] * fTN[6]) + z*(fTN[0] * fTN[5] - fTN[1] * fTN[4]);
789}
790
791////////////////////////////////////////////////////////////////////////////////
792/// Find critical PHI sectors.
793///
794/// Input:
795/// - IOPT - options:
796/// 1. from BACK to FRONT 'BF'
797/// 2. from FRONT to BACK 'FB'
798/// - KPHI - number of phi sectors
799/// - APHI(*) - PHI separators (modified internally)
800///
801/// Output:
802/// - IPHI1 - initial sector
803/// - IPHI2 - final sector
804
805void TView3D::FindPhiSectors(Int_t iopt, Int_t &kphi, Double_t *aphi, Int_t &iphi1, Int_t &iphi2)
806{
807 Int_t iphi[2], i, k;
808 Double_t dphi;
809 Double_t x1, x2, z1, z2, phi1, phi2;
810
811 // Parameter adjustments
812 --aphi;
813
814 if (aphi[kphi + 1] == aphi[1]) aphi[kphi + 1] += 360;
815 dphi = TMath::Abs(aphi[kphi + 1] - aphi[1]);
816 if (dphi != 360) {
817 aphi[kphi + 2] = (aphi[1] + aphi[kphi + 1]) / (float)2. + 180;
818 aphi[kphi + 3] = aphi[1] + 360;
819 kphi += 2;
820 }
821
822 //*-*- F I N D C R I T I C A L S E C T O R S
823 k = 0;
824 for (i = 1; i <= kphi; ++i) {
825 phi1 = kRad*aphi[i];
826 phi2 = kRad*aphi[i + 1];
827 x1 = fTN[0]*TMath::Cos(phi1) + fTN[1]*TMath::Sin(phi1);
828 x2 = fTN[0]*TMath::Cos(phi2) + fTN[1]*TMath::Sin(phi2);
829 if (x1 >= 0 && x2 > 0) continue;
830 if (x1 <= 0 && x2 < 0) continue;
831 ++k;
832 if (k == 3) break;
833 iphi[k - 1] = i;
834 }
835 if (k != 2) {
836 Error("FindPhiSectors", "something strange: num. of critical sector not equal 2");
837 iphi1 = 1;
838 iphi2 = 2;
839 return;
840 }
841
842 //*-*- F I N D O R D E R O F C R I T I C A L S E C T O R S
843 phi1 = kRad*(aphi[iphi[0]] + aphi[iphi[0] + 1]) / (float)2.;
844 phi2 = kRad*(aphi[iphi[1]] + aphi[iphi[1] + 1]) / (float)2.;
845 z1 = fTN[8]*TMath::Cos(phi1) + fTN[9]*TMath::Sin(phi1);
846 z2 = fTN[8]*TMath::Cos(phi2) + fTN[9]*TMath::Sin(phi2);
847 if ((z1 <= z2 && iopt == 1) || (z1 > z2 && iopt == 2)) {
848 iphi1 = iphi[0];
849 iphi2 = iphi[1];
850 } else {
851 iphi1 = iphi[1];
852 iphi2 = iphi[0];
853 }
854}
855
856////////////////////////////////////////////////////////////////////////////////
857/// Find critical THETA sectors for given PHI sector.
858///
859/// Input:
860/// - IOPT - options:
861/// 1. from BACK to FRONT 'BF'
862/// 2. from FRONT to BACK 'FB'
863/// - PHI - PHI sector
864/// - KTH - number of THETA sectors
865/// - ATH(*) - THETA separators (modified internally)
866///
867/// Output:
868/// - ITH1 - initial sector
869/// - ITH2 - final sector
870
871void TView3D::FindThetaSectors(Int_t iopt, Double_t phi, Int_t &kth, Double_t *ath, Int_t &ith1, Int_t &ith2)
872{
873 Int_t i, k, ith[2];
874 Double_t z1, z2, cosphi, sinphi, tncons, th1, th2, dth;
875
876 // Parameter adjustments
877 --ath;
878
879 // Function Body
880 dth = TMath::Abs(ath[kth + 1] - ath[1]);
881 if (dth != 360) {
882 ath[kth + 2] = 0.5*(ath[1] + ath[kth + 1]) + 180;
883 ath[kth + 3] = ath[1] + 360;
884 kth += 2;
885 }
886
887 //*-*- F I N D C R I T I C A L S E C T O R S
888 cosphi = TMath::Cos(phi*kRad);
889 sinphi = TMath::Sin(phi*kRad);
890 k = 0;
891 for (i = 1; i <= kth; ++i) {
892 th1 = kRad*ath[i];
893 th2 = kRad*ath[i + 1];
894 FindNormal(TMath::Cos(th1)*cosphi, TMath::Cos(th1)*sinphi, -TMath::Sin(th1), z1);
895 FindNormal(TMath::Cos(th2)*cosphi, TMath::Cos(th2)*sinphi, -TMath::Sin(th2), z2);
896 if (z1 >= 0 && z2 > 0) continue;
897 if (z1 <= 0 && z2 < 0) continue;
898 ++k;
899 if (k == 3) break;
900 ith[k - 1] = i;
901 }
902 if (k != 2) {
903 Error("FindThetaSectors", "Something strange: num. of critical sectors not equal 2");
904 ith1 = 1;
905 ith2 = 2;
906 return;
907 }
908
909 //*-*- F I N D O R D E R O F C R I T I C A L S E C T O R S
910 tncons = fTN[8]*TMath::Cos(phi*kRad) + fTN[9]*TMath::Sin(phi*kRad);
911 th1 = kRad*(ath[ith[0]] + ath[ith[0] + 1]) / (float)2.;
912 th2 = kRad*(ath[ith[1]] + ath[ith[1] + 1]) / (float)2.;
913 z1 = tncons*TMath::Sin(th1) + fTN[10]*TMath::Cos(th1);
914 z2 = tncons*TMath::Sin(th2) + fTN[10]*TMath::Cos(th2);
915 if ((z1 <= z2 && iopt == 1) || (z1 > z2 && iopt == 2)) {
916 ith1 = ith[0];
917 ith2 = ith[1];
918 } else {
919 ith1 = ith[1];
920 ith2 = ith[0];
921 }
922}
923
924////////////////////////////////////////////////////////////////////////////////
925/// Find centre of a MIN-MAX scope and scale factors
926///
927/// Output:
928/// - SCALE(3) - scale factors
929/// - CENTER(3) - centre
930/// - IREP - reply (-1 if error in min-max)
931
932void TView3D::FindScope(Double_t *scale, Double_t *center, Int_t &irep)
933{
934 irep = 0;
935 Double_t sqrt3 = 0.5*TMath::Sqrt(3.0);
936
937 for (Int_t i = 0; i < 3; i++) {
938 if (fRmin[i] >= fRmax[i]) { irep = -1; return;}
939 scale[i] = sqrt3*(fRmax[i] - fRmin[i]);
940 center[i] = 0.5*(fRmax[i] + fRmin[i]);
941 }
942}
943
944////////////////////////////////////////////////////////////////////////////////
945/// Return distance to axis from point px,py.
946///
947/// Algorithm:
948///
949/// ~~~ {.cpp}
950/// A(x1,y1) P B(x2,y2)
951/// ------------------------------------------------
952/// I
953/// I
954/// I
955/// I
956/// M(x,y)
957///
958/// Let us call a = distance AM A=a**2
959/// b = distance BM B=b**2
960/// c = distance AB C=c**2
961/// d = distance PM D=d**2
962/// u = distance AP U=u**2
963/// v = distance BP V=v**2 c = u + v
964///
965/// D = A - U
966/// D = B - V = B -(c-u)**2
967/// ==> u = (A -B +C)/2c
968/// ~~~
969
971{
972 Double_t x1,y1,x2,y2;
973 Double_t x = px;
974 Double_t y = py;
975 ratio = 0;
976
977 if (fSystem != kCARTESIAN) return 9998; // only implemented for Cartesian coordinates
978 if (axis == 1) {
979 x1 = gPad->XtoAbsPixel(fX1[0]);
980 y1 = gPad->YtoAbsPixel(fX1[1]);
981 x2 = gPad->XtoAbsPixel(fX2[0]);
982 y2 = gPad->YtoAbsPixel(fX2[1]);
983 } else if (axis == 2) {
984 x1 = gPad->XtoAbsPixel(fY1[0]);
985 y1 = gPad->YtoAbsPixel(fY1[1]);
986 x2 = gPad->XtoAbsPixel(fY2[0]);
987 y2 = gPad->YtoAbsPixel(fY2[1]);
988 } else {
989 x1 = gPad->XtoAbsPixel(fZ1[0]);
990 y1 = gPad->YtoAbsPixel(fZ1[1]);
991 x2 = gPad->XtoAbsPixel(fZ2[0]);
992 y2 = gPad->YtoAbsPixel(fZ2[1]);
993 }
994 Double_t xx1 = x - x1;
995 Double_t xx2 = x - x2;
996 Double_t x1x2 = x1 - x2;
997 Double_t yy1 = y - y1;
998 Double_t yy2 = y - y2;
999 Double_t y1y2 = y1 - y2;
1000 Double_t a = xx1*xx1 + yy1*yy1;
1001 Double_t b = xx2*xx2 + yy2*yy2;
1002 Double_t c = x1x2*x1x2 + y1y2*y1y2;
1003 if (c <= 0) return 9999;
1005 Double_t u = (a - b + c)/(2*v);
1006 Double_t d = TMath::Abs(a - u*u);
1007
1008 Int_t dist = Int_t(TMath::Sqrt(d) - 0.5);
1009 ratio = u/v;
1010 return dist;
1011}
1012
1013////////////////////////////////////////////////////////////////////////////////
1014/// Get maximum view extent.
1015
1017{
1018 Double_t dx = 0.5*(fRmax[0]-fRmin[0]);
1019 Double_t dy = 0.5*(fRmax[1]-fRmin[1]);
1020 Double_t dz = 0.5*(fRmax[2]-fRmin[2]);
1021 Double_t extent = TMath::Sqrt(dx*dx+dy*dy+dz*dz);
1022 return extent;
1023}
1024
1025////////////////////////////////////////////////////////////////////////////////
1026/// Get Range function.
1027
1029{
1030 for (Int_t i = 0; i < 3; max[i] = fRmax[i], min[i] = fRmin[i], i++) { }
1031}
1032
1033////////////////////////////////////////////////////////////////////////////////
1034/// Get Range function.
1035
1037{
1038 for (Int_t i = 0; i < 3; max[i] = fRmax[i], min[i] = fRmin[i], i++) { }
1039}
1040
1041////////////////////////////////////////////////////////////////////////////////
1042/// Get current window extent.
1043
1045{
1046 u0 = fUVcoord[0];
1047 v0 = fUVcoord[1];
1048 du = fUVcoord[2];
1049 dv = fUVcoord[3];
1050}
1051
1052////////////////////////////////////////////////////////////////////////////////
1053/// Check if point is clipped in perspective view.
1054
1056{
1057 if (TMath::Abs(p[0])>p[2]) return kTRUE;
1058 if (TMath::Abs(p[1])>p[2]) return kTRUE;
1059 return kFALSE;
1060}
1061
1062////////////////////////////////////////////////////////////////////////////////
1063/// Transfer point from normalized to world coordinates.
1064///
1065/// Input:
1066/// - PN(3) - point in world coordinate system
1067/// - PW(3) - point in normalized coordinate system
1068
1070{
1071 Float_t x = pn[0], y = pn[1], z = pn[2];
1072 pw[0] = fTback[0]*x + fTback[1]*y + fTback[2]*z + fTback[3];
1073 pw[1] = fTback[4]*x + fTback[5]*y + fTback[6]*z + fTback[7];
1074 pw[2] = fTback[8]*x + fTback[9]*y + fTback[10]*z + fTback[11];
1075}
1076
1077////////////////////////////////////////////////////////////////////////////////
1078/// Transfer point from normalized to world coordinates.
1079///
1080/// Input:
1081/// - PN(3) - point in world coordinate system
1082/// - PW(3) - point in normalized coordinate system
1083
1085{
1086 Double_t x = pn[0], y = pn[1], z = pn[2];
1087 pw[0] = fTback[0]*x + fTback[1]*y + fTback[2]*z + fTback[3];
1088 pw[1] = fTback[4]*x + fTback[5]*y + fTback[6]*z + fTback[7];
1089 pw[2] = fTback[8]*x + fTback[9]*y + fTback[10]*z + fTback[11];
1090}
1091
1092////////////////////////////////////////////////////////////////////////////////
1093/// Transfer vector of NORMAL from word to normalized coordinates.
1094///
1095/// Input:
1096/// - PW(3) - vector of NORMAL in word coordinate system
1097/// - PN(3) - vector of NORMAL in normalized coordinate system
1098
1100{
1101 Double_t x, y, z, a1, a2, a3, b1, b2, b3, c1, c2, c3;
1102
1103 x = pw[0];
1104 y = pw[1];
1105 z = pw[2];
1106 a1 = fTnorm[0];
1107 a2 = fTnorm[1];
1108 a3 = fTnorm[2];
1109 b1 = fTnorm[4];
1110 b2 = fTnorm[5];
1111 b3 = fTnorm[6];
1112 c1 = fTnorm[8];
1113 c2 = fTnorm[9];
1114 c3 = fTnorm[10];
1115 pn[0] = x*(b2*c3 - b3*c2) + y*(b3*c1 - b1*c3) + z*(b1*c2 - b2*c1);
1116 pn[1] = x*(c2*a3 - c3*a2) + y*(c3*a1 - c1*a3) + z*(c1*a2 - c2*a1);
1117 pn[2] = x*(a2*b3 - a3*b2) + y*(a3*b1 - a1*b3) + z*(a1*b2 - a2*b1);
1118}
1119
1120////////////////////////////////////////////////////////////////////////////////
1121/// Transfer vector of NORMAL from word to normalized coordinates.
1122///
1123/// Input:
1124/// - PW(3) - vector of NORMAL in word coordinate system
1125/// - PN(3) - vector of NORMAL in normalized coordinate system
1126
1128{
1129 Double_t x, y, z, a1, a2, a3, b1, b2, b3, c1, c2, c3;
1130
1131 x = pw[0];
1132 y = pw[1];
1133 z = pw[2];
1134 a1 = fTnorm[0];
1135 a2 = fTnorm[1];
1136 a3 = fTnorm[2];
1137 b1 = fTnorm[4];
1138 b2 = fTnorm[5];
1139 b3 = fTnorm[6];
1140 c1 = fTnorm[8];
1141 c2 = fTnorm[9];
1142 c3 = fTnorm[10];
1143 pn[0] = x*(b2*c3 - b3*c2) + y*(b3*c1 - b1*c3) + z*(b1*c2 - b2*c1);
1144 pn[1] = x*(c2*a3 - c3*a2) + y*(c3*a1 - c1*a3) + z*(c1*a2 - c2*a1);
1145 pn[2] = x*(a2*b3 - a3*b2) + y*(a3*b1 - a1*b3) + z*(a1*b2 - a2*b1);
1146}
1147
1148////////////////////////////////////////////////////////////////////////////////
1149/// Set the correct window size for lego and surface plots.
1150///
1151/// Set the correct window size for lego and surface plots.
1152/// And draw the background if necessary.
1153///
1154/// Input parameters:
1155/// - RBACK : Background colour
1156
1158{
1159 Int_t i, k;
1160 Double_t x, y, z, r1, r2, r3, xx, yy, smax[2];
1161 Double_t xgraf[6], ygraf[6];
1162
1163 for (i = 1; i <= 2; ++i) {
1164 smax[i - 1] = fTnorm[(i << 2) - 1];
1165 for (k = 1; k <= 3; ++k) {
1166 if (fTnorm[k + (i << 2) - 5] < 0) {
1167 smax[i - 1] += fTnorm[k + (i << 2) - 5]*fRmin[k-1];
1168 } else {
1169 smax[i - 1] += fTnorm[k + (i << 2) - 5]*fRmax[k-1];
1170 }
1171 }
1172 }
1173
1174 //*-*- Compute x,y range
1175 Double_t xmin = -smax[0];
1176 Double_t xmax = smax[0];
1177 Double_t ymin = -smax[1];
1178 Double_t ymax = smax[1];
1179 Double_t dx = xmax-xmin;
1180 Double_t dy = ymax-ymin;
1181 Double_t dxr = dx/(1 - gPad->GetLeftMargin() - gPad->GetRightMargin());
1182 Double_t dyr = dy/(1 - gPad->GetBottomMargin() - gPad->GetTopMargin());
1183
1184 // Range() could change the size of the pad pixmap and therefore should
1185 // be called before the other paint routines
1186 gPad->Range(xmin - dxr*gPad->GetLeftMargin(),
1187 ymin - dyr*gPad->GetBottomMargin(),
1188 xmax + dxr*gPad->GetRightMargin(),
1189 ymax + dyr*gPad->GetTopMargin());
1190 gPad->RangeAxis(xmin, ymin, xmax, ymax);
1191
1192 //*-*- Draw the background if necessary
1193 if (rback > 0) {
1194 r1 = -1;
1195 r2 = -1;
1196 r3 = -1;
1197 xgraf[0] = -smax[0];
1198 xgraf[1] = -smax[0];
1199 xgraf[2] = -smax[0];
1200 xgraf[3] = -smax[0];
1201 xgraf[4] = smax[0];
1202 xgraf[5] = smax[0];
1203 ygraf[0] = -smax[1];
1204 ygraf[1] = smax[1];
1205 ygraf[2] = -smax[1];
1206 ygraf[3] = smax[1];
1207 ygraf[5] = smax[1];
1208 ygraf[4] = -smax[1];
1209 for (i = 1; i <= 8; ++i) {
1210 x = 0.5*((1 - r1)*fRmin[0] + (r1 + 1)*fRmax[0]);
1211 y = 0.5*((1 - r2)*fRmin[1] + (r2 + 1)*fRmax[1]);
1212 z = 0.5*((1 - r3)*fRmin[2] + (r3 + 1)*fRmax[2]);
1213 xx = fTnorm[0]*x + fTnorm[1]*y + fTnorm[2]*z + fTnorm[3];
1214 yy = fTnorm[4]*x + fTnorm[5]*y + fTnorm[6]*z + fTnorm[7];
1215 if (TMath::Abs(xx - xgraf[1]) <= 1e-4) {
1216 if (ygraf[1] >= yy) ygraf[1] = yy;
1217 if (ygraf[2] <= yy) ygraf[2] = yy;
1218 }
1219 if (TMath::Abs(xx - xgraf[5]) <= 1e-4) {
1220 if (ygraf[5] >= yy) ygraf[5] = yy;
1221 if (ygraf[4] <= yy) ygraf[4] = yy;
1222 }
1223 if (TMath::Abs(yy - ygraf[0]) <= 1e-4) xgraf[0] = xx;
1224 if (TMath::Abs(yy - ygraf[3]) <= 1e-4) xgraf[3] = xx;
1225 r1 = -r1;
1226 if (i % 2 == 0) r2 = -r2;
1227 if (i >= 4) r3 = 1;
1228 }
1229 gPad->PaintFillArea(6, xgraf, ygraf);
1230 }
1231}
1232
1233////////////////////////////////////////////////////////////////////////////////
1234/// Store axis coordinates in the NDC system.
1235
1236void TView3D::SetAxisNDC(const Double_t *x1, const Double_t *x2, const Double_t *y1, const Double_t *y2, const Double_t *z1, const Double_t *z2)
1237{
1238 for (Int_t i=0;i<3;i++) {
1239 fX1[i] = x1[i];
1240 fX2[i] = x2[i];
1241 fY1[i] = y1[i];
1242 fY2[i] = y2[i];
1243 fZ1[i] = z1[i];
1244 fZ2[i] = z2[i];
1245 }
1246}
1247
1248////////////////////////////////////////////////////////////////////////////////
1249/// Set default viewing window.
1250
1252{
1253 if (!gPad) return;
1254 Double_t screen_factor = 1.;
1255 Double_t du, dv;
1256 Double_t extent = GetExtent();
1257 fDview = 3*extent;
1258 fDproj = 0.5*extent;
1259
1260 // width in pixels
1261 fUpix = gPad->GetWw()*gPad->GetAbsWNDC();
1262
1263 // height in pixels
1264 fVpix = gPad->GetWh()*gPad->GetAbsHNDC();
1265 du = 0.5*screen_factor*fDproj;
1266 dv = du*fVpix/fUpix; // keep aspect ratio
1267 SetWindow(0, 0, du, dv);
1268}
1269
1270////////////////////////////////////////////////////////////////////////////////
1271/// This is a function which creates default outline.
1272///
1273/// ~~~ {.cpp}
1274/// x = fRmin[0] X = fRmax[0]
1275/// y = fRmin[1] Y = fRmax[1]
1276/// z = fRmin[2] Z = fRmax[2]
1277///
1278/// (x,Y,Z) +---------+ (X,Y,Z)
1279/// / /|
1280/// / / |
1281/// / / |
1282/// (x,y,Z) +---------+ |
1283/// | | + (X,Y,z)
1284/// | | /
1285/// | | /
1286/// | |/
1287/// +---------+
1288/// (x,y,z) (X,y,z)
1289/// ~~~
1290
1292{
1293 if (!fOutline) {
1295 fOutline = new TList();
1296 }
1298}
1299
1300////////////////////////////////////////////////////////////////////////////////
1301/// Set the parallel option (default).
1302
1304{
1305 if (!IsPerspective()) return;
1307 Int_t irep;
1309}
1310
1311////////////////////////////////////////////////////////////////////////////////
1312/// Set perspective option.
1313
1315{
1316 if (IsPerspective()) return;
1318 Int_t irep;
1321}
1322
1323////////////////////////////////////////////////////////////////////////////////
1324/// Set Range function.
1325
1326void TView3D::SetRange(const Double_t *min, const Double_t *max)
1327{
1328 Int_t irep;
1329 for (Int_t i = 0; i < 3; fRmax[i] = max[i], fRmin[i] = min[i], i++) { }
1332 if(irep < 0)
1333 Error("SetRange", "problem setting view");
1335}
1336
1337////////////////////////////////////////////////////////////////////////////////
1338/// Set 3-D View range.
1339///
1340/// Input:
1341/// - x0, y0, z0 are minimum coordinates
1342/// - x1, y1, z1 are maximum coordinates
1343///
1344/// - flag values are:
1345/// - 0 (set always) <- default
1346/// - 1 (shrink view)
1347/// - 2 (expand view)
1348
1350{
1351 Double_t rmax[3], rmin[3];
1352
1353 switch (flag) {
1354 case 2: // expand view
1355 GetRange(rmin, rmax);
1356 rmin[0] = x0 < rmin[0] ? x0 : rmin[0];
1357 rmin[1] = y0 < rmin[1] ? y0 : rmin[1];
1358 rmin[2] = z0 < rmin[2] ? z0 : rmin[2];
1359 rmax[0] = x1 > rmax[0] ? x1 : rmax[0];
1360 rmax[1] = y1 > rmax[1] ? y1 : rmax[1];
1361 rmax[2] = z1 > rmax[2] ? z1 : rmax[2];
1362 break;
1363
1364 case 1: // shrink view
1365 GetRange(rmin, rmax);
1366 rmin[0] = x0 > rmin[0] ? x0 : rmin[0];
1367 rmin[1] = y0 > rmin[1] ? y0 : rmin[1];
1368 rmin[2] = z0 > rmin[2] ? z0 : rmin[2];
1369 rmax[0] = x1 < rmax[0] ? x1 : rmax[0];
1370 rmax[1] = y1 < rmax[1] ? y1 : rmax[1];
1371 rmax[2] = z1 < rmax[2] ? z1 : rmax[2];
1372 break;
1373
1374 default:
1375 rmin[0] = x0; rmax[0] = x1;
1376 rmin[1] = y0; rmax[1] = y1;
1377 rmin[2] = z0; rmax[2] = z1;
1378 }
1379 SetRange(rmin, rmax);
1380}
1381
1382////////////////////////////////////////////////////////////////////////////////
1383/// Set viewing window.
1384
1386{
1387 fUVcoord[0] = u0;
1388 fUVcoord[1] = v0;
1389 fUVcoord[2] = du;
1390 fUVcoord[3] = dv;
1391}
1392
1393////////////////////////////////////////////////////////////////////////////////
1394/// Set view parameters.
1395
1396void TView3D::SetView(Double_t longitude, Double_t latitude, Double_t psi, Int_t &irep)
1397{
1398 ResetView(longitude, latitude, psi, irep);
1399}
1400
1401////////////////////////////////////////////////////////////////////////////////
1402/// Recompute window for perspective view.
1403
1405{
1406 if (!IsPerspective()) return;
1407 Double_t upix = fUpix;
1408 Double_t vpix = fVpix;
1409
1410 // width in pixels
1411 fUpix = gPad->GetWw()*gPad->GetAbsWNDC();
1412
1413 // height in pixels
1414 fVpix = gPad->GetWh()*gPad->GetAbsHNDC();
1415 Double_t u0 = fUVcoord[0]*fUpix/upix;
1416 Double_t v0 = fUVcoord[1]*fVpix/vpix;
1417 Double_t du = fUVcoord[2]*fUpix/upix;
1418 Double_t dv = fUVcoord[3]*fVpix/vpix;
1419 SetWindow(u0, v0, du, dv);
1421}
1422
1423////////////////////////////////////////////////////////////////////////////////
1424/// Set view direction (in spherical coordinates).
1425///
1426/// Input
1427/// - PHI - longitude
1428/// - THETA - latitude (angle between +Z and view direction)
1429/// - PSI - rotation in screen plane
1430///
1431/// Output:
1432/// - IREP - reply (-1 if error in min-max)
1433///
1434/// Errors: error in min-max scope
1435
1436void TView3D::ResetView(Double_t longitude, Double_t latitude, Double_t psi, Int_t &irep)
1437{
1438 Double_t scale[3], centre[3];
1439 Double_t c1, c2, c3, s1, s2, s3;
1440
1441 //*-*- F I N D C E N T E R O F S C O P E A N D
1442 //*-*- S C A L E F A C T O R S
1443 FindScope(scale, centre, irep);
1444 if (irep < 0) {
1445 Error("ResetView", "Error in min-max scope");
1446 return;
1447 }
1448
1449 //*-*- S E T T R A N S F O R M A T I O N M A T R I C E S
1450 fLongitude = longitude;
1451 fPsi = psi;
1452 fLatitude = latitude;
1453
1454 if (IsPerspective()) {
1456 return;
1457 }
1458
1459 c1 = TMath::Cos(longitude*kRad);
1460 s1 = TMath::Sin(longitude*kRad);
1461 c2 = TMath::Cos(latitude*kRad);
1462 s2 = TMath::Sin(latitude*kRad);
1463 c3 = TMath::Cos(psi*kRad);
1464 s3 = TMath::Sin(psi*kRad);
1465 DefineViewDirection(scale, centre, c1, s1, c2, s2, c3, s3, fTnorm, fTback);
1466 c3 = 1;
1467 s3 = 0;
1468 DefineViewDirection(scale, centre, c1, s1, c2, s2, c3, s3, fTN, fTB);
1469}
1470
1471////////////////////////////////////////////////////////////////////////////////
1472/// Transfer point from world to normalized coordinates.
1473///
1474/// Input:
1475/// - PW(3) - point in world coordinate system
1476/// - PN(3) - point in normalized coordinate system
1477
1479{
1480 Float_t x = pw[0], y = pw[1], z = pw[2];
1481
1482 // perspective view
1483 if (IsPerspective()) {
1484 for (Int_t i=0; i<3; i++) {
1485 pn[i] = fTnorm[i]*x + fTnorm[i+4]*y + fTnorm[i+8]*z + fTnorm[i+12];
1486 }
1487 if (pn[2]>0) {
1488 pn[0] /= pn[2];
1489 pn[1] /= pn[2];
1490 } else {
1491 pn[0] *= 1000.;
1492 pn[1] *= 1000.;
1493 }
1494 return;
1495 }
1496
1497 // parallel view
1498 pn[0] = fTnorm[0]*x + fTnorm[1]*y + fTnorm[2]*z + fTnorm[3];
1499 pn[1] = fTnorm[4]*x + fTnorm[5]*y + fTnorm[6]*z + fTnorm[7];
1500 pn[2] = fTnorm[8]*x + fTnorm[9]*y + fTnorm[10]*z + fTnorm[11];
1501}
1502
1503////////////////////////////////////////////////////////////////////////////////
1504/// Transfer point from world to normalized coordinates.
1505///
1506/// Input:
1507/// - PW(3) - point in world coordinate system
1508/// - PN(3) - point in normalized coordinate system
1509
1511{
1512 Double_t x = pw[0], y = pw[1], z = pw[2];
1513
1514 // perspective view
1515 if (IsPerspective()) {
1516 for (Int_t i=0; i<3; i++) {
1517 pn[i] = fTnorm[i]*x + fTnorm[i+4]*y + fTnorm[i+8]*z + fTnorm[i+12];
1518 }
1519 if (pn[2]>0) {
1520 pn[0] /= pn[2];
1521 pn[1] /= pn[2];
1522 } else {
1523 pn[0] *= 1000.;
1524 pn[1] *= 1000.;
1525 }
1526 return;
1527 }
1528
1529 // parallel view
1530 pn[0] = fTnorm[0]*x + fTnorm[1]*y + fTnorm[2]*z + fTnorm[3];
1531 pn[1] = fTnorm[4]*x + fTnorm[5]*y + fTnorm[6]*z + fTnorm[7];
1532 pn[2] = fTnorm[8]*x + fTnorm[9]*y + fTnorm[10]*z + fTnorm[11];
1533}
1534
1535////////////////////////////////////////////////////////////////////////////////
1536/// Force the current pad to be updated.
1537
1539{
1540 TVirtualPad *thisPad = pad;
1541 if (!thisPad) thisPad = gPad;
1542 if (thisPad) {
1543#ifdef R__HAS_COCOA
1544 thisPad->AbsCoordinates(kFALSE);
1545#endif
1546 thisPad->Modified();
1547 thisPad->Update();
1548 }
1549}
1550
1551////////////////////////////////////////////////////////////////////////////////
1552/// API to rotate view and adjust the pad provided it the current one.
1553
1555{
1556 Int_t iret;
1557 Double_t p = phi;
1558 Double_t t = theta;
1559 SetView(p, t, 0, iret);
1560
1561 // Adjust current pad too
1562 TVirtualPad *thisPad = pad;
1563 if (!thisPad) thisPad = gPad;
1564 if (thisPad) {
1565 thisPad->SetPhi(-90-p);
1566 thisPad->SetTheta(90-t);
1567 thisPad->Modified();
1568 thisPad->Update();
1569 }
1570}
1571
1572////////////////////////////////////////////////////////////////////////////////
1573/// Set to side view.
1574
1576{
1577 RotateView(0,90.0,pad);
1578}
1579
1580////////////////////////////////////////////////////////////////////////////////
1581/// Set to front view.
1582
1584{
1585 RotateView(270.0,90.0,pad);
1586}
1587
1588////////////////////////////////////////////////////////////////////////////////
1589/// Set to top view.
1590
1592{
1593 RotateView(270.0,0.0,pad);
1594}
1595
1596////////////////////////////////////////////////////////////////////////////////
1597/// Turn on /off 3D axis.
1598
1600{
1602}
1603
1604////////////////////////////////////////////////////////////////////////////////
1605/// Turn on /off the interactive option to
1606/// Zoom / Move / Change attributes of 3D axis correspond this view.
1607
1609{
1611}
1612
1613////////////////////////////////////////////////////////////////////////////////
1614/// Adjust all sides of view in respect of the biggest one.
1615
1617{
1618 Double_t min[3],max[3];
1619 GetRange(min,max);
1620 int i;
1621 Double_t maxSide = 0;
1622 // Find the largest side
1623 for (i=0;i<3; i++) maxSide = TMath::Max(maxSide,max[i]-min[i]);
1624 //Adjust scales:
1625 for (i=0;i<3; i++) max[i] += maxSide - (max[i]-min[i]);
1626 SetRange(min,max);
1627
1628 AdjustPad(pad);
1629}
1630
1631////////////////////////////////////////////////////////////////////////////////
1632/// Move view into the center of the scene.
1633
1635{
1636 Double_t min[3],max[3];
1637 GetRange(min,max);
1638 int i;
1639 for (i=0;i<3; i++) {
1640 if (max[i] > 0) min[i] = -max[i];
1641 else max[i] = -min[i];
1642 }
1643 SetRange(min,max);
1644 AdjustPad(pad);
1645}
1646
1647////////////////////////////////////////////////////////////////////////////////
1648/// unZOOM this view.
1649
1651{
1652 if (TMath::Abs(unZoomFactor) < 0.001) return;
1653 ZoomView(pad,1./unZoomFactor);
1654}
1655
1656////////////////////////////////////////////////////////////////////////////////
1657/// ZOOM this view.
1658
1660{
1661 if (TMath::Abs(zoomFactor) < 0.001) return;
1662 Double_t min[3],max[3];
1663 GetRange(min,max);
1664 int i;
1665 for (i=0;i<3; i++) {
1666 // Find center
1667 Double_t c = (max[i]+min[i])/2;
1668 // Find a new size
1669 Double_t s = (max[i]-min[i])/(2*zoomFactor);
1670 // Set a new size
1671 max[i] = c + s;
1672 min[i] = c - s;
1673 }
1674 SetRange(min,max);
1675 AdjustPad(pad);
1676}
1677
1678////////////////////////////////////////////////////////////////////////////////
1679/// Move focus to a different box position and extent in nsteps. Perform
1680/// rotation with dlat,dlong,dpsi at each step.
1681
1683 Double_t dlong, Double_t dlat, Double_t dpsi)
1684{
1685 if (!IsPerspective()) return;
1686 if (nsteps<1) return;
1687 Double_t fc = 1./Double_t(nsteps);
1688 Double_t oc[3], od[3], dir[3];
1689 dir[0] = 0;
1690 dir[1] = 0;
1691 dir[2] = 1.;
1692 Int_t i, j;
1693 for (i=0; i<3; i++) {
1694 oc[i] = 0.5*(fRmin[i]+fRmax[i]);
1695 od[i] = 0.5*(fRmax[i]-fRmin[i]);
1696 }
1697 Double_t dox = cov[0]-oc[0];
1698 Double_t doy = cov[1]-oc[1];
1699 Double_t doz = cov[2]-oc[2];
1700
1701 Double_t dd = TMath::Sqrt(dox*dox+doy*doy+doz*doz);
1702 if (dd!=0) {;
1703 dir[0] = dox/dd;
1704 dir[1] = doy/dd;
1705 dir[2] = doz/dd;
1706 }
1707 dd *= fc;
1708 dox = fc*(dx-od[0]);
1709 doy = fc*(dy-od[1]);
1710 doz = fc*(dz-od[2]);
1711 for (i=0; i<nsteps; i++) {
1712 oc[0] += dd*dir[0];
1713 oc[1] += dd*dir[1];
1714 oc[2] += dd*dir[2];
1715 od[0] += dox;
1716 od[1] += doy;
1717 od[2] += doz;
1718 for (j=0; j<3; j++) {
1719 fRmin[j] = oc[j]-od[j];
1720 fRmax[j] = oc[j]+od[j];
1721 }
1723 fLatitude += dlat;
1724 fLongitude += dlong;
1725 fPsi += dpsi;
1727 if (gPad) {
1728 gPad->Modified();
1729 gPad->Update();
1730 }
1731 }
1732}
1733
1734////////////////////////////////////////////////////////////////////////////////
1735/// - 'a' increase scale factor (clip cube borders)
1736/// - 's' decrease scale factor (clip cube borders)
1737
1739{
1740 if (count <= 0) count = 1;
1741 switch (option) {
1742 case '+':
1743 ZoomView();
1744 break;
1745 case '-':
1746 UnzoomView();
1747 break;
1748 case 's':
1749 case 'S':
1750 UnzoomView();
1751 break;
1752 case 'a':
1753 case 'A':
1754 ZoomView();
1755 break;
1756 case 'l':
1757 case 'L':
1758 case 'h':
1759 case 'H':
1760 case 'u':
1761 case 'U':
1762 case 'i':
1763 case 'I':
1764 MoveWindow(option);
1765 break;
1766 case 'j':
1767 case 'J':
1768 ZoomIn();
1769 break;
1770 case 'k':
1771 case 'K':
1772 ZoomOut();
1773 break;
1774 default:
1775 break;
1776 }
1777}
1778
1779////////////////////////////////////////////////////////////////////////////////
1780/// Move view window :
1781/// - l,L - left
1782/// - h,H - right
1783/// - u,U - down
1784/// - i,I - up
1785
1787{
1788 if (!IsPerspective()) return;
1789 Double_t shiftu = 0.1*fUVcoord[2];
1790 Double_t shiftv = 0.1*fUVcoord[3];
1791 switch (option) {
1792 case 'l':
1793 case 'L':
1794 fUVcoord[0] += shiftu;
1795 break;
1796 case 'h':
1797 case 'H':
1798 fUVcoord[0] -= shiftu;
1799 break;
1800 case 'u':
1801 case 'U':
1802 fUVcoord[1] += shiftv;
1803 break;
1804 case 'i':
1805 case 'I':
1806 fUVcoord[1] -= shiftv;
1807 break;
1808 default:
1809 return;
1810 }
1812 if (gPad) {
1813 gPad->Modified();
1814 gPad->Update();
1815 }
1816}
1817
1818////////////////////////////////////////////////////////////////////////////////
1819/// Zoom in.
1820
1822{
1823 if (!IsPerspective()) return;
1824 Double_t extent = GetExtent();
1825 Double_t fc = 0.1;
1826 if (fDview<extent) {
1827 fDview -= fc*extent;
1828 } else {
1829 fDview /= 1.25;
1830 }
1832 if (gPad) {
1833 gPad->Modified();
1834 gPad->Update();
1835 }
1836}
1837
1838////////////////////////////////////////////////////////////////////////////////
1839/// Zoom out.
1840
1842{
1843 if (!IsPerspective()) return;
1844 Double_t extent = GetExtent();
1845 Double_t fc = 0.1;
1846 if (fDview<extent) {
1847 fDview += fc*extent;
1848 } else {
1849 fDview *= 1.25;
1850 }
1852 if (gPad) {
1853 gPad->Modified();
1854 gPad->Update();
1855 }
1856}
1857
1858////////////////////////////////////////////////////////////////////////////////
1859/// Stream an object of class TView3D.
1860
1861void TView3D::Streamer(TBuffer &R__b)
1862{
1863 if (R__b.IsReading()) {
1864 UInt_t R__s, R__c;
1865 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
1866 if (R__v > 1) {
1867 R__b.ReadClassBuffer(TView3D::Class(), this, R__v, R__s, R__c);
1868 return;
1869 }
1870 //====process old versions before automatic schema evolution
1871 //unfortunately we forgot to increment the TView3D version number
1872 //when the class was upgraded to double precision in version 2.25.
1873 //we are forced to use the file version number to recognize old files.
1874 if (R__b.GetParent() && R__b.GetVersionOwner() < 22500) { //old version in single precision
1875 TObject::Streamer(R__b);
1876 TAttLine::Streamer(R__b);
1877 Float_t single, sa[12];
1878 Int_t i;
1879 R__b >> fSystem;
1880 R__b >> single; fLatitude = single;
1881 R__b >> single; fLongitude = single;
1882 R__b >> single; fPsi = single;
1883 R__b.ReadStaticArray(sa); for (i=0;i<12;i++) fTN[i] = sa[i];
1884 R__b.ReadStaticArray(sa); for (i=0;i<12;i++) fTB[i] = sa[i];
1885 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fRmax[i] = sa[i];
1886 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fRmin[i] = sa[i];
1887 R__b.ReadStaticArray(sa); for (i=0;i<12;i++) fTnorm[i] = sa[i];
1888 R__b.ReadStaticArray(sa); for (i=0;i<12;i++) fTback[i] = sa[i];
1889 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fX1[i] = sa[i];
1890 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fX2[i] = sa[i];
1891 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fY1[i] = sa[i];
1892 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fY2[i] = sa[i];
1893 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fZ1[i] = sa[i];
1894 R__b.ReadStaticArray(sa); for (i=0;i<3;i++) fZ2[i] = sa[i];
1895 R__b >> fOutline;
1896 R__b >> fDefaultOutline;
1897 R__b >> fAutoRange;
1898 } else {
1899 TObject::Streamer(R__b);
1900 TAttLine::Streamer(R__b);
1901 R__b >> fLatitude;
1902 R__b >> fLongitude;
1903 R__b >> fPsi;
1904 R__b.ReadStaticArray(fTN);
1905 R__b.ReadStaticArray(fTB);
1906 R__b.ReadStaticArray(fRmax);
1907 R__b.ReadStaticArray(fRmin);
1908 R__b.ReadStaticArray(fTnorm);
1909 R__b.ReadStaticArray(fTback);
1910 R__b.ReadStaticArray(fX1);
1911 R__b.ReadStaticArray(fX2);
1912 R__b.ReadStaticArray(fY1);
1913 R__b.ReadStaticArray(fY2);
1914 R__b.ReadStaticArray(fZ1);
1915 R__b.ReadStaticArray(fZ2);
1916 R__b >> fSystem;
1917 R__b >> fOutline;
1918 R__b >> fDefaultOutline;
1919 R__b >> fAutoRange;
1920 }
1921 //====end of old versions
1922
1923 } else {
1924 R__b.WriteClassBuffer(TView3D::Class(),this);
1925 }
1926}
1927
1928// Shortcuts for menus
1937
1938
@ kMouseMotion
Definition: Buttons.h:23
@ kKeyPress
Definition: Buttons.h:20
@ kButton1Motion
Definition: Buttons.h:20
@ kButton1Up
Definition: Buttons.h:19
@ kButton1Down
Definition: Buttons.h:17
void Class()
Definition: Class.C:29
SVector< double, 2 > v
Definition: Dict.h:5
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
#define c(i)
Definition: RSha256.hxx:101
#define s1(x)
Definition: RSha256.hxx:91
#define e(i)
Definition: RSha256.hxx:103
static const double x2[5]
static const double x1[5]
int Int_t
Definition: RtypesCore.h:41
short Version_t
Definition: RtypesCore.h:61
char Char_t
Definition: RtypesCore.h:29
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
float xmin
Definition: THbookFile.cxx:93
float ymin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
float ymax
Definition: THbookFile.cxx:93
#define gROOT
Definition: TROOT.h:410
const Double_t kRad
Definition: TView3D.cxx:34
const Int_t kCARTESIAN
Definition: TView3D.cxx:32
const Int_t kPOLAR
Definition: TView3D.cxx:33
#define gPad
Definition: TVirtualPad.h:286
#define gVirtualX
Definition: TVirtualX.h:345
@ kRotate
Definition: TVirtualX.h:46
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:3728
static TAxis3D * ToggleZoom(TVirtualPad *pad=0)
Turn ON / OFF the "Ruler" and "zoom mode" of the TAxis3D object attached to the current pad (if pad =...
Definition: TAxis3D.cxx:766
static TAxis3D * ToggleRulers(TVirtualPad *pad=0)
Turn ON / OFF the "Ruler", TAxis3D object attached to the current pad.
Definition: TAxis3D.cxx:738
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
TObject * GetParent() const
Return pointer to parent of this buffer.
Definition: TBuffer.cxx:241
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t GetVersionOwner() const =0
Bool_t IsReading() const
Definition: TBuffer.h:83
virtual Int_t ReadStaticArray(Bool_t *b)=0
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual void Delete(Option_t *option="")=0
Delete this object.
virtual void Paint(Option_t *option="")
Paint all objects in this collection.
A doubly linked list.
Definition: TList.h:44
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.h:271
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition: TObject.h:60
static void DrawOutlineCube(TList *outline, Double_t *rmin, Double_t *rmax)
Draw cube outline with 3d polylines.
The 3D view class.
Definition: TView3D.h:29
@ kPerspective
Definition: TView3D.h:67
virtual void Zoom()
Definition: TView3D.cxx:1935
Double_t fDproj
Definition: TView3D.h:36
Bool_t fDefaultOutline
Definition: TView3D.h:54
virtual void UnzoomView(TVirtualPad *pad=0, Double_t unZoomFactor=1.25)
unZOOM this view.
Definition: TView3D.cxx:1650
Double_t fTB[16]
Definition: TView3D.h:40
TView3D()
Default constructor.
Definition: TView3D.cxx:108
Double_t fUpix
Definition: TView3D.h:37
virtual void Side()
Definition: TView3D.cxx:1932
Double_t fX2[3]
Definition: TView3D.h:47
Bool_t fChanged
Definition: TView3D.h:56
static void AdjustPad(TVirtualPad *pad=0)
Force the current pad to be updated.
Definition: TView3D.cxx:1538
virtual void ZoomView(TVirtualPad *pad=0, Double_t zoomFactor=1.25)
ZOOM this view.
Definition: TView3D.cxx:1659
virtual void SetPerspective()
Set perspective option.
Definition: TView3D.cxx:1314
Double_t fPsi
Definition: TView3D.h:34
virtual ~TView3D()
TView3D default destructor.
Definition: TView3D.cxx:275
virtual void Centered3DImages(TVirtualPad *pad=0)
Move view into the center of the scene.
Definition: TView3D.cxx:1634
virtual Int_t GetSystem()
Definition: TView3D.h:101
virtual void SetParallel()
Set the parallel option (default).
Definition: TView3D.cxx:1303
Double_t fLongitude
Definition: TView3D.h:33
virtual void RotateView(Double_t phi, Double_t theta, TVirtualPad *pad=0)
API to rotate view and adjust the pad provided it the current one.
Definition: TView3D.cxx:1554
Double_t fZ1[3]
Definition: TView3D.h:50
Double_t fUVcoord[4]
Definition: TView3D.h:43
virtual void FindPhiSectors(Int_t iopt, Int_t &kphi, Double_t *aphi, Int_t &iphi1, Int_t &iphi2)
Find critical PHI sectors.
Definition: TView3D.cxx:805
Double_t fZ2[3]
Definition: TView3D.h:51
TView3D & operator=(const TView3D &)
Assignment operator.
Definition: TView3D.cxx:234
Double_t fY1[3]
Definition: TView3D.h:48
virtual void PadRange(Int_t rback)
Set the correct window size for lego and surface plots.
Definition: TView3D.cxx:1157
virtual Int_t GetDistancetoAxis(Int_t axis, Int_t px, Int_t py, Double_t &ratio)
Return distance to axis from point px,py.
Definition: TView3D.cxx:970
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)
Transfer point from world to normalized coordinates.
Definition: TView3D.cxx:1478
Double_t fTback[16]
Definition: TView3D.h:45
Double_t fTnorm[16]
Definition: TView3D.h:44
void ResetView(Double_t longitude, Double_t latitude, Double_t psi, Int_t &irep)
Set view direction (in spherical coordinates).
Definition: TView3D.cxx:1436
virtual void DrawOutlineCube(TList *outline, Double_t *rmin, Double_t *rmax)
Draw the outline of a cube while rotating a 3-d object in the pad.
Definition: TView3D.cxx:619
Double_t fLatitude
Definition: TView3D.h:32
virtual void NormalWCtoNDC(const Float_t *pw, Float_t *pn)
Transfer vector of NORMAL from word to normalized coordinates.
Definition: TView3D.cxx:1099
virtual void MoveFocus(Double_t *center, Double_t dx, Double_t dy, Double_t dz, Int_t nsteps=10, Double_t dlong=0, Double_t dlat=0, Double_t dpsi=0)
Move focus to a different box position and extent in nsteps.
Definition: TView3D.cxx:1682
virtual void DefinePerspectiveView()
Define perspective view.
Definition: TView3D.cxx:399
virtual Double_t GetPsi()
Definition: TView3D.h:92
virtual void AdjustScales(TVirtualPad *pad=0)
Adjust all sides of view in respect of the biggest one.
Definition: TView3D.cxx:1616
virtual void MoveWindow(Char_t option)
Move view window :
Definition: TView3D.cxx:1786
virtual void SetAxisNDC(const Double_t *x1, const Double_t *x2, const Double_t *y1, const Double_t *y2, const Double_t *z1, const Double_t *z2)
Store axis coordinates in the NDC system.
Definition: TView3D.cxx:1236
virtual void SetWindow(Double_t u0, Double_t v0, Double_t du, Double_t dv)
Set viewing window.
Definition: TView3D.cxx:1385
virtual void SideView(TVirtualPad *pad=0)
Set to side view.
Definition: TView3D.cxx:1575
virtual void SetView(Double_t longitude, Double_t latitude, Double_t psi, Int_t &irep)
Set view parameters.
Definition: TView3D.cxx:1396
virtual void FindScope(Double_t *scale, Double_t *center, Int_t &irep)
Find centre of a MIN-MAX scope and scale factors.
Definition: TView3D.cxx:932
virtual void SetRange(const Double_t *min, const Double_t *max)
Set Range function.
Definition: TView3D.cxx:1326
virtual Bool_t IsClippedNDC(Double_t *p) const
Check if point is clipped in perspective view.
Definition: TView3D.cxx:1055
virtual void TopView(TVirtualPad *pad=0)
Set to top view.
Definition: TView3D.cxx:1591
virtual void Top()
Definition: TView3D.cxx:1933
virtual void NDCtoWC(const Float_t *pn, Float_t *pw)
Transfer point from normalized to world coordinates.
Definition: TView3D.cxx:1069
virtual void ResizePad()
Recompute window for perspective view.
Definition: TView3D.cxx:1404
Double_t fY2[3]
Definition: TView3D.h:49
virtual void SetDefaultWindow()
Set default viewing window.
Definition: TView3D.cxx:1251
virtual void Front()
Definition: TView3D.cxx:1930
Double_t fX1[3]
Definition: TView3D.h:46
virtual void GetWindow(Double_t &u0, Double_t &v0, Double_t &du, Double_t &dv) const
Get current window extent.
Definition: TView3D.cxx:1044
virtual void ExecuteRotateView(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TView3D.cxx:641
virtual void FindNormal(Double_t x, Double_t y, Double_t z, Double_t &zn)
Find Z component of NORMAL in normalized coordinates.
Definition: TView3D.cxx:785
Double_t fTN[16]
Definition: TView3D.h:39
virtual void UnZoom()
Definition: TView3D.cxx:1936
virtual Bool_t IsPerspective() const
Definition: TView3D.h:109
Bool_t fAutoRange
Definition: TView3D.h:55
virtual void ZoomMove()
Definition: TView3D.cxx:1934
virtual void ToggleZoom(TVirtualPad *pad=0)
Turn on /off the interactive option to Zoom / Move / Change attributes of 3D axis correspond this vie...
Definition: TView3D.cxx:1608
Double_t fRmin[3]
Definition: TView3D.h:42
Double_t fRmax[3]
Definition: TView3D.h:41
virtual void MoveViewCommand(Char_t chCode, Int_t count=1)
Definition: TView3D.cxx:1738
virtual void ToggleRulers(TVirtualPad *pad=0)
Turn on /off 3D axis.
Definition: TView3D.cxx:1599
virtual void ZoomOut()
Zoom out.
Definition: TView3D.cxx:1841
Double_t fDview
Definition: TView3D.h:35
virtual void Centered()
Definition: TView3D.cxx:1929
virtual void FrontView(TVirtualPad *pad=0)
Set to front view.
Definition: TView3D.cxx:1583
virtual Double_t GetExtent() const
Get maximum view extent.
Definition: TView3D.cxx:1016
TSeqCollection * fOutline
Definition: TView3D.h:53
virtual void DefineViewDirection(const Double_t *s, const Double_t *c, Double_t cosphi, Double_t sinphi, Double_t costhe, Double_t sinthe, Double_t cospsi, Double_t sinpsi, Double_t *tnorm, Double_t *tback)
Define view direction (in spherical coordinates)
Definition: TView3D.cxx:513
virtual void GetRange(Float_t *min, Float_t *max)
Get Range function.
Definition: TView3D.cxx:1028
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TView3D.cxx:627
virtual void ZoomIn()
Zoom in.
Definition: TView3D.cxx:1821
virtual void FindThetaSectors(Int_t iopt, Double_t phi, Int_t &kth, Double_t *ath, Int_t &ith1, Int_t &ith2)
Find critical THETA sectors for given PHI sector.
Definition: TView3D.cxx:871
Int_t fSystem
Definition: TView3D.h:52
virtual void SetOutlineToCube()
This is a function which creates default outline.
Definition: TView3D.cxx:1291
Double_t fVpix
Definition: TView3D.h:38
virtual void AxisVertex(Double_t ang, Double_t *av, Int_t &ix1, Int_t &ix2, Int_t &iy1, Int_t &iy2, Int_t &iz1, Int_t &iz2)
Define axis vertices.
Definition: TView3D.cxx:297
virtual void ShowAxis()
Definition: TView3D.cxx:1931
See TView3D.
Definition: TView.h:25
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:50
virtual void Modified(Bool_t flag=1)=0
virtual void AbsCoordinates(Bool_t set)=0
virtual void Update()=0
virtual void SetTheta(Double_t theta=30)=0
virtual void SetPhi(Double_t phi=30)=0
Abstract 3D shapes viewer.
return c1
Definition: legend1.C:41
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
return c2
Definition: legend2.C:14
return c3
Definition: legend3.C:15
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
static constexpr double s
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
Double_t Sqrt(Double_t x)
Definition: TMath.h:679
Double_t Cos(Double_t)
Definition: TMath.h:629
Double_t Sin(Double_t)
Definition: TMath.h:625
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * th2
Definition: textalign.C:17
auto * th1
Definition: textalign.C:13
auto * a
Definition: textangle.C:12