Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoPainter.cxx
Go to the documentation of this file.
1// @(#)root/geompainter:$Id: 58726ead32989b65bb2cbff2af4235fe9c6b12ae $
2// Author: Andrei Gheata 05/03/02
3/*************************************************************************
4 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11/** \class TGeoPainter
12\ingroup Geometry_painter
13
14Class implementing all draw interfaces for a generic 3D viewer
15using TBuffer3D mechanism.
16*/
17
18#include <map>
19#include "TROOT.h"
20#include "TClass.h"
21#include "TColor.h"
22#include "TPoint.h"
23#include "TView.h"
24#include "TAttLine.h"
25#include "TAttFill.h"
26#include "TVirtualPad.h"
27#include "TStopwatch.h"
28#include "TCanvas.h"
29#include "TCanvasImp.h"
30#include "TF1.h"
31#include "TGraph.h"
32#include "TPluginManager.h"
33#include "TVirtualPadEditor.h"
34
35#include "TPolyMarker3D.h"
36
37#include "TGeoAtt.h"
38#include "TGeoVolume.h"
39#include "TGeoNode.h"
40#include "TGeoElement.h"
41#include "TGeoManager.h"
42#include "TGeoTrack.h"
43#include "TGeoOverlap.h"
44#include "TGeoPhysicalNode.h"
45#include "TGeoPolygon.h"
46#include "TGeoCompositeShape.h"
47#include "TGeoShapeAssembly.h"
48#include "TGeoPainter.h"
49#include "TMath.h"
50#include "TVirtualGeoChecker.h"
51
52#include "X3DBuffer.h"
53
54#include "TBuffer3D.h"
55#include "TBuffer3DTypes.h"
56#include "TVirtualViewer3D.h"
57#include "TVirtualX.h"
58
59
60////////////////////////////////////////////////////////////////////////////////
61/// Default constructor.
62
64{
66 if (manager)
68 else {
69 Error("ctor", "No geometry loaded");
70 return;
71 }
73 fNVisNodes = 0;
74 fBombX = 1.3;
75 fBombY = 1.3;
76 fBombZ = 1.3;
77 fBombR = 1.3;
81 fVisBranch = "";
82 fVolInfo = "";
87 fPlugin = nullptr;
88 fVisVolumes = new TObjArray();
89 fOverlap = nullptr;
90 fGlobal = new TGeoHMatrix();
91 fBuffer = new TBuffer3D(TBuffer3DTypes::kGeneric, 20, 3 * 20, 0, 0, 0, 0);
92 fClippingShape = nullptr;
93 fLastVolume = nullptr;
94 fTopVolume = nullptr;
96 memset(&fCheckedBox[0], 0, 6 * sizeof(Double_t));
97
100 DefineColors();
101}
102////////////////////////////////////////////////////////////////////////////////
103/// Default destructor.
104
106{
107 delete fVisVolumes;
108 delete fGlobal;
109 delete fBuffer;
110 if (fPlugin)
111 delete fPlugin;
112}
113////////////////////////////////////////////////////////////////////////////////
114/// Add numpoints, numsegs, numpolys to the global 3D size.
115
117{
118 gSize3D.numPoints += numpoints;
119 gSize3D.numSegs += numsegs;
120 gSize3D.numPolys += numpolys;
121}
122////////////////////////////////////////////////////////////////////////////////
123/// Create a primary TGeoTrack.
124
129
130////////////////////////////////////////////////////////////////////////////////
131/// Average center of view of all painted tracklets and compute view box.
132
134{
135 static Int_t npoints = 0;
136 static Double_t xmin[3] = {0, 0, 0};
137 static Double_t xmax[3] = {0, 0, 0};
138 Int_t i;
139 if (reset) {
140 memset(box, 0, 6 * sizeof(Double_t));
141 memset(xmin, 0, 3 * sizeof(Double_t));
142 memset(xmax, 0, 3 * sizeof(Double_t));
143 npoints = 0;
144 return;
145 }
146 if (npoints == 0) {
147 for (i = 0; i < 3; i++)
148 xmin[i] = xmax[i] = 0;
149 npoints++;
150 }
151 npoints++;
153 for (i = 0; i < 3; i++) {
154 box[i] += ninv * (point[i] - box[i]);
155 if (point[i] < xmin[i])
156 xmin[i] = point[i];
157 if (point[i] > xmax[i])
158 xmax[i] = point[i];
159 box[i + 3] = 0.5 * (xmax[i] - xmin[i]);
160 }
161}
162
163////////////////////////////////////////////////////////////////////////////////
164/// Get the new 'bombed' translation vector according current exploded view mode.
165
167{
168 memcpy(bombtr, tr, 3 * sizeof(Double_t));
169 switch (fExplodedView) {
170 case kGeoNoBomb: return;
171 case kGeoBombXYZ:
172 bombtr[0] *= fBombX;
173 bombtr[1] *= fBombY;
174 bombtr[2] *= fBombZ;
175 return;
176 case kGeoBombCyl:
177 bombtr[0] *= fBombR;
178 bombtr[1] *= fBombR;
179 bombtr[2] *= fBombZ;
180 return;
181 case kGeoBombSph:
182 bombtr[0] *= fBombR;
183 bombtr[1] *= fBombR;
184 bombtr[2] *= fBombR;
185 return;
186 default: return;
187 }
188}
189
190////////////////////////////////////////////////////////////////////////////////
191/// Clear the list of visible volumes
192/// reset the kVisOnScreen bit for volumes previously in the list
193
195{
196 if (!fVisVolumes)
197 return;
198 TIter next(fVisVolumes);
199 TGeoVolume *vol;
200 while ((vol = (TGeoVolume *)next())) {
202 }
204}
205
206////////////////////////////////////////////////////////////////////////////////
207/// Define 100 colors with increasing light intensities for each basic color (1-7)
208/// Register these colors at indexes starting with 1000.
209
211{
212 static Int_t color = 0;
213 if (!color) {
215 for (auto icol = 1; icol < 10; ++icol)
216 color = GetColor(icol, 0.5);
217 }
218}
219
220////////////////////////////////////////////////////////////////////////////////
221/// Get index of a base color with given light intensity (0,1)
222
224{
225 using IntMap_t = std::map<Int_t, Int_t>;
226 constexpr Int_t ncolors = 100;
227 constexpr Float_t lmin = 0.25;
228 constexpr Float_t lmax = 0.75;
229 static IntMap_t colmap;
230 Int_t color = base;
231 // Search color in the map
232 auto it = colmap.find(base);
233 if (it != colmap.end())
234 return (it->second + light * (ncolors - 1));
235 // Get color pointer if stored
236 TColor *col_base = gROOT->GetColor(base);
237 if (!col_base) {
238 // If color not defined, use gray palette
239 it = colmap.find(kBlack);
240 if (it != colmap.end())
241 return (it->second + light * (ncolors - 1));
242 col_base = gROOT->GetColor(kBlack);
243 color = 1;
244 }
245 // Create a color palette for col_base
246 Float_t r = 0., g = 0., b = 0., h = 0., l = 0., s = 0.;
247 Double_t red[2], green[2], blue[2];
248 Double_t stop[] = {0., 1.0};
249
250 if (col_base)
251 col_base->GetRGB(r, g, b);
252 TColor::RGB2HLS(r, g, b, h, l, s);
253 TColor::HLS2RGB(h, lmin, s, r, g, b);
254 red[0] = r;
255 green[0] = g;
256 blue[0] = b;
257 TColor::HLS2RGB(h, lmax, s, r, g, b);
258 red[1] = r;
259 green[1] = g;
260 blue[1] = b;
262 colmap[color] = color_map_idx;
263 return (color_map_idx + light * (ncolors - 1));
264}
265
266////////////////////////////////////////////////////////////////////////////////
267/// Get currently drawn volume.
268
270{
271 if (!gPad)
272 return nullptr;
273 return fTopVolume;
274}
275
276////////////////////////////////////////////////////////////////////////////////
277/// Compute the closest distance of approach from point px,py to a volume.
278
280{
281 const Int_t big = 9999;
282 const Int_t inaxis = 7;
283 const Int_t maxdist = 5;
284
285 if (fTopVolume != volume)
286 fTopVolume = volume;
287 TView *view = gPad->GetView();
288 if (!view)
289 return big;
290 TGeoBBox *box;
291 fGlobal->Clear();
293
294 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
295 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
296 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
297 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
298 // return if point not in user area
299 if (px < puxmin - inaxis)
300 return big;
301 if (py > puymin + inaxis)
302 return big;
303 if (px > puxmax + inaxis)
304 return big;
305 if (py < puymax - inaxis)
306 return big;
307
309 gPad->SetSelected(view);
310 Int_t dist = big;
311 // Int_t id;
312
313 if (fPaintingOverlaps) {
317 dist = crt->GetShape()->DistancetoPrimitive(px, py);
318 if (dist < maxdist) {
319 gPad->SetSelected(crt);
320 box = (TGeoBBox *)crt->GetShape();
321 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
322 fCheckedBox[3] = box->GetDX();
323 fCheckedBox[4] = box->GetDY();
324 fCheckedBox[5] = box->GetDZ();
325 return 0;
326 }
329 dist = crt->GetShape()->DistancetoPrimitive(px, py);
330 if (dist < maxdist) {
331 gPad->SetSelected(crt);
332 box = (TGeoBBox *)crt->GetShape();
333 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
334 fCheckedBox[3] = box->GetDX();
335 fCheckedBox[4] = box->GetDY();
336 fCheckedBox[5] = box->GetDZ();
337 return 0;
338 }
339 return big;
340 }
341 // Compute distance to the right edge
342 if ((puxmax + inaxis - px) < 40) {
343 if ((py - puymax + inaxis) < 40) {
344 // when the mouse points to the (40x40) right corner of the pad, the manager class is selected
345 gPad->SetSelected(fGeoManager);
347 box = (TGeoBBox *)volume->GetShape();
348 memcpy(fCheckedBox, box->GetOrigin(), 3 * sizeof(Double_t));
349 fCheckedBox[3] = box->GetDX();
350 fCheckedBox[4] = box->GetDY();
351 fCheckedBox[5] = box->GetDZ();
352 return 0;
353 }
354 // when the mouse points to the (40 pix) right edge of the pad, the top volume is selected
355 gPad->SetSelected(volume);
356 fVolInfo = volume->GetName();
357 box = (TGeoBBox *)volume->GetShape();
358 memcpy(fCheckedBox, box->GetOrigin(), 3 * sizeof(Double_t));
359 fCheckedBox[3] = box->GetDX();
360 fCheckedBox[4] = box->GetDY();
361 fCheckedBox[5] = box->GetDZ();
362 return 0;
363 }
364
365 TGeoVolume *vol = volume;
366 Bool_t vis = vol->IsVisible();
367 // Bool_t drawDaughters = kTRUE;
368 // Do we need to check a branch only?
369 if (volume->IsVisBranch()) {
370 if (!fGeoManager->IsClosed())
371 return big;
374 while (fGeoManager->GetLevel()) {
377 dist = vol->GetShape()->DistancetoPrimitive(px, py);
378 if (dist < maxdist) {
380 box = (TGeoBBox *)vol->GetShape();
381 fGeoManager->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
384 gPad->SetSelected(fCheckedNode);
385 else
386 gPad->SetSelected(vol);
387 fCheckedBox[3] = box->GetDX();
388 fCheckedBox[4] = box->GetDY();
389 fCheckedBox[5] = box->GetDZ();
391 return 0;
392 }
393 fGeoManager->CdUp();
394 }
396 return dist;
397 }
398
399 // Do I need to look for the top volume ?
400 if ((fTopVisible && vis) || !vol->GetNdaughters() || !vol->IsVisDaughters() || vol->IsVisOnly()) {
401 dist = vol->GetShape()->DistancetoPrimitive(px, py);
402 if (dist < maxdist) {
403 fVolInfo = vol->GetName();
404 gPad->SetSelected(vol);
405 box = (TGeoBBox *)vol->GetShape();
406 memcpy(fCheckedBox, box->GetOrigin(), 3 * sizeof(Double_t));
407 fCheckedBox[3] = box->GetDX();
408 fCheckedBox[4] = box->GetDY();
409 fCheckedBox[5] = box->GetDZ();
410 return 0;
411 }
412 if (vol->IsVisOnly() || !vol->GetNdaughters() || !vol->IsVisDaughters())
413 return dist;
414 }
415
416 // Iterate the volume content
417 TGeoIterator next(vol);
418 next.SetTopName(TString::Format("%s_1", vol->GetName()));
420
421 Int_t level, nd;
422 Bool_t last;
423
424 while ((daughter = next())) {
425 vol = daughter->GetVolume();
426 level = next.GetLevel();
427 nd = daughter->GetNdaughters();
428 vis = daughter->IsVisible();
429 if (volume->IsVisContainers()) {
430 if (vis && level <= fVisLevel) {
431 *fGlobal = next.GetCurrentMatrix();
432 dist = vol->GetShape()->DistancetoPrimitive(px, py);
433 if (dist < maxdist) {
434 next.GetPath(fVolInfo);
435 box = (TGeoBBox *)vol->GetShape();
436 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
439 gPad->SetSelected(fCheckedNode);
440 else
441 gPad->SetSelected(vol);
442 fCheckedBox[3] = box->GetDX();
443 fCheckedBox[4] = box->GetDY();
444 fCheckedBox[5] = box->GetDZ();
446 return 0;
447 }
448 }
449 // Check if we have to skip this branch
450 if (level == fVisLevel || !daughter->IsVisDaughters()) {
451 next.Skip();
452 continue;
453 }
454 } else if (volume->IsVisLeaves()) {
455 last = ((nd == 0) || (level == fVisLevel) || (!daughter->IsVisDaughters())) ? kTRUE : kFALSE;
456 if (vis && last) {
457 *fGlobal = next.GetCurrentMatrix();
458 dist = vol->GetShape()->DistancetoPrimitive(px, py);
459 if (dist < maxdist) {
460 next.GetPath(fVolInfo);
461 box = (TGeoBBox *)vol->GetShape();
462 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
465 gPad->SetSelected(fCheckedNode);
466 else
467 gPad->SetSelected(vol);
468 fCheckedBox[3] = box->GetDX();
469 fCheckedBox[4] = box->GetDY();
470 fCheckedBox[5] = box->GetDZ();
472 return 0;
473 }
474 }
475 // Check if we have to skip the branch
476 if (last || !daughter->IsVisDaughters())
477 next.Skip();
478 }
479 }
480 return dist;
481}
482
483////////////////////////////////////////////////////////////////////////////////
484/// Set default angles for the current view.
485
487{
488 if (gPad) {
489 Int_t irep;
490 TView *view = gPad->GetView();
491 if (!view)
492 return;
493 view->SetView(-206, 126, 75, irep);
494 ModifiedPad();
495 }
496}
497
498////////////////////////////////////////////////////////////////////////////////
499/// Set default volume colors according to tracking media
500
502{
504 TGeoVolume *vol;
505 while ((vol = (TGeoVolume *)next()))
507 ModifiedPad();
508}
509
510////////////////////////////////////////////////////////////////////////////////
511/// Count number of visible nodes down to a given level.
512
514{
515 TGeoVolume *vol = volume;
516 Int_t count = 0;
517 Bool_t vis = vol->IsVisible();
518 // Do I need to look for the top volume ?
519 if ((fTopVisible && vis) || !vol->GetNdaughters() || !vol->IsVisDaughters() || vol->IsVisOnly())
520 count++;
521 // Is this the only volume?
522 if (volume->IsVisOnly())
523 return count;
524
525 // Do we need to check a branch only?
526 if (volume->IsVisBranch()) {
529 count = fGeoManager->GetLevel() + 1;
531 return count;
532 }
533 // Iterate the volume content
534 TGeoIterator next(vol);
536 Int_t level, nd;
537 Bool_t last;
538
539 while ((daughter = next())) {
540 // vol = daughter->GetVolume();
541 level = next.GetLevel();
542 nd = daughter->GetNdaughters();
543 vis = daughter->IsVisible();
544 if (volume->IsVisContainers()) {
545 if (vis && level <= rlevel)
546 count++;
547 // Check if we have to skip this branch
548 if (level == rlevel || !daughter->IsVisDaughters()) {
549 next.Skip();
550 continue;
551 }
552 } else if (volume->IsVisLeaves()) {
553 last = ((nd == 0) || (level == rlevel) || (!daughter->IsVisDaughters())) ? kTRUE : kFALSE;
554 if (vis && last)
555 count++;
556 // Check if we have to skip the branch
557 if (last)
558 next.Skip();
559 }
560 }
561 return count;
562}
563
564////////////////////////////////////////////////////////////////////////////////
565/// Count total number of visible nodes.
566
568{
570 Int_t vislevel = fGeoManager->GetVisLevel();
571 // TGeoVolume *top = fGeoManager->GetTopVolume();
572 TGeoVolume *top = fTopVolume;
573 if (maxnodes <= 0 && top) {
574 fNVisNodes = CountNodes(top, vislevel);
575 SetVisLevel(vislevel);
576 return fNVisNodes;
577 }
578 // if (the total number of nodes of the top volume is less than maxnodes
579 // we can visualize everything.
580 // recompute the best visibility level
581 if (!top) {
582 SetVisLevel(vislevel);
583 return 0;
584 }
585 fNVisNodes = -1;
587 for (Int_t level = 1; level < 20; level++) {
588 vislevel = level;
589 Int_t nnodes = CountNodes(top, level);
590 if (top->IsVisOnly() || top->IsVisBranch()) {
591 vislevel = fVisLevel;
593 break;
594 }
595 if (nnodes > maxnodes) {
596 vislevel--;
597 break;
598 }
599 if (nnodes == fNVisNodes) {
600 if (again)
601 break;
602 again = kTRUE;
603 }
605 }
606 SetVisLevel(vislevel);
607 return fNVisNodes;
608}
609
610////////////////////////////////////////////////////////////////////////////////
611/// Check if Ged library is loaded and load geometry editor classe.
612
614{
615 if (fIsEditable)
616 return;
617 if (!TClass::GetClass("TGedEditor"))
618 return;
620 if ((h = gROOT->GetPluginManager()->FindHandler("TGeoManagerEditor"))) {
621 if (h->LoadPlugin() == -1)
622 return;
623 h->ExecPlugin(0);
624 }
626}
627
628////////////////////////////////////////////////////////////////////////////////
629/// Start the geometry editor.
630
632{
633 if (!gPad)
634 return;
635 if (!fIsEditable) {
636 if (!option[0])
637 gPad->GetCanvas()->GetCanvasImp()->ShowEditor();
638 else
640 CheckEdit();
641 }
642 gPad->SetSelected(fGeoManager);
643 gPad->GetCanvas()->Selected(gPad, fGeoManager, kButton1Down);
644}
645
646////////////////////////////////////////////////////////////////////////////////
647/// Draw method.
648
653
654////////////////////////////////////////////////////////////////////////////////
655/// Draw the time evolution of a radionuclide.
656
658{
659 Int_t ncoeff = sol->GetNcoeff();
660 if (!ncoeff)
661 return;
662 Double_t tlo = 0., thi = 0.;
663 Double_t cn = 0., lambda = 0.;
664 Int_t i;
665 sol->GetRange(tlo, thi);
666 Bool_t autorange = (thi == 0.) ? kTRUE : kFALSE;
667
668 // Try to find the optimum range in time.
669 if (autorange)
670 tlo = 0.;
671 sol->GetCoeff(0, cn, lambda);
672 Double_t lambdamin = lambda;
673 TString formula = "";
674 for (i = 0; i < ncoeff; i++) {
675 sol->GetCoeff(i, cn, lambda);
676 formula += TString::Format("%g*exp(-%g*x)", cn, lambda);
677 if (i < ncoeff - 1)
678 formula += "+";
680 lambdamin = lambda;
681 }
682 if (autorange)
683 thi = 10. / lambdamin;
684 // Create a function
685 TF1 *func = new TF1(TString::Format("conc%s", sol->GetElement()->GetName()), formula.Data(), tlo, thi);
686 func->SetTitle(formula + ";time[s]" + TString::Format(";Concentration_of_%s", sol->GetElement()->GetName()));
687 func->SetMinimum(1.e-3);
688 func->SetMaximum(1.25 * TMath::Max(sol->Concentration(tlo), sol->Concentration(thi)));
689 func->SetLineColor(sol->GetLineColor());
690 func->SetLineStyle(sol->GetLineStyle());
691 func->SetLineWidth(sol->GetLineWidth());
692 func->SetMarkerColor(sol->GetMarkerColor());
693 func->SetMarkerStyle(sol->GetMarkerStyle());
694 func->SetMarkerSize(sol->GetMarkerSize());
695 func->Draw(option);
696}
697
698////////////////////////////////////////////////////////////////////////////////
699/// Draw a polygon in 3D.
700
702{
703 Int_t nvert = poly->GetNvert();
704 if (!nvert) {
705 Error("DrawPolygon", "No vertices defined");
706 return;
707 }
708 Int_t nconv = poly->GetNconvex();
709 Double_t *x = new Double_t[nvert + 1];
710 Double_t *y = new Double_t[nvert + 1];
711 poly->GetVertices(x, y);
712 x[nvert] = x[0];
713 y[nvert] = y[0];
714 TGraph *g1 = new TGraph(nvert + 1, x, y);
715 g1->SetTitle(Form("Polygon with %d vertices (outscribed %d)", nvert, nconv));
716 g1->SetLineColor(kRed);
717 g1->SetMarkerColor(kRed);
718 g1->SetMarkerStyle(4);
719 g1->SetMarkerSize(0.8);
720 delete[] x;
721 delete[] y;
722 Double_t *xc = nullptr;
723 Double_t *yc = nullptr;
724 TGraph *g2 = nullptr;
725 if (nconv && !poly->IsConvex()) {
726 xc = new Double_t[nconv + 1];
727 yc = new Double_t[nconv + 1];
728 poly->GetConvexVertices(xc, yc);
729 xc[nconv] = xc[0];
730 yc[nconv] = yc[0];
731 g2 = new TGraph(nconv + 1, xc, yc);
732 g2->SetLineColor(kBlue);
733 g2->SetLineColor(kBlue);
734 g2->SetMarkerColor(kBlue);
735 g2->SetMarkerStyle(21);
736 g2->SetMarkerSize(0.4);
737 delete[] xc;
738 delete[] yc;
739 }
740 if (!gPad) {
741 gROOT->MakeDefCanvas();
742 }
743 g1->Draw("ALP");
744 if (g2)
745 g2->Draw("LP");
746}
747
748////////////////////////////////////////////////////////////////////////////////
749/// Draw method.
750
752{
753 fTopVolume = vol;
754 fLastVolume = nullptr;
756 // if (fVisOption==kGeoVisOnly ||
757 // fVisOption==kGeoVisBranch) fGeoManager->SetVisOption(kGeoVisLeaves);
759 TString opt = option;
760 opt.ToLower();
762 fOverlap = nullptr;
763
764 if (fVisLock) {
767 }
768 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
769 // Clear pad if option "same" not given
770 if (!gPad) {
771 gROOT->MakeDefCanvas();
772 }
773 if (!opt.Contains("same"))
774 gPad->Clear();
775 // append this volume to pad
777
778 // Create a 3-D view
779 TView *view = gPad->GetView();
780 if (!view) {
781 view = TView::CreateView(11, nullptr, nullptr);
782 // Set the view to perform a first autorange (frame) draw.
783 // TViewer3DPad will revert view to normal painting after this
784 view->SetAutoRange(kTRUE);
785 if (has_pad)
786 gPad->Update();
787 }
788 if (!opt.Contains("same"))
789 Paint("range");
790 else
791 Paint(opt);
792 view->SetAutoRange(kFALSE);
793 // If we are drawing into the pad, then the view needs to be
794 // set to perspective
795 // if (!view->IsPerspective()) view->SetPerspective();
796
798
799 // Create a 3D viewer to paint us
800 gPad->GetViewer3D(option);
801}
802
803////////////////////////////////////////////////////////////////////////////////
804/// Draw a shape.
805
807{
808 TString opt = option;
809 opt.ToLower();
811 fOverlap = nullptr;
813
814 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
815 // Clear pad if option "same" not given
816 if (!gPad) {
817 gROOT->MakeDefCanvas();
818 }
819 if (!opt.Contains("same"))
820 gPad->Clear();
821 // append this shape to pad
822 shape->AppendPad(option);
823
824 // Create a 3-D view
825 TView *view = gPad->GetView();
826 if (!view) {
827 view = TView::CreateView(11, nullptr, nullptr);
828 // Set the view to perform a first autorange (frame) draw.
829 // TViewer3DPad will revert view to normal painting after this
830 view->SetAutoRange(kTRUE);
831 if (has_pad)
832 gPad->Update();
833 }
834 PaintShape(shape, "range");
835 view->SetAutoRange(kFALSE);
836 view->SetPerspective();
837 // Create a 3D viewer to paint us
838 gPad->GetViewer3D(option);
839}
840
841////////////////////////////////////////////////////////////////////////////////
842/// Draw an overlap.
843
845{
846 TString opt = option;
849 if (!overlap)
850 return;
851
854 opt.ToLower();
855 if (fVisLock) {
858 }
859 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
860 // Clear pad if option "same" not given
861 if (!gPad) {
862 gROOT->MakeDefCanvas();
863 }
864 if (!opt.Contains("same"))
865 gPad->Clear();
866 // append this volume to pad
867 overlap->AppendPad(option);
868
869 // Create a 3-D view
870 // Create a 3D viewer to paint us
871 gPad->GetViewer3D(option);
872 TView *view = gPad->GetView();
873 if (!view) {
874 view = TView::CreateView(11, nullptr, nullptr);
875 // Set the view to perform a first autorange (frame) draw.
876 // TViewer3DPad will revert view to normal painting after this
877 view->SetAutoRange(kTRUE);
878 PaintOverlap(ovlp, "range");
879 overlap->GetPolyMarker()->Draw("SAME");
880 if (has_pad)
881 gPad->Update();
882 }
883
884 // If we are drawing into the pad, then the view needs to be
885 // set to perspective
886 // if (!view->IsPerspective()) view->SetPerspective();
887 fVisLock = kTRUE;
888}
889
890////////////////////////////////////////////////////////////////////////////////
891/// Draw only one volume.
892
894{
895 TString opt = option;
896 opt.ToLower();
897 if (fVisLock) {
900 }
903 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
904 // Clear pad if option "same" not given
905 if (!gPad) {
906 gROOT->MakeDefCanvas();
907 }
908 if (!opt.Contains("same"))
909 gPad->Clear();
910 // append this volume to pad
913
914 // Create a 3-D view
915 TView *view = gPad->GetView();
916 if (!view) {
917 view = TView::CreateView(11, nullptr, nullptr);
918 // Set the view to perform a first autorange (frame) draw.
919 // TViewer3DPad will revert view to normal painting after this
920 view->SetAutoRange(kTRUE);
922 if (has_pad)
923 gPad->Update();
924 }
925
926 // If we are drawing into the pad, then the view needs to be
927 // set to perspective
928 // if (!view->IsPerspective()) view->SetPerspective();
929 fVisLock = kTRUE;
930}
931
932////////////////////////////////////////////////////////////////////////////////
933/// Draw current point in the same view.
934
936{
937 if (!gPad)
938 return;
939 if (!gPad->GetView())
940 return;
942 pm->SetMarkerColor(color);
943 const Double_t *point = fGeoManager->GetCurrentPoint();
944 pm->SetNextPoint(point[0], point[1], point[2]);
945 pm->SetMarkerStyle(8);
946 pm->SetMarkerSize(0.5);
947 pm->Draw("SAME");
948}
949
950////////////////////////////////////////////////////////////////////////////////
951
953
954////////////////////////////////////////////////////////////////////////////////
955/// Draw all volumes for a given path.
956
966
967////////////////////////////////////////////////////////////////////////////////
968/// Estimate camera movement between tmin and tmax for best track display
969
971{
972 if (!gPad)
973 return;
974 TIter next(gPad->GetListOfPrimitives());
976 TObject *obj;
977 Int_t ntracks = 0;
978 Double_t *point = nullptr;
979 AddTrackPoint(point, start, kTRUE);
980 while ((obj = next())) {
981 if (strcmp(obj->ClassName(), "TGeoTrack"))
982 continue;
983 track = (TVirtualGeoTrack *)obj;
984 ntracks++;
985 track->PaintCollect(tmin, start);
986 }
987
988 if (!ntracks)
989 return;
990 next.Reset();
991 AddTrackPoint(point, end, kTRUE);
992 while ((obj = next())) {
993 if (strcmp(obj->ClassName(), "TGeoTrack"))
994 continue;
995 track = (TVirtualGeoTrack *)obj;
996 if (!track)
997 continue;
998 track->PaintCollect(tmax, end);
999 }
1000}
1001
1002////////////////////////////////////////////////////////////////////////////////
1003/// Execute mouse actions on a given volume.
1004
1005void TGeoPainter::ExecuteManagerEvent(TGeoManager * /*geom*/, Int_t event, Int_t /*px*/, Int_t /*py*/)
1006{
1007 if (!gPad)
1008 return;
1009 gPad->SetCursor(kPointer);
1010 switch (event) {
1011 case kButton1Down:
1012 if (!fIsEditable)
1013 CheckEdit();
1014 }
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018/// Execute mouse actions on a given shape.
1019
1020void TGeoPainter::ExecuteShapeEvent(TGeoShape * /*shape*/, Int_t event, Int_t /*px*/, Int_t /*py*/)
1021{
1022 if (!gPad)
1023 return;
1024 gPad->SetCursor(kHand);
1025 switch (event) {
1026 case kButton1Down:
1027 if (!fIsEditable)
1028 CheckEdit();
1029 }
1030}
1031
1032////////////////////////////////////////////////////////////////////////////////
1033/// Execute mouse actions on a given volume.
1034
1035void TGeoPainter::ExecuteVolumeEvent(TGeoVolume * /*volume*/, Int_t event, Int_t /*px*/, Int_t /*py*/)
1036{
1037 if (!gPad)
1038 return;
1039 if (!fIsEditable)
1040 CheckEdit();
1041 // if (fIsRaytracing) return;
1042 // Bool_t istop = (volume==fTopVolume)?kTRUE:kFALSE;
1043 // if (istop) gPad->SetCursor(kHand);
1044 // else gPad->SetCursor(kPointer);
1045 gPad->SetCursor(kHand);
1046 // static Int_t width, color;
1047 switch (event) {
1048 case kMouseEnter:
1049 // width = volume->GetLineWidth();
1050 // color = volume->GetLineColor();
1051 break;
1052
1053 case kMouseLeave:
1054 // volume->SetLineWidth(width);
1055 // volume->SetLineColor(color);
1056 break;
1057
1058 case kButton1Down:
1059 // volume->SetLineWidth(3);
1060 // volume->SetLineColor(2);
1061 // gPad->Modified();
1062 // gPad->Update();
1063 break;
1064
1065 case kButton1Up:
1066 // volume->SetLineWidth(width);
1067 // volume->SetLineColor(color);
1068 // gPad->Modified();
1069 // gPad->Update();
1070 break;
1071
1072 case kButton1Double:
1073 gPad->SetCursor(kWatch);
1074 GrabFocus();
1075 break;
1076 }
1077}
1078
1079////////////////////////////////////////////////////////////////////////////////
1080/// Get some info about the current selected volume.
1081
1082const char *TGeoPainter::GetVolumeInfo(const TGeoVolume *volume, Int_t /*px*/, Int_t /*py*/) const
1083{
1084 static TString info;
1085 info = "";
1086 if (!gPad)
1087 return info;
1088 if (fPaintingOverlaps) {
1089 if (!fOverlap) {
1090 info = "wrong overlapping flag";
1091 return info;
1092 }
1094 if (fOverlap->IsExtrusion())
1095 ovtype = "EXTRUSION";
1096 else
1097 ovtype = "OVERLAP";
1098 if (volume == fOverlap->GetFirstVolume())
1099 name = volume->GetName();
1100 else
1102 info = TString::Format("%s: %s of %g", name.Data(), ovtype.Data(), fOverlap->GetOverlap());
1103 return info;
1104 } else
1105 info = TString::Format("%s, shape=%s", fVolInfo.Data(), volume->GetShape()->ClassName());
1106 return info;
1107}
1108
1109////////////////////////////////////////////////////////////////////////////////
1110/// Get the current view angles.
1111
1113{
1114 if (!gPad)
1115 return;
1116 TView *view = gPad->GetView();
1117 if (!view)
1118 return;
1119 longitude = view->GetLongitude();
1120 latitude = view->GetLatitude();
1121 psi = view->GetPsi();
1122}
1123
1124////////////////////////////////////////////////////////////////////////////////
1125/// Move focus to current volume
1126
1128{
1129 if (!gPad)
1130 return;
1131 TView *view = gPad->GetView();
1132 if (!view)
1133 return;
1135 printf("Woops!!!\n");
1137 memcpy(&fCheckedBox[0], box->GetOrigin(), 3 * sizeof(Double_t));
1138 fCheckedBox[3] = box->GetDX();
1139 fCheckedBox[4] = box->GetDY();
1140 fCheckedBox[5] = box->GetDZ();
1141 }
1142 view->SetPerspective();
1144 Int_t nframes = nfr;
1145 if (nfr == 0) {
1146 nframes = 1;
1147 if (nvols < 1500)
1148 nframes = 10;
1149 if (nvols < 1000)
1150 nframes = 20;
1151 if (nvols < 200)
1152 nframes = 50;
1153 if (nvols < 100)
1154 nframes = 100;
1155 }
1157}
1158
1159////////////////////////////////////////////////////////////////////////////////
1160/// Convert a local vector according view rotation matrix
1161
1163{
1164 for (Int_t i = 0; i < 3; i++)
1165 master[i] = -local[0] * fMat[i] - local[1] * fMat[i + 3] - local[2] * fMat[i + 6];
1166}
1167
1168////////////////////////////////////////////////////////////////////////////////
1169/// Check if a pad and view are present and send signal "Modified" to pad.
1170
1172{
1173 if (!gPad)
1174 return;
1175 if (update) {
1176 gPad->Update();
1177 return;
1178 }
1179 TView *view = gPad->GetView();
1180 if (!view)
1181 return;
1182 view->SetViewChanged();
1183 gPad->Modified();
1184 if (gROOT->FromPopUp())
1185 gPad->Update();
1186}
1187
1188////////////////////////////////////////////////////////////////////////////////
1189/// Paint current geometry according to option.
1190
1192{
1193 if (!fGeoManager || !fTopVolume)
1194 return;
1196 if (gPad)
1197 is_padviewer = (!strcmp(gPad->GetViewer3D()->ClassName(), "TViewer3DPad")) ? kTRUE : kFALSE;
1198
1202 else if (fTopVolume->IsVisLeaves())
1204 else if (fTopVolume->IsVisOnly())
1206 else if (fTopVolume->IsVisBranch())
1208
1209 if (!fIsRaytracing || !is_padviewer) {
1210 if (fGeoManager->IsDrawingExtra()) {
1211 // loop the list of physical volumes
1212 fGeoManager->CdTop();
1214 Int_t nnodes = nodeList->GetEntriesFast();
1215 Int_t inode;
1216 TGeoPhysicalNode *node;
1217 for (inode = 0; inode < nnodes; inode++) {
1218 node = (TGeoPhysicalNode *)nodeList->UncheckedAt(inode);
1220 }
1221 } else {
1223 }
1224 fVisLock = kTRUE;
1225 }
1226 // Check if we have to raytrace (only in pad)
1228 Raytrace();
1229}
1230
1231////////////////////////////////////////////////////////////////////////////////
1232/// Paint an overlap.
1233
1235{
1236 if (!fGeoManager)
1237 return;
1239 if (!overlap)
1240 return;
1241 Int_t color, transparency;
1242 if (fOverlap != overlap)
1243 fOverlap = overlap;
1246 TGeoVolume *vol;
1247 TGeoVolume *vol1 = overlap->GetFirstVolume();
1248 TGeoVolume *vol2 = overlap->GetSecondVolume();
1249 TGeoHMatrix *matrix1 = overlap->GetFirstMatrix();
1250 TGeoHMatrix *matrix2 = overlap->GetSecondMatrix();
1251 //
1252 vol = vol1;
1253 *hmat = matrix1;
1254 fGeoManager->SetMatrixReflection(matrix1->IsReflection());
1255 if (!fVisLock)
1256 fVisVolumes->Add(vol);
1258 color = vol->GetLineColor();
1260 vol->SetLineColor(kGreen);
1261 vol->SetTransparency(40);
1262 if (!strstr(option, "range"))
1263 ((TAttLine *)vol)->Modify();
1264 PaintShape(*(vol->GetShape()), option);
1265 vol->SetLineColor(color);
1267 vol = vol2;
1268 *hmat = matrix2;
1269 fGeoManager->SetMatrixReflection(matrix2->IsReflection());
1270 if (!fVisLock)
1271 fVisVolumes->Add(vol);
1273 color = vol->GetLineColor();
1275 vol->SetLineColor(kBlue);
1276 vol->SetTransparency(40);
1277 if (!strstr(option, "range"))
1278 ((TAttLine *)vol)->Modify();
1279 PaintShape(*(vol->GetShape()), option);
1280 vol->SetLineColor(color);
1283 fVisLock = kTRUE;
1284}
1285
1286////////////////////////////////////////////////////////////////////////////////
1287/// Paint recursively a node and its content according to visualization options.
1288
1293
1294////////////////////////////////////////////////////////////////////////////////
1295/// Paint recursively a node and its content according to visualization options.
1296
1298{
1299 if (fTopVolume != top) {
1301 fVisLock = kFALSE;
1302 }
1303 fTopVolume = top;
1304 if (!fVisLevel)
1305 return;
1306 TGeoVolume *vol = top;
1307 if (global)
1308 *fGlobal = *global;
1309 else
1310 fGlobal->Clear();
1313 Bool_t vis = (top->IsVisible() && !top->IsAssembly());
1314 Int_t transparency = 0;
1315
1316 // Update pad attributes in case we need to paint VOL
1317 if (!strstr(option, "range"))
1318 ((TAttLine *)vol)->Modify();
1319
1320 // Do we need to draw a branch ?
1321 if (top->IsVisBranch()) {
1324 // while (fGeoManager->GetLevel()) {
1326 if (!fVisLock) {
1327 fVisVolumes->Add(vol);
1329 }
1332 vol->SetTransparency(40);
1333 if (!strstr(option, "range"))
1334 ((TAttLine *)vol)->Modify();
1335 if (global) {
1336 *fGlobal = *global;
1338 } else {
1340 }
1342 PaintShape(*(vol->GetShape()), option);
1344 fGeoManager->CdUp();
1345 // }
1346 fVisLock = kTRUE;
1349 return;
1350 }
1351
1352 // Do I need to draw the top volume ?
1353 if ((fTopVisible && vis) || !top->GetNdaughters() || !top->IsVisDaughters() || top->IsVisOnly()) {
1356 PaintShape(*(vol->GetShape()), option);
1357 if (!fVisLock && !vol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1358 fVisVolumes->Add(vol);
1360 }
1361 if (top->IsVisOnly() || !top->GetNdaughters() || !top->IsVisDaughters()) {
1362 fVisLock = kTRUE;
1363 return;
1364 }
1365 }
1366
1367 // Iterate the volume content
1368 TGeoIterator next(vol);
1369 if (fPlugin)
1370 next.SetUserPlugin(fPlugin);
1372 // TGeoMatrix *glmat;
1373 Int_t level, nd;
1374 Bool_t last;
1375 Int_t line_color = 0, line_width = 0, line_style = 0;
1376 while ((daughter = next())) {
1377 vol = daughter->GetVolume();
1379 level = next.GetLevel();
1380 nd = daughter->GetNdaughters();
1381 vis = daughter->IsVisible();
1383 if (top->IsVisContainers()) {
1384 if (vis && level <= fVisLevel) {
1385 if (fPlugin) {
1386 line_color = vol->GetLineColor();
1387 line_width = vol->GetLineWidth();
1388 line_style = vol->GetLineStyle();
1391 }
1392 if (!strstr(option, "range"))
1393 ((TAttLine *)vol)->Modify();
1394 if (global) {
1395 *fGlobal = *global;
1396 *fGlobal *= *next.GetCurrentMatrix();
1397 } else {
1398 *fGlobal = next.GetCurrentMatrix();
1399 }
1402 if (fPlugin) {
1407 }
1408 if (!fVisLock && !daughter->IsOnScreen()) {
1409 fVisVolumes->Add(vol);
1411 }
1412 }
1413 // Check if we have to skip this branch
1414 if (!drawDaughters || level == fVisLevel || !daughter->IsVisDaughters()) {
1415 next.Skip();
1416 continue;
1417 }
1418 } else if (top->IsVisLeaves()) {
1419 last = ((nd == 0) || (level == fVisLevel) || (!daughter->IsVisDaughters())) ? kTRUE : kFALSE;
1420 if (vis && last) {
1421 if (fPlugin) {
1422 line_color = vol->GetLineColor();
1423 line_width = vol->GetLineWidth();
1424 line_style = vol->GetLineStyle();
1427 }
1428 if (!strstr(option, "range"))
1429 ((TAttLine *)vol)->Modify();
1430 if (global) {
1431 *fGlobal = *global;
1432 *fGlobal *= *next.GetCurrentMatrix();
1433 } else {
1434 *fGlobal = next.GetCurrentMatrix();
1435 }
1438 if (fPlugin) {
1443 }
1444 if (!fVisLock && !daughter->IsOnScreen()) {
1445 fVisVolumes->Add(vol);
1447 }
1448 }
1449 // Check if we have to skip the branch
1450 if (!drawDaughters || last || !daughter->IsVisDaughters())
1451 next.Skip();
1452 }
1453 }
1454 if (fPlugin)
1455 fPlugin->SetIterator(nullptr);
1457 fVisLock = kTRUE;
1458}
1459
1460////////////////////////////////////////////////////////////////////////////////
1461/// Paint the supplied shape into the current 3D viewer
1462
1464{
1466
1467 TVirtualViewer3D *viewer = gPad->GetViewer3D();
1468
1469 if (!viewer || shape.IsA() == TGeoShapeAssembly::Class()) {
1470 return addDaughters;
1471 }
1472
1473 // For non-composite shapes we are the main paint method & perform the negotiation
1474 // with the viewer here
1475 if (!shape.IsComposite()) {
1476 // Does viewer prefer local frame positions?
1477 Bool_t localFrame = viewer->PreferLocalFrame();
1478 // Perform first fetch of buffer from the shape and try adding it
1479 // to the viewer
1480 const TBuffer3D &buffer =
1482 Int_t reqSections = viewer->AddObject(buffer, &addDaughters);
1483
1484 // If the viewer requires additional sections fetch from the shape (if possible)
1485 // and add again
1488 viewer->AddObject(buffer, &addDaughters);
1489 }
1490 }
1491 // Composite shapes have their own internal hierarchy of shapes, each
1492 // of which generate a filled TBuffer3D. Therefore we can't pass up a
1493 // single buffer to here. So as a special case the TGeoCompositeShape
1494 // performs it's own painting & negotiation with the viewer.
1495 else {
1496 const TGeoCompositeShape *composite = static_cast<const TGeoCompositeShape *>(&shape);
1497
1498 // We need the addDaughters flag returned from the viewer from paint
1499 // so can't use the normal TObject::Paint()
1500 // TGeoHMatrix *matrix = (TGeoHMatrix*)TGeoShape::GetTransform();
1501 // if (viewer->PreferLocalFrame()) matrix->Clear();
1502 addDaughters = composite->PaintComposite(option);
1503 }
1504
1505 return addDaughters;
1506}
1507
1508////////////////////////////////////////////////////////////////////////////////
1509/// Paint an overlap.
1510
1518
1519////////////////////////////////////////////////////////////////////////////////
1520/// Paints a physical node associated with a path.
1521
1523{
1524 if (!node->IsVisible())
1525 return;
1526 Int_t level = node->GetLevel();
1527 Int_t i, col, wid, sty;
1528 TGeoShape *shape;
1532 if (!node->IsVisibleFull()) {
1533 // Paint only last node in the branch
1534 vcrt = node->GetVolume();
1535 if (!strstr(option, "range"))
1536 ((TAttLine *)vcrt)->Modify();
1537 shape = vcrt->GetShape();
1538 *matrix = node->GetMatrix();
1539 fGeoManager->SetMatrixReflection(matrix->IsReflection());
1541 if (!node->IsVolAttributes() && !strstr(option, "range")) {
1542 col = vcrt->GetLineColor();
1543 wid = vcrt->GetLineWidth();
1544 sty = vcrt->GetLineStyle();
1545 vcrt->SetLineColor(node->GetLineColor());
1546 vcrt->SetLineWidth(node->GetLineWidth());
1547 vcrt->SetLineStyle(node->GetLineStyle());
1548 ((TAttLine *)vcrt)->Modify();
1549 PaintShape(*shape, option);
1550 vcrt->SetLineColor(col);
1551 vcrt->SetLineWidth(wid);
1552 vcrt->SetLineStyle(sty);
1553 } else {
1554 PaintShape(*shape, option);
1555 }
1556 } else {
1557 // Paint full branch, except top node
1558 for (i = 1; i <= level; i++) {
1559 vcrt = node->GetVolume(i);
1560 if (!strstr(option, "range"))
1561 ((TAttLine *)vcrt)->Modify();
1562 shape = vcrt->GetShape();
1563 *matrix = node->GetMatrix(i);
1564 fGeoManager->SetMatrixReflection(matrix->IsReflection());
1566 if (!node->IsVolAttributes() && !strstr(option, "range")) {
1567 col = vcrt->GetLineColor();
1568 wid = vcrt->GetLineWidth();
1569 sty = vcrt->GetLineStyle();
1570 vcrt->SetLineColor(node->GetLineColor());
1571 vcrt->SetLineWidth(node->GetLineWidth());
1572 vcrt->SetLineStyle(node->GetLineStyle());
1573 ((TAttLine *)vcrt)->Modify();
1574 PaintShape(*shape, option);
1575 vcrt->SetLineColor(col);
1576 vcrt->SetLineWidth(wid);
1577 vcrt->SetLineStyle(sty);
1578 } else {
1579 PaintShape(*shape, option);
1580 }
1581 }
1582 }
1584}
1585
1586////////////////////////////////////////////////////////////////////////////////
1587/// Raytrace current drawn geometry
1588
1590{
1591 if (!gPad || gPad->IsBatch())
1592 return;
1593 TView *view = gPad->GetView();
1594 if (!view)
1595 return;
1598 if (top != fTopVolume)
1600 if (!view->IsPerspective())
1601 view->SetPerspective();
1602 gVirtualX->SetMarkerSize(1);
1603 gVirtualX->SetMarkerStyle(1);
1606 Double_t lat = view->GetLatitude();
1607 Double_t longit = view->GetLongitude();
1608 Double_t psi = view->GetPsi();
1615 fMat[0] = c1 * c3 - s1 * c2 * s3;
1616 fMat[1] = c1 * s3 + s1 * c2 * c3;
1617 fMat[2] = s1 * s2;
1618
1619 fMat[3] = -s1 * c3 - c1 * c2 * s3;
1620 fMat[4] = -s1 * s3 + c1 * c2 * c3;
1621 fMat[5] = c1 * s2;
1622
1623 fMat[6] = s2 * s3;
1624 fMat[7] = -s2 * c3;
1625 fMat[8] = c2;
1626 Double_t u0, v0, du, dv;
1627 view->GetWindow(u0, v0, du, dv);
1628 Double_t dview = view->GetDview();
1629 Double_t dproj = view->GetDproj();
1630 Double_t local[3] = {0, 0, 1};
1631 Double_t dir[3], normal[3];
1633 Double_t min[3], max[3];
1634 view->GetRange(min, max);
1635 Double_t cov[3];
1636 for (Int_t i = 0; i < 3; i++)
1637 cov[i] = 0.5 * (min[i] + max[i]);
1638 Double_t cop[3];
1639 for (Int_t i = 0; i < 3; i++)
1640 cop[i] = cov[i] - dir[i] * dview;
1641 fGeoManager->InitTrack(cop, dir);
1644 if (fClippingShape)
1646 Int_t px, py;
1649 pxmin = gPad->UtoAbsPixel(0);
1650 pxmax = gPad->UtoAbsPixel(1);
1651 pymin = gPad->VtoAbsPixel(1);
1652 pymax = gPad->VtoAbsPixel(0);
1653 TGeoNode *next = nullptr;
1654 TGeoNode *nextnode = nullptr;
1655 Double_t step, steptot;
1656 Double_t *norm;
1657 const Double_t *point = fGeoManager->GetCurrentPoint();
1658 Double_t *ppoint = (Double_t *)point;
1659 Double_t tosource[3];
1660 Double_t calf;
1661 Double_t phi = 45. * krad;
1662 tosource[0] = -dir[0] * TMath::Cos(phi) + dir[1] * TMath::Sin(phi);
1663 tosource[1] = -dir[0] * TMath::Sin(phi) - dir[1] * TMath::Cos(phi);
1664 tosource[2] = -dir[2];
1665
1667
1668 Bool_t done;
1669 // Int_t istep;
1670 Int_t base_color, color;
1673 TPoint *pxy = new TPoint[1];
1675 Int_t up;
1676 Int_t ntotal = pxmax * pymax;
1677 Int_t nrays = 0;
1678 TStopwatch *timer = new TStopwatch();
1679 timer->Start();
1680 for (px = pxmin; px < pxmax; px++) {
1681 for (py = pymin; py < pymax; py++) {
1682 if ((nrays % 100) == 0)
1683 checker->OpProgress("Raytracing", nrays, ntotal, timer, kFALSE);
1684 nrays++;
1685 base_color = 1;
1686 steptot = 0;
1688 xloc = gPad->AbsPixeltoX(pxmin + pxmax - px);
1689 xloc = xloc * du - u0;
1690 yloc = gPad->AbsPixeltoY(pymin + pymax - py);
1691 yloc = yloc * dv - v0;
1693 local[0] = xloc / modloc;
1694 local[1] = yloc / modloc;
1695 local[2] = dproj / modloc;
1701 // fGeoManager->InitTrack(cop,dir);
1702 // current ray pointing to pixel (px,py)
1703 done = kFALSE;
1704 norm = nullptr;
1705 // propagate to the clipping shape if any
1706 if (fClippingShape) {
1707 if (inclip) {
1710 } else {
1712 stemin = 0;
1713 }
1714 }
1715
1716 while (!done) {
1717 if (fClippingShape) {
1718 if (stemin > 1E10)
1719 break;
1720 if (stemin > 0) {
1721 // we are inside clipping shape
1723 next = fGeoManager->Step();
1724 steptot = 0;
1725 stemin = 0;
1726 up = 0;
1727 while (next) {
1728 // we found something after clipping region
1729 nextvol = next->GetVolume();
1730 if (nextvol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1731 done = kTRUE;
1732 base_color = nextvol->GetLineColor();
1734 norm = normal;
1735 break;
1736 }
1737 up++;
1738 next = fGeoManager->GetMother(up);
1739 }
1740 if (done)
1741 break;
1743 fGeoManager->SetStep(1E-3);
1744 while (inclip) {
1745 fGeoManager->Step();
1747 }
1749 }
1750 }
1752 step = fGeoManager->GetStep();
1753 if (step > 1E10)
1754 break;
1755 steptot += step;
1756 next = nextnode;
1757 // Check the step
1758 if (fClippingShape) {
1759 if (steptot > stemax) {
1760 steptot = 0;
1762 if (inclip) {
1765 continue;
1766 } else {
1767 stemin = 0;
1769 }
1770 }
1771 }
1772 // Check if next node is visible
1773 if (!nextnode)
1774 continue;
1775 nextvol = nextnode->GetVolume();
1776 if (nextvol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1777 done = kTRUE;
1778 base_color = nextvol->GetLineColor();
1779 next = nextnode;
1780 break;
1781 }
1782 }
1783 if (!done)
1784 continue;
1785 // current ray intersect a visible volume having color=base_color
1786 if (rtMode > 0) {
1789 for (Int_t i = 0; i < 3; ++i)
1790 local[i] += 1.E-8 * dir[i];
1791 step = next->GetVolume()->GetShape()->DistFromInside(local, dir, 3);
1792 for (Int_t i = 0; i < 3; ++i)
1793 local[i] += step * dir[i];
1794 next->GetVolume()->GetShape()->ComputeNormal(local, dir, normal);
1795 norm = normal;
1796 } else {
1797 if (!norm)
1799 if (!norm)
1800 continue;
1801 }
1802 calf = norm[0] * tosource[0] + norm[1] * tosource[1] + norm[2] * tosource[2];
1804 color = GetColor(base_color, light);
1805 // Now we know the color of the pixel, just draw it
1806 gVirtualX->SetMarkerColor(color);
1807 pxy[0].fX = px;
1808 pxy[0].fY = py;
1809 gVirtualX->DrawPolyMarker(1, pxy);
1810 }
1811 }
1812 delete[] pxy;
1813 timer->Stop();
1814 checker->OpProgress("Raytracing", nrays, ntotal, timer, kTRUE);
1815 delete timer;
1816}
1817
1818////////////////////////////////////////////////////////////////////////////////
1819/// Set cartesian and radial bomb factors for translations.
1820
1830
1831////////////////////////////////////////////////////////////////////////////////
1832/// Set type of exploding view.
1833
1835{
1836 if ((ibomb < 0) || (ibomb > 3)) {
1837 Warning("SetExplodedView", "exploded view can be 0-3");
1838 return;
1839 }
1840 if ((Int_t)ibomb == fExplodedView)
1841 return;
1842 Bool_t change = (gPad == nullptr) ? kFALSE : kTRUE;
1843
1844 if (ibomb == kGeoNoBomb) {
1846 }
1847 if (ibomb == kGeoBombXYZ) {
1849 }
1850 if (ibomb == kGeoBombCyl) {
1852 }
1853 if (ibomb == kGeoBombSph) {
1855 }
1857 if (change)
1858 ModifiedPad();
1859}
1860
1861////////////////////////////////////////////////////////////////////////////////
1862/// Set number of segments to approximate circles.
1863
1865{
1866 if (nseg < 3) {
1867 Warning("SetNsegments", "number of segments should be > 2");
1868 return;
1869 }
1870 if (fNsegments == nseg)
1871 return;
1872 fNsegments = nseg;
1873 ModifiedPad();
1874}
1875
1876////////////////////////////////////////////////////////////////////////////////
1877/// Set default level down to which visualization is performed
1878
1880{
1881 if (level == fVisLevel && fLastVolume == fTopVolume)
1882 return;
1883 fVisLevel = level;
1884 if (!fTopVolume)
1885 return;
1886 if (fVisLock) {
1888 fVisLock = kFALSE;
1889 }
1890 if (!fLastVolume) {
1891 // printf("--- Drawing %6d nodes with %d visible levels\n",fNVisNodes,fVisLevel);
1892 return;
1893 }
1894 if (!gPad)
1895 return;
1896 if (gPad->GetView()) {
1897 // printf("--- Drawing %6d nodes with %d visible levels\n",fNVisNodes,fVisLevel);
1898 ModifiedPad();
1899 }
1900}
1901
1902////////////////////////////////////////////////////////////////////////////////
1903/// Set top geometry volume as visible.
1904
1906{
1907 if (fTopVisible == vis)
1908 return;
1909 fTopVisible = vis;
1910 ModifiedPad();
1911}
1912
1913////////////////////////////////////////////////////////////////////////////////
1914/// Set drawing mode :
1915/// - option=0 (default) all nodes drawn down to vislevel
1916/// - option=1 leaves and nodes at vislevel drawn
1917/// - option=2 path is drawn
1918
1920{
1921 if ((fVisOption < 0) || (fVisOption > 4)) {
1922 Warning("SetVisOption", "wrong visualization option");
1923 return;
1924 }
1925
1926 if (option == kGeoVisChanged) {
1927 if (fVisLock) {
1929 fVisLock = kFALSE;
1930 }
1931 ModifiedPad();
1932 return;
1933 }
1934
1935 if (fTopVolume) {
1937 att->SetAttBit(TGeoAtt::kVisBranch, kFALSE);
1938 att->SetAttBit(TGeoAtt::kVisContainers, kFALSE);
1939 att->SetAttBit(TGeoAtt::kVisOnly, kFALSE);
1940 switch (option) {
1941 case kGeoVisDefault: att->SetAttBit(TGeoAtt::kVisContainers, kTRUE); break;
1942 case kGeoVisLeaves: break;
1943 case kGeoVisOnly: att->SetAttBit(TGeoAtt::kVisOnly, kTRUE); break;
1944 }
1945 }
1946
1947 if (fVisOption == option)
1948 return;
1950 if (fVisLock) {
1952 fVisLock = kFALSE;
1953 }
1954 ModifiedPad();
1955}
1956
1957////////////////////////////////////////////////////////////////////////////////
1958/// Returns distance between point px,py on the pad an a shape.
1959
1961{
1962 const Int_t inaxis = 7;
1963 const Int_t maxdist = 5;
1964 const Int_t big = 9999;
1965 Int_t dist = big;
1966 if (!gPad)
1967 return dist;
1968 TView *view = gPad->GetView();
1969 if (!(numpoints && view))
1970 return dist;
1971 if (shape->IsA() == TGeoShapeAssembly::Class())
1972 return dist;
1973
1974 if (fIsPaintingShape) {
1975 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
1976 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
1977 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
1978 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
1979 // return if point not in user area
1980 if (px < puxmin - inaxis)
1981 return big;
1982 if (py > puymin + inaxis)
1983 return big;
1984 if (px > puxmax + inaxis)
1985 return big;
1986 if (py < puymax - inaxis)
1987 return big;
1988 if ((puxmax + inaxis - px) < 40) {
1989 // when the mouse points to the (40 pix) right edge of the pad, the manager class is selected
1990 gPad->SetSelected(fGeoManager);
1991 return 0;
1992 }
1993 }
1994
1995 fBuffer->SetRawSizes(numpoints, 3 * numpoints, 0, 0, 0, 0);
1997 shape->SetPoints(points);
1998 Double_t dpoint2, x1, y1, xndc[3];
1999 Double_t dmaster[3];
2000 Int_t j;
2001 for (Int_t i = 0; i < numpoints; i++) {
2002 j = 3 * i;
2003 TGeoShape::GetTransform()->LocalToMaster(&points[j], dmaster);
2004 points[j] = dmaster[0];
2005 points[j + 1] = dmaster[1];
2006 points[j + 2] = dmaster[2];
2007 view->WCtoNDC(&points[j], xndc);
2008 x1 = gPad->XtoAbsPixel(xndc[0]);
2009 y1 = gPad->YtoAbsPixel(xndc[1]);
2010 dpoint2 = (px - x1) * (px - x1) + (py - y1) * (py - y1);
2011 if (dpoint2 < dist)
2012 dist = (Int_t)dpoint2;
2013 }
2014 if (dist > 100)
2015 return dist;
2016 dist = Int_t(TMath::Sqrt(Double_t(dist)));
2017 if (dist < maxdist && fIsPaintingShape)
2018 gPad->SetSelected((TObject *)shape);
2019 return dist;
2020}
2021
2022////////////////////////////////////////////////////////////////////////////////
2023/// Get the new 'unbombed' translation vector according current exploded view mode.
2024
2026{
2027 memcpy(bombtr, tr, 3 * sizeof(Double_t));
2028 switch (fExplodedView) {
2029 case kGeoNoBomb: return;
2030 case kGeoBombXYZ:
2031 bombtr[0] /= fBombX;
2032 bombtr[1] /= fBombY;
2033 bombtr[2] /= fBombZ;
2034 return;
2035 case kGeoBombCyl:
2036 bombtr[0] /= fBombR;
2037 bombtr[1] /= fBombR;
2038 bombtr[2] /= fBombZ;
2039 return;
2040 case kGeoBombSph:
2041 bombtr[0] /= fBombR;
2042 bombtr[1] /= fBombR;
2043 bombtr[2] /= fBombR;
2044 return;
2045 default: return;
2046 }
2047}
@ kButton1Double
Definition Buttons.h:24
@ kButton1Up
Definition Buttons.h:19
@ kMouseLeave
Definition Buttons.h:23
@ kButton1Down
Definition Buttons.h:17
@ kMouseEnter
Definition Buttons.h:23
@ kWatch
Definition GuiTypes.h:375
@ kHand
Definition GuiTypes.h:374
@ kPointer
Definition GuiTypes.h:375
#define b(i)
Definition RSha256.hxx:100
#define g(i)
Definition RSha256.hxx:105
#define s1(x)
Definition RSha256.hxx:91
#define h(i)
Definition RSha256.hxx:106
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
@ kRed
Definition Rtypes.h:67
@ kBlack
Definition Rtypes.h:66
@ kGreen
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize wid
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x1
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 y1
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float xmax
#define gROOT
Definition TROOT.h:411
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
#define gPad
#define gVirtualX
Definition TVirtualX.h:337
#define gSize3D
Definition X3DBuffer.h:40
const_iterator end() const
Line Attributes class.
Definition TAttLine.h:20
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:37
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:36
virtual void Modify()
Change current line attributes if necessary.
Definition TAttLine.cxx:246
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:46
Generic 3D primitive description class.
Definition TBuffer3D.h:18
@ kBoundingBox
Definition TBuffer3D.h:51
@ kShapeSpecific
Definition TBuffer3D.h:52
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
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition TClass.cxx:2973
The color creation and management class.
Definition TColor.h:22
static void HLS2RGB(Float_t h, Float_t l, Float_t s, Float_t &r, Float_t &g, Float_t &b)
Static method to compute RGB from HLS.
Definition TColor.cxx:1581
static void InitializeColors()
Initialize colors used by the TCanvas based graphics (via TColor objects).
Definition TColor.cxx:1172
static Int_t CreateGradientColorTable(UInt_t Number, Double_t *Stops, Double_t *Red, Double_t *Green, Double_t *Blue, UInt_t NColors, Float_t alpha=1., Bool_t setPalette=kTRUE)
Static function creating a color table with several connected linear gradients.
Definition TColor.cxx:2742
static void RGB2HLS(Float_t r, Float_t g, Float_t b, Float_t &h, Float_t &l, Float_t &s)
Static method to compute HLS from RGB.
Definition TColor.cxx:1721
1-Dim function class
Definition TF1.h:182
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum value along Y for this function In case the function is already drawn,...
Definition TF1.cxx:3425
void SetTitle(const char *title="") override
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition TF1.cxx:3589
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1340
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum value along Y for this function In case the function is already drawn,...
Definition TF1.cxx:3438
Visualization and tracking attributes for volumes and nodes.
Definition TGeoAtt.h:17
Bool_t TestAttBit(UInt_t f) const
Definition TGeoAtt.h:64
Bool_t IsVisBranch() const
Definition TGeoAtt.h:85
@ kVisContainers
Definition TGeoAtt.h:32
@ kVisOnly
Definition TGeoAtt.h:33
@ kVisOnScreen
Definition TGeoAtt.h:31
@ kVisBranch
Definition TGeoAtt.h:34
void ResetAttBit(UInt_t f)
Definition TGeoAtt.h:63
void SetVisRaytrace(Bool_t flag=kTRUE)
Definition TGeoAtt.h:66
Bool_t IsVisDaughters() const
Definition TGeoAtt.h:84
void SetAttBit(UInt_t f)
Definition TGeoAtt.h:61
Box class.
Definition TGeoBBox.h:17
Class representing the Bateman solution for a decay branch.
Composite shapes are Boolean combinations of two or more shape components.
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:458
void Clear(Option_t *option="") override
clear the data for this matrix
void SetIterator(const TGeoIterator *iter)
Definition TGeoNode.h:237
virtual void ProcessNode()=0
A geometry iterator.
Definition TGeoNode.h:248
const TGeoMatrix * GetCurrentMatrix() const
Returns global matrix for current node.
void SetTopName(const char *name)
Set the top name for path.
Int_t GetLevel() const
Definition TGeoNode.h:294
void GetPath(TString &path) const
Returns the path for the current node.
void SetUserPlugin(TGeoIteratorPlugin *plugin)
Set a plugin.
void Skip()
Stop iterating the current branch.
The manager class for any TGeo geometry.
Definition TGeoManager.h:45
TGeoNode * GetMother(Int_t up=1) const
Double_t * FindNormalFast()
Computes fast normal to next crossed boundary, assuming that the current point is close enough to the...
void DoRestoreState()
Restore a backed-up state without affecting the cache stack.
const Double_t * GetCurrentDirection() const
TVirtualGeoChecker * GetGeomChecker()
Make a default checker if none present. Returns pointer to it.
void CdUp()
Go one level up in geometry.
void DoBackupState()
Backup the current state without affecting the cache stack.
TObjArray * GetListOfVolumes() const
void SetMatrixReflection(Bool_t flag=kTRUE)
virtual Bool_t cd(const char *path="")
Browse the tree of nodes starting from fTopNode according to pathname.
void LocalToMaster(const Double_t *local, Double_t *master) const
Int_t GetRTmode() const
Bool_t IsClosed() const
TGeoNode * GetCurrentNode() const
void SetCurrentDirection(Double_t *dir)
TGeoNode * Step(Bool_t is_geom=kTRUE, Bool_t cross=kTRUE)
Make a rectilinear step of length fStep from current point (fPoint) on current direction (fDirection)...
Bool_t IsDrawingExtra() const
void SetOutside(Bool_t flag=kTRUE)
TGeoNode * FindNextBoundaryAndStep(Double_t stepmax=TGeoShape::Big(), Bool_t compsafe=kFALSE)
Compute distance to next boundary within STEPMAX.
Int_t GetMaxVisNodes() const
void SetCurrentPoint(Double_t *point)
Int_t GetVisOption() const
Returns current depth to which geometry is drawn.
const Double_t * GetCurrentPoint() const
Bool_t IsOutside() const
TGeoNode * InitTrack(const Double_t *point, const Double_t *dir)
Initialize current point and current direction vector (normalized) in MARS.
void SetTopVolume(TGeoVolume *vol)
Set the top volume and corresponding node as starting point of the geometry.
Int_t GetLevel() const
Double_t GetStep() const
TGeoHMatrix * GetCurrentMatrix() const
TGeoNode * GetTopNode() const
void MasterToLocalVect(const Double_t *master, Double_t *local) const
void SetPaintVolume(TGeoVolume *vol)
void SetStep(Double_t step)
TGeoVolume * GetCurrentVolume() const
Int_t GetVisLevel() const
Returns current depth to which geometry is drawn.
Int_t GetBombMode() const
void CdTop()
Make top level node the current node.
void MasterToLocal(const Double_t *master, Double_t *local) const
Int_t PushPath(Int_t startlevel=0)
Int_t GetNsegments() const
Get number of segments approximating circles.
Bool_t IsNodeSelectable() const
TGeoVolume * GetTopVolume() const
TObjArray * GetListOfPhysicalNodes()
Bool_t PopPath()
virtual Int_t GetDefaultColor() const
Get some default color related to this material.
Geometrical transformation package.
Definition TGeoMatrix.h:38
Bool_t IsReflection() const
Definition TGeoMatrix.h:66
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:39
TGeoVolume * GetVolume() const
Definition TGeoNode.h:99
Base class describing geometry overlaps.
Definition TGeoOverlap.h:37
TGeoVolume * GetSecondVolume() const
Definition TGeoOverlap.h:66
TGeoHMatrix * GetFirstMatrix() const
Definition TGeoOverlap.h:67
Bool_t IsExtrusion() const
Definition TGeoOverlap.h:70
Double_t GetOverlap() const
Definition TGeoOverlap.h:69
TGeoHMatrix * GetSecondMatrix() const
Definition TGeoOverlap.h:68
TGeoVolume * GetFirstVolume() const
Definition TGeoOverlap.h:65
TVirtualGeoTrack * AddTrack(Int_t id, Int_t pdgcode, TObject *part) override
Create a primary TGeoTrack.
void PaintOverlap(void *ovlp, Option_t *option="") override
Paint an overlap.
Double_t fBombZ
Definition TGeoPainter.h:41
void EstimateCameraMove(Double_t tmin, Double_t tmax, Double_t *start, Double_t *end) override
Estimate camera movement between tmin and tmax for best track display.
TObjArray * fVisVolumes
Definition TGeoPainter.h:66
TGeoIteratorPlugin * fPlugin
Definition TGeoPainter.h:65
Double_t fMat[9]
Definition TGeoPainter.h:44
TBuffer3D * fBuffer
Definition TGeoPainter.h:60
void BombTranslation(const Double_t *tr, Double_t *bombtr) override
Get the new 'bombed' translation vector according current exploded view mode.
void Paint(Option_t *option="") override
Paint current geometry according to option.
void EditGeometry(Option_t *option="") override
Start the geometry editor.
void GetViewAngles(Double_t &longitude, Double_t &latitude, Double_t &psi) override
Get the current view angles.
TGeoManager * fGeoManager
Definition TGeoPainter.h:61
void Raytrace(Option_t *option="") override
Raytrace current drawn geometry.
TString fVisBranch
Definition TGeoPainter.h:55
TGeoVolume * fTopVolume
Definition TGeoPainter.h:63
~TGeoPainter() override
Default destructor.
void SetVisOption(Int_t option=0) override
Set drawing mode :
Bool_t IsExplodedView() const override
Bool_t fIsPaintingShape
Definition TGeoPainter.h:54
Double_t fBombR
Definition TGeoPainter.h:42
Int_t ShapeDistancetoPrimitive(const TGeoShape *shape, Int_t numpoints, Int_t px, Int_t py) const override
Returns distance between point px,py on the pad an a shape.
void Draw(Option_t *option="") override
Draw method.
Int_t CountVisibleNodes() override
Count total number of visible nodes.
void DefineColors() const
Define 100 colors with increasing light intensities for each basic color (1-7) Register these colors ...
void SetExplodedView(Int_t iopt=0) override
Set type of exploding view.
void ExecuteManagerEvent(TGeoManager *geom, Int_t event, Int_t px, Int_t py) override
Execute mouse actions on a given volume.
TGeoVolume * fLastVolume
Definition TGeoPainter.h:64
Int_t DistanceToPrimitiveVol(TGeoVolume *vol, Int_t px, Int_t py) override
Compute the closest distance of approach from point px,py to a volume.
void UnbombTranslation(const Double_t *tr, Double_t *bombtr) override
Get the new 'unbombed' translation vector according current exploded view mode.
Bool_t fVisLock
Definition TGeoPainter.h:50
void PaintNode(TGeoNode *node, Option_t *option="", TGeoMatrix *global=nullptr) override
Paint recursively a node and its content according to visualization options.
void DefaultAngles() override
Set default angles for the current view.
TGeoShape * fClippingShape
Definition TGeoPainter.h:62
Bool_t PaintShape(const TGeoShape &shape, Option_t *option) const
Paint the supplied shape into the current 3D viewer.
Bool_t fIsEditable
Definition TGeoPainter.h:67
Int_t CountNodes(TGeoVolume *vol, Int_t level) const
Count number of visible nodes down to a given level.
Double_t fCheckedBox[6]
Definition TGeoPainter.h:43
void DrawShape(TGeoShape *shape, Option_t *option="") override
Draw a shape.
void DrawOnly(Option_t *option="") override
Draw only one volume.
TGeoOverlap * fOverlap
Definition TGeoPainter.h:58
Double_t fBombY
Definition TGeoPainter.h:40
Bool_t fPaintingOverlaps
Definition TGeoPainter.h:52
void SetBombFactors(Double_t bombx=1.3, Double_t bomby=1.3, Double_t bombz=1.3, Double_t bombr=1.3) override
Set cartesian and radial bomb factors for translations.
void DefaultColors() override
Set default volume colors according to tracking media.
void CheckEdit()
Check if Ged library is loaded and load geometry editor classe.
void ClearVisibleVolumes()
Clear the list of visible volumes reset the kVisOnScreen bit for volumes previously in the list.
Int_t GetColor(Int_t base, Float_t light) const override
Get index of a base color with given light intensity (0,1)
void LocalToMasterVect(const Double_t *local, Double_t *master) const
Convert a local vector according view rotation matrix.
void DrawVolume(TGeoVolume *vol, Option_t *option="") override
Draw method.
TGeoVolume * GetDrawnVolume() const override
Get currently drawn volume.
void DrawPolygon(const TGeoPolygon *poly) override
Draw a polygon in 3D.
Int_t fNsegments
Definition TGeoPainter.h:45
void PaintVolume(TGeoVolume *vol, Option_t *option="", TGeoMatrix *global=nullptr) override
Paint recursively a node and its content according to visualization options.
TString fVolInfo
Definition TGeoPainter.h:56
Bool_t fTopVisible
Definition TGeoPainter.h:51
TGeoHMatrix * fGlobal
Definition TGeoPainter.h:59
Int_t fVisLevel
Definition TGeoPainter.h:47
void SetNsegments(Int_t nseg=20) override
Set number of segments to approximate circles.
void AddTrackPoint(Double_t *point, Double_t *box, Bool_t reset=kFALSE) override
Average center of view of all painted tracklets and compute view box.
void GrabFocus(Int_t nfr=0, Double_t dlong=0, Double_t dlat=0, Double_t dpsi=0) override
Move focus to current volume.
void DrawOverlap(void *ovlp, Option_t *option="") override
Draw an overlap.
void SetVisLevel(Int_t level=3) override
Set default level down to which visualization is performed.
void AddSize3D(Int_t numpoints, Int_t numsegs, Int_t numpolys) override
Add numpoints, numsegs, numpolys to the global 3D size.
TGeoNode * fCheckedNode
Definition TGeoPainter.h:57
void PaintPhysicalNode(TGeoPhysicalNode *node, Option_t *option="")
Paints a physical node associated with a path.
Int_t fNVisNodes
Definition TGeoPainter.h:46
Int_t fExplodedView
Definition TGeoPainter.h:49
const char * GetVolumeInfo(const TGeoVolume *volume, Int_t px, Int_t py) const override
Get some info about the current selected volume.
void DrawCurrentPoint(Int_t color) override
Draw current point in the same view.
void DrawPath(const char *path, Option_t *option="") override
Draw all volumes for a given path.
void DrawBatemanSol(TGeoBatemanSol *sol, Option_t *option="") override
Draw the time evolution of a radionuclide.
Bool_t fIsRaytracing
Definition TGeoPainter.h:53
void ExecuteShapeEvent(TGeoShape *shape, Int_t event, Int_t px, Int_t py) override
Execute mouse actions on a given shape.
void ExecuteVolumeEvent(TGeoVolume *volume, Int_t event, Int_t px, Int_t py) override
Execute mouse actions on a given volume.
TGeoPainter(TGeoManager *manager)
Default constructor.
Double_t fBombX
Definition TGeoPainter.h:39
void ModifiedPad(Bool_t update=kFALSE) const override
Check if a pad and view are present and send signal "Modified" to pad.
void SetTopVisible(Bool_t vis=kTRUE) override
Set top geometry volume as visible.
void DrawPanel() override
Int_t fVisOption
Definition TGeoPainter.h:48
Physical nodes are the actual 'touchable' objects in the geometry, representing a path of positioned ...
Int_t GetLevel() const
Bool_t IsVisible() const
Bool_t IsVisibleFull() const
TGeoHMatrix * GetMatrix(Int_t level=-1) const
Return global matrix for node at LEVEL.
Bool_t IsVolAttributes() const
TGeoVolume * GetVolume(Int_t level=-1) const
Return volume associated with node at LEVEL in the branch.
An arbitrary polygon defined by vertices.
Definition TGeoPolygon.h:19
static TClass * Class()
Base abstract class for all shapes.
Definition TGeoShape.h:25
Int_t DistancetoPrimitive(Int_t px, Int_t py) override=0
Computes distance from point (px,py) to the object.
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Stub implementation to avoid forcing implementation at this stage.
static Double_t Big()
Definition TGeoShape.h:94
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const =0
static void SetTransform(TGeoMatrix *matrix)
Set current transformation matrix that applies to shape.
virtual Bool_t IsComposite() const
Definition TGeoShape.h:138
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const =0
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const =0
static TGeoMatrix * GetTransform()
Returns current transformation matrix that applies to shape.
virtual Bool_t Contains(const Double_t *point) const =0
virtual void SetPoints(Double_t *points) const =0
TClass * IsA() const override
Definition TGeoShape.h:179
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
Bool_t IsVisContainers() const
Definition TGeoVolume.h:157
void SetLineWidth(Width_t lwidth) override
Set the line width.
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:174
Int_t GetNdaughters() const
Definition TGeoVolume.h:362
void SetTransparency(Char_t transparency=0)
Definition TGeoVolume.h:376
void SetLineColor(Color_t lcolor) override
Set the line color.
TGeoShape * GetShape() const
Definition TGeoVolume.h:190
Bool_t IsRaytracing() const
Check if the painter is currently ray-tracing the content of this volume.
Bool_t IsVisLeaves() const
Definition TGeoVolume.h:158
void SetLineStyle(Style_t lstyle) override
Set the line style.
Bool_t IsVisOnly() const
Definition TGeoVolume.h:159
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
Char_t GetTransparency() const
Definition TGeoVolume.h:369
virtual Bool_t IsVisible() const
Definition TGeoVolume.h:155
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
void Reset()
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
void Clear(Option_t *option="") override
Remove all objects from the array.
void Add(TObject *obj) override
Definition TObjArray.h:68
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1057
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:203
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
Longptr_t ExecPlugin(int nargs)
A 3D polymarker.
Stopwatch class.
Definition TStopwatch.h:28
Basic string class.
Definition TString.h:138
void ToLower()
Change string to lower-case.
Definition TString.cxx:1189
const char * Data() const
Definition TString.h:384
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2384
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:640
See TView3D.
Definition TView.h:25
virtual void SetPerspective()=0
virtual Double_t GetPsi()=0
virtual void GetWindow(Double_t &u0, Double_t &v0, Double_t &du, Double_t &dv) const =0
virtual Double_t GetLongitude()=0
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
virtual Double_t GetDview() const =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:26
virtual void SetAutoRange(Bool_t autorange=kTRUE)=0
virtual Double_t GetDproj() const =0
virtual Bool_t IsPerspective() const =0
virtual void SetViewChanged(Bool_t flag=kTRUE)=0
virtual void GetRange(Float_t *min, Float_t *max)=0
virtual Double_t GetLatitude()=0
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)=0
virtual void SetView(Double_t longitude, Double_t latitude, Double_t psi, Int_t &irep)=0
Abstract class for geometry painters.
static void SetPainter(const TVirtualGeoPainter *painter)
Static function to set an alternative histogram painter.
Base class for user-defined tracks attached to a geometry.
static void ShowEditor()
Show the global pad editor. Static method.
Abstract 3D shapes viewer.
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
return c2
Definition legend2.C:14
return c3
Definition legend3.C:15
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:251
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:82
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
TLine l
Definition textangle.C:4