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 "TBuffer3D.h"
53#include "TBuffer3DTypes.h"
54#include "TVirtualViewer3D.h"
55#include "TVirtualX.h"
56
57#include <cstring>
58
59namespace {
60
62{
63 return viewer && viewer->IsA() && std::strcmp(viewer->IsA()->GetName(), "TViewerX3D") == 0;
64}
65
66} // namespace
67
68////////////////////////////////////////////////////////////////////////////////
69/// Default constructor.
70
72{
74 if (manager)
76 else {
77 Error("ctor", "No geometry loaded");
78 return;
79 }
81 fNVisNodes = 0;
82 fBombX = 1.3;
83 fBombY = 1.3;
84 fBombZ = 1.3;
85 fBombR = 1.3;
89 fVisBranch = "";
90 fVolInfo = "";
95 fPlugin = nullptr;
96 fVisVolumes = new TObjArray();
97 fOverlap = nullptr;
98 fGlobal = new TGeoHMatrix();
99 fBuffer = new TBuffer3D(TBuffer3DTypes::kGeneric, 20, 3 * 20, 0, 0, 0, 0);
100 fClippingShape = nullptr;
101 fLastVolume = nullptr;
102 fTopVolume = nullptr;
104 memset(&fCheckedBox[0], 0, 6 * sizeof(Double_t));
105
108 DefineColors();
109}
110////////////////////////////////////////////////////////////////////////////////
111/// Default destructor.
112
114{
115 delete fVisVolumes;
116 delete fGlobal;
117 delete fBuffer;
118 if (fPlugin)
119 delete fPlugin;
120}
121////////////////////////////////////////////////////////////////////////////////
122/// Add numpoints, numsegs, numpolys to the global 3D size.
123
125{
126 // Legacy no-op kept for the obsolete x3d sizing interface.
127 (void)numpoints;
128 (void)numsegs;
129 (void)numpolys;
130}
131////////////////////////////////////////////////////////////////////////////////
132/// Create a primary TGeoTrack.
133
138
139////////////////////////////////////////////////////////////////////////////////
140/// Average center of view of all painted tracklets and compute view box.
141
143{
144 static Int_t npoints = 0;
145 static Double_t xmin[3] = {0, 0, 0};
146 static Double_t xmax[3] = {0, 0, 0};
147 Int_t i;
148 if (reset) {
149 memset(box, 0, 6 * sizeof(Double_t));
150 memset(xmin, 0, 3 * sizeof(Double_t));
151 memset(xmax, 0, 3 * sizeof(Double_t));
152 npoints = 0;
153 return;
154 }
155 if (npoints == 0) {
156 for (i = 0; i < 3; i++)
157 xmin[i] = xmax[i] = 0;
158 npoints++;
159 }
160 npoints++;
162 for (i = 0; i < 3; i++) {
163 box[i] += ninv * (point[i] - box[i]);
164 if (point[i] < xmin[i])
165 xmin[i] = point[i];
166 if (point[i] > xmax[i])
167 xmax[i] = point[i];
168 box[i + 3] = 0.5 * (xmax[i] - xmin[i]);
169 }
170}
171
172////////////////////////////////////////////////////////////////////////////////
173/// Get the new 'bombed' translation vector according current exploded view mode.
174
176{
177 memcpy(bombtr, tr, 3 * sizeof(Double_t));
178 switch (fExplodedView) {
179 case kGeoNoBomb: return;
180 case kGeoBombXYZ:
181 bombtr[0] *= fBombX;
182 bombtr[1] *= fBombY;
183 bombtr[2] *= fBombZ;
184 return;
185 case kGeoBombCyl:
186 bombtr[0] *= fBombR;
187 bombtr[1] *= fBombR;
188 bombtr[2] *= fBombZ;
189 return;
190 case kGeoBombSph:
191 bombtr[0] *= fBombR;
192 bombtr[1] *= fBombR;
193 bombtr[2] *= fBombR;
194 return;
195 default: return;
196 }
197}
198
199////////////////////////////////////////////////////////////////////////////////
200/// Clear the list of visible volumes
201/// reset the kVisOnScreen bit for volumes previously in the list
202
204{
205 if (!fVisVolumes)
206 return;
207 TIter next(fVisVolumes);
208 TGeoVolume *vol;
209 while ((vol = (TGeoVolume *)next())) {
211 }
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Define 100 colors with increasing light intensities for each basic color (1-7)
217/// Register these colors at indexes starting with 1000.
218
220{
221 static Int_t color = 0;
222 if (!color) {
224 for (auto icol = 1; icol < 10; ++icol)
225 color = GetColor(icol, 0.5);
226 }
227}
228
229////////////////////////////////////////////////////////////////////////////////
230/// Get index of a base color with given light intensity (0,1)
231
233{
234 using IntMap_t = std::map<Int_t, Int_t>;
235 constexpr Int_t ncolors = 100;
236 constexpr Float_t lmin = 0.25;
237 constexpr Float_t lmax = 0.75;
238 static IntMap_t colmap;
239 Int_t color = base;
240 // Search color in the map
241 auto it = colmap.find(base);
242 if (it != colmap.end())
243 return (it->second + light * (ncolors - 1));
244 // Get color pointer if stored
245 TColor *col_base = gROOT->GetColor(base);
246 if (!col_base) {
247 // If color not defined, use gray palette
248 it = colmap.find(kBlack);
249 if (it != colmap.end())
250 return (it->second + light * (ncolors - 1));
251 col_base = gROOT->GetColor(kBlack);
252 color = 1;
253 }
254 // Create a color palette for col_base
255 Float_t r = 0., g = 0., b = 0., h = 0., l = 0., s = 0.;
256 Double_t red[2], green[2], blue[2];
257 Double_t stop[] = {0., 1.0};
258
259 if (col_base)
260 col_base->GetRGB(r, g, b);
261 TColor::RGB2HLS(r, g, b, h, l, s);
262 TColor::HLS2RGB(h, lmin, s, r, g, b);
263 red[0] = r;
264 green[0] = g;
265 blue[0] = b;
266 TColor::HLS2RGB(h, lmax, s, r, g, b);
267 red[1] = r;
268 green[1] = g;
269 blue[1] = b;
271 colmap[color] = color_map_idx;
272 return (color_map_idx + light * (ncolors - 1));
273}
274
275////////////////////////////////////////////////////////////////////////////////
276/// Get currently drawn volume.
277
279{
280 if (!gPad)
281 return nullptr;
282 return fTopVolume;
283}
284
285////////////////////////////////////////////////////////////////////////////////
286/// Compute the closest distance of approach from point px,py to a volume.
287
289{
290 const Int_t big = 9999;
291 const Int_t inaxis = 7;
292 const Int_t maxdist = 5;
293
294 if (fTopVolume != volume)
295 fTopVolume = volume;
296 TView *view = gPad->GetView();
297 if (!view)
298 return big;
299 TGeoBBox *box;
300 fGlobal->Clear();
302
303 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
304 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
305 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
306 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
307 // return if point not in user area
308 if (px < puxmin - inaxis)
309 return big;
310 if (py > puymin + inaxis)
311 return big;
312 if (px > puxmax + inaxis)
313 return big;
314 if (py < puymax - inaxis)
315 return big;
316
318 gPad->SetSelected(view);
319 Int_t dist = big;
320 // Int_t id;
321
322 if (fPaintingOverlaps) {
326 dist = crt->GetShape()->DistancetoPrimitive(px, py);
327 if (dist < maxdist) {
328 gPad->SetSelected(crt);
329 box = (TGeoBBox *)crt->GetShape();
330 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
331 fCheckedBox[3] = box->GetDX();
332 fCheckedBox[4] = box->GetDY();
333 fCheckedBox[5] = box->GetDZ();
334 return 0;
335 }
338 dist = crt->GetShape()->DistancetoPrimitive(px, py);
339 if (dist < maxdist) {
340 gPad->SetSelected(crt);
341 box = (TGeoBBox *)crt->GetShape();
342 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
343 fCheckedBox[3] = box->GetDX();
344 fCheckedBox[4] = box->GetDY();
345 fCheckedBox[5] = box->GetDZ();
346 return 0;
347 }
348 return big;
349 }
350 // Compute distance to the right edge
351 if ((puxmax + inaxis - px) < 40) {
352 if ((py - puymax + inaxis) < 40) {
353 // when the mouse points to the (40x40) right corner of the pad, the manager class is selected
354 gPad->SetSelected(fGeoManager);
356 box = (TGeoBBox *)volume->GetShape();
357 memcpy(fCheckedBox, box->GetOrigin(), 3 * sizeof(Double_t));
358 fCheckedBox[3] = box->GetDX();
359 fCheckedBox[4] = box->GetDY();
360 fCheckedBox[5] = box->GetDZ();
361 return 0;
362 }
363 // when the mouse points to the (40 pix) right edge of the pad, the top volume is selected
364 gPad->SetSelected(volume);
365 fVolInfo = volume->GetName();
366 box = (TGeoBBox *)volume->GetShape();
367 memcpy(fCheckedBox, box->GetOrigin(), 3 * sizeof(Double_t));
368 fCheckedBox[3] = box->GetDX();
369 fCheckedBox[4] = box->GetDY();
370 fCheckedBox[5] = box->GetDZ();
371 return 0;
372 }
373
374 TGeoVolume *vol = volume;
375 Bool_t vis = vol->IsVisible();
376 // Bool_t drawDaughters = kTRUE;
377 // Do we need to check a branch only?
378 if (volume->IsVisBranch()) {
379 if (!fGeoManager->IsClosed())
380 return big;
383 while (fGeoManager->GetLevel()) {
386 dist = vol->GetShape()->DistancetoPrimitive(px, py);
387 if (dist < maxdist) {
389 box = (TGeoBBox *)vol->GetShape();
390 fGeoManager->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
393 gPad->SetSelected(fCheckedNode);
394 else
395 gPad->SetSelected(vol);
396 fCheckedBox[3] = box->GetDX();
397 fCheckedBox[4] = box->GetDY();
398 fCheckedBox[5] = box->GetDZ();
400 return 0;
401 }
402 fGeoManager->CdUp();
403 }
405 return dist;
406 }
407
408 // Do I need to look for the top volume ?
409 if ((fTopVisible && vis) || !vol->GetNdaughters() || !vol->IsVisDaughters() || vol->IsVisOnly()) {
410 dist = vol->GetShape()->DistancetoPrimitive(px, py);
411 if (dist < maxdist) {
412 fVolInfo = vol->GetName();
413 gPad->SetSelected(vol);
414 box = (TGeoBBox *)vol->GetShape();
415 memcpy(fCheckedBox, box->GetOrigin(), 3 * sizeof(Double_t));
416 fCheckedBox[3] = box->GetDX();
417 fCheckedBox[4] = box->GetDY();
418 fCheckedBox[5] = box->GetDZ();
419 return 0;
420 }
421 if (vol->IsVisOnly() || !vol->GetNdaughters() || !vol->IsVisDaughters())
422 return dist;
423 }
424
425 // Iterate the volume content
426 TGeoIterator next(vol);
427 next.SetTopName(TString::Format("%s_1", vol->GetName()));
429
430 Int_t level, nd;
431 Bool_t last;
432
433 while ((daughter = next())) {
434 vol = daughter->GetVolume();
435 level = next.GetLevel();
436 nd = daughter->GetNdaughters();
437 vis = daughter->IsVisible();
438 if (volume->IsVisContainers()) {
439 if (vis && level <= fVisLevel) {
440 *fGlobal = next.GetCurrentMatrix();
441 dist = vol->GetShape()->DistancetoPrimitive(px, py);
442 if (dist < maxdist) {
443 next.GetPath(fVolInfo);
444 box = (TGeoBBox *)vol->GetShape();
445 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
448 gPad->SetSelected(fCheckedNode);
449 else
450 gPad->SetSelected(vol);
451 fCheckedBox[3] = box->GetDX();
452 fCheckedBox[4] = box->GetDY();
453 fCheckedBox[5] = box->GetDZ();
455 return 0;
456 }
457 }
458 // Check if we have to skip this branch
459 if (level == fVisLevel || !daughter->IsVisDaughters()) {
460 next.Skip();
461 continue;
462 }
463 } else if (volume->IsVisLeaves()) {
464 last = ((nd == 0) || (level == fVisLevel) || (!daughter->IsVisDaughters())) ? kTRUE : kFALSE;
465 if (vis && last) {
466 *fGlobal = next.GetCurrentMatrix();
467 dist = vol->GetShape()->DistancetoPrimitive(px, py);
468 if (dist < maxdist) {
469 next.GetPath(fVolInfo);
470 box = (TGeoBBox *)vol->GetShape();
471 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
474 gPad->SetSelected(fCheckedNode);
475 else
476 gPad->SetSelected(vol);
477 fCheckedBox[3] = box->GetDX();
478 fCheckedBox[4] = box->GetDY();
479 fCheckedBox[5] = box->GetDZ();
481 return 0;
482 }
483 }
484 // Check if we have to skip the branch
485 if (last || !daughter->IsVisDaughters())
486 next.Skip();
487 }
488 }
489 return dist;
490}
491
492////////////////////////////////////////////////////////////////////////////////
493/// Set default angles for the current view.
494
496{
497 if (gPad) {
498 Int_t irep;
499 TView *view = gPad->GetView();
500 if (!view)
501 return;
502 view->SetView(-206, 126, 75, irep);
503 ModifiedPad();
504 }
505}
506
507////////////////////////////////////////////////////////////////////////////////
508/// Set default volume colors according to tracking media
509
511{
513 TGeoVolume *vol;
514 while ((vol = (TGeoVolume *)next()))
516 ModifiedPad();
517}
518
519////////////////////////////////////////////////////////////////////////////////
520/// Count number of visible nodes down to a given level.
521
523{
524 TGeoVolume *vol = volume;
525 Int_t count = 0;
526 Bool_t vis = vol->IsVisible();
527 // Do I need to look for the top volume ?
528 if ((fTopVisible && vis) || !vol->GetNdaughters() || !vol->IsVisDaughters() || vol->IsVisOnly())
529 count++;
530 // Is this the only volume?
531 if (volume->IsVisOnly())
532 return count;
533
534 // Do we need to check a branch only?
535 if (volume->IsVisBranch()) {
538 count = fGeoManager->GetLevel() + 1;
540 return count;
541 }
542 // Iterate the volume content
543 TGeoIterator next(vol);
545 Int_t level, nd;
546 Bool_t last;
547
548 while ((daughter = next())) {
549 // vol = daughter->GetVolume();
550 level = next.GetLevel();
551 nd = daughter->GetNdaughters();
552 vis = daughter->IsVisible();
553 if (volume->IsVisContainers()) {
554 if (vis && level <= rlevel)
555 count++;
556 // Check if we have to skip this branch
557 if (level == rlevel || !daughter->IsVisDaughters()) {
558 next.Skip();
559 continue;
560 }
561 } else if (volume->IsVisLeaves()) {
562 last = ((nd == 0) || (level == rlevel) || (!daughter->IsVisDaughters())) ? kTRUE : kFALSE;
563 if (vis && last)
564 count++;
565 // Check if we have to skip the branch
566 if (last)
567 next.Skip();
568 }
569 }
570 return count;
571}
572
573////////////////////////////////////////////////////////////////////////////////
574/// Count total number of visible nodes.
575
577{
579 Int_t vislevel = fGeoManager->GetVisLevel();
580 // TGeoVolume *top = fGeoManager->GetTopVolume();
581 TGeoVolume *top = fTopVolume;
582 if (maxnodes <= 0 && top) {
583 fNVisNodes = CountNodes(top, vislevel);
584 SetVisLevel(vislevel);
585 return fNVisNodes;
586 }
587 // if (the total number of nodes of the top volume is less than maxnodes
588 // we can visualize everything.
589 // recompute the best visibility level
590 if (!top) {
591 SetVisLevel(vislevel);
592 return 0;
593 }
594 fNVisNodes = -1;
596 for (Int_t level = 1; level < 20; level++) {
597 vislevel = level;
598 Int_t nnodes = CountNodes(top, level);
599 if (top->IsVisOnly() || top->IsVisBranch()) {
600 vislevel = fVisLevel;
602 break;
603 }
604 if (nnodes > maxnodes) {
605 vislevel--;
606 break;
607 }
608 if (nnodes == fNVisNodes) {
609 if (again)
610 break;
611 again = kTRUE;
612 }
614 }
615 SetVisLevel(vislevel);
616 return fNVisNodes;
617}
618
619////////////////////////////////////////////////////////////////////////////////
620/// Check if Ged library is loaded and load geometry editor classe.
621
623{
624 if (fIsEditable)
625 return;
626 if (!TClass::GetClass("TGedEditor"))
627 return;
629 if ((h = gROOT->GetPluginManager()->FindHandler("TGeoManagerEditor"))) {
630 if (h->LoadPlugin() == -1)
631 return;
632 h->ExecPlugin(0);
633 }
635}
636
637////////////////////////////////////////////////////////////////////////////////
638/// Start the geometry editor.
639
641{
642 if (!gPad)
643 return;
644 if (!fIsEditable) {
645 if (!option[0])
646 gPad->GetCanvas()->GetCanvasImp()->ShowEditor();
647 else
649 CheckEdit();
650 }
651 gPad->SetSelected(fGeoManager);
652 gPad->GetCanvas()->Selected(gPad, fGeoManager, kButton1Down);
653}
654
655////////////////////////////////////////////////////////////////////////////////
656/// Draw method.
657
662
663////////////////////////////////////////////////////////////////////////////////
664/// Draw the time evolution of a radionuclide.
665
667{
668 Int_t ncoeff = sol->GetNcoeff();
669 if (!ncoeff)
670 return;
671 Double_t tlo = 0., thi = 0.;
672 Double_t cn = 0., lambda = 0.;
673 Int_t i;
674 sol->GetRange(tlo, thi);
675 Bool_t autorange = (thi == 0.) ? kTRUE : kFALSE;
676
677 // Try to find the optimum range in time.
678 if (autorange)
679 tlo = 0.;
680 sol->GetCoeff(0, cn, lambda);
681 Double_t lambdamin = lambda;
682 TString formula = "";
683 for (i = 0; i < ncoeff; i++) {
684 sol->GetCoeff(i, cn, lambda);
685 formula += TString::Format("%g*exp(-%g*x)", cn, lambda);
686 if (i < ncoeff - 1)
687 formula += "+";
689 lambdamin = lambda;
690 }
691 if (autorange)
692 thi = 10. / lambdamin;
693 // Create a function
694 TF1 *func = new TF1(TString::Format("conc%s", sol->GetElement()->GetName()), formula.Data(), tlo, thi);
695 func->SetTitle(formula + ";time[s]" + TString::Format(";Concentration_of_%s", sol->GetElement()->GetName()));
696 func->SetMinimum(1.e-3);
697 func->SetMaximum(1.25 * TMath::Max(sol->Concentration(tlo), sol->Concentration(thi)));
698 func->SetLineColor(sol->GetLineColor());
699 func->SetLineStyle(sol->GetLineStyle());
700 func->SetLineWidth(sol->GetLineWidth());
701 func->SetMarkerColor(sol->GetMarkerColor());
702 func->SetMarkerStyle(sol->GetMarkerStyle());
703 func->SetMarkerSize(sol->GetMarkerSize());
704 func->Draw(option);
705}
706
707////////////////////////////////////////////////////////////////////////////////
708/// Draw a polygon in 3D.
709
711{
712 Int_t nvert = poly->GetNvert();
713 if (!nvert) {
714 Error("DrawPolygon", "No vertices defined");
715 return;
716 }
717 Int_t nconv = poly->GetNconvex();
718 Double_t *x = new Double_t[nvert + 1];
719 Double_t *y = new Double_t[nvert + 1];
720 poly->GetVertices(x, y);
721 x[nvert] = x[0];
722 y[nvert] = y[0];
723 TGraph *g1 = new TGraph(nvert + 1, x, y);
724 g1->SetTitle(Form("Polygon with %d vertices (outscribed %d)", nvert, nconv));
725 g1->SetLineColor(kRed);
726 g1->SetMarkerColor(kRed);
727 g1->SetMarkerStyle(4);
728 g1->SetMarkerSize(0.8);
729 delete[] x;
730 delete[] y;
731 Double_t *xc = nullptr;
732 Double_t *yc = nullptr;
733 TGraph *g2 = nullptr;
734 if (nconv && !poly->IsConvex()) {
735 xc = new Double_t[nconv + 1];
736 yc = new Double_t[nconv + 1];
737 poly->GetConvexVertices(xc, yc);
738 xc[nconv] = xc[0];
739 yc[nconv] = yc[0];
740 g2 = new TGraph(nconv + 1, xc, yc);
741 g2->SetLineColor(kBlue);
742 g2->SetLineColor(kBlue);
743 g2->SetMarkerColor(kBlue);
744 g2->SetMarkerStyle(21);
745 g2->SetMarkerSize(0.4);
746 delete[] xc;
747 delete[] yc;
748 }
749 if (!gPad) {
750 gROOT->MakeDefCanvas();
751 }
752 g1->Draw("ALP");
753 if (g2)
754 g2->Draw("LP");
755}
756
757////////////////////////////////////////////////////////////////////////////////
758/// Draw method.
759
761{
762 fTopVolume = vol;
763 fLastVolume = nullptr;
765 // if (fVisOption==kGeoVisOnly ||
766 // fVisOption==kGeoVisBranch) fGeoManager->SetVisOption(kGeoVisLeaves);
768 TString opt = option;
769 opt.ToLower();
771 fOverlap = nullptr;
772
773 if (fVisLock) {
776 }
777 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
778 // Clear pad if option "same" not given
779 if (!gPad) {
780 gROOT->MakeDefCanvas();
781 }
782 if (!opt.Contains("same"))
783 gPad->Clear();
784 // append this volume to pad
786
787 // Create a 3-D view
788 TView *view = gPad->GetView();
789 if (!view) {
790 view = TView::CreateView(11, nullptr, nullptr);
791 // Set the view to perform a first autorange (frame) draw.
792 // TViewer3DPad will revert view to normal painting after this
793 view->SetAutoRange(kTRUE);
794 if (has_pad)
795 gPad->Update();
796 }
797 if (!opt.Contains("same"))
798 Paint("range");
799 else
800 Paint(opt);
801 view->SetAutoRange(kFALSE);
802 // If we are drawing into the pad, then the view needs to be
803 // set to perspective
804 // if (!view->IsPerspective()) view->SetPerspective();
805
807
808 // Create a 3D viewer to paint us
809 gPad->GetViewer3D(option);
810}
811
812////////////////////////////////////////////////////////////////////////////////
813/// Draw a shape.
814
816{
817 TString opt = option;
818 opt.ToLower();
820 fOverlap = nullptr;
822
823 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
824 // Clear pad if option "same" not given
825 if (!gPad) {
826 gROOT->MakeDefCanvas();
827 }
828 if (!opt.Contains("same"))
829 gPad->Clear();
830 // append this shape to pad
831 shape->AppendPad(option);
832
833 // Create a 3-D view
834 TView *view = gPad->GetView();
835 if (!view) {
836 view = TView::CreateView(11, nullptr, nullptr);
837 // Set the view to perform a first autorange (frame) draw.
838 // TViewer3DPad will revert view to normal painting after this
839 view->SetAutoRange(kTRUE);
840 if (has_pad)
841 gPad->Update();
842 }
843 PaintShape(shape, "range");
844 view->SetAutoRange(kFALSE);
845 view->SetPerspective();
846 // Create a 3D viewer to paint us
847 gPad->GetViewer3D(option);
848}
849
850////////////////////////////////////////////////////////////////////////////////
851/// Draw an overlap.
852
854{
855 TString opt = option;
858 if (!overlap)
859 return;
860
863 opt.ToLower();
864 if (fVisLock) {
867 }
868 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
869 // Clear pad if option "same" not given
870 if (!gPad) {
871 gROOT->MakeDefCanvas();
872 }
873 if (!opt.Contains("same"))
874 gPad->Clear();
875 // append this volume to pad
876 overlap->AppendPad(option);
877
878 // Create a 3-D view
879 // Create a 3D viewer to paint us
880 gPad->GetViewer3D(option);
881 TView *view = gPad->GetView();
882 if (!view) {
883 view = TView::CreateView(11, nullptr, nullptr);
884 // Set the view to perform a first autorange (frame) draw.
885 // TViewer3DPad will revert view to normal painting after this
886 view->SetAutoRange(kTRUE);
887 PaintOverlap(ovlp, "range");
888 overlap->GetPolyMarker()->Draw("SAME");
889 if (has_pad)
890 gPad->Update();
891 }
892
893 // If we are drawing into the pad, then the view needs to be
894 // set to perspective
895 // if (!view->IsPerspective()) view->SetPerspective();
896 fVisLock = kTRUE;
897}
898
899////////////////////////////////////////////////////////////////////////////////
900/// Draw only one volume.
901
903{
904 TString opt = option;
905 opt.ToLower();
906 if (fVisLock) {
909 }
912 Bool_t has_pad = (gPad == nullptr) ? kFALSE : kTRUE;
913 // Clear pad if option "same" not given
914 if (!gPad) {
915 gROOT->MakeDefCanvas();
916 }
917 if (!opt.Contains("same"))
918 gPad->Clear();
919 // append this volume to pad
922
923 // Create a 3-D view
924 TView *view = gPad->GetView();
925 if (!view) {
926 view = TView::CreateView(11, nullptr, nullptr);
927 // Set the view to perform a first autorange (frame) draw.
928 // TViewer3DPad will revert view to normal painting after this
929 view->SetAutoRange(kTRUE);
931 if (has_pad)
932 gPad->Update();
933 }
934
935 // If we are drawing into the pad, then the view needs to be
936 // set to perspective
937 // if (!view->IsPerspective()) view->SetPerspective();
938 fVisLock = kTRUE;
939}
940
941////////////////////////////////////////////////////////////////////////////////
942/// Draw current point in the same view.
943
945{
946 if (!gPad)
947 return;
948 if (!gPad->GetView())
949 return;
950
951 TVirtualViewer3D *viewer = gPad->GetViewer3D();
952 if (!viewer || IsX3DViewer(viewer))
953 return;
954
956 pm->SetMarkerColor(color);
957 const Double_t *point = fGeoManager->GetCurrentPoint();
958 pm->SetNextPoint(point[0], point[1], point[2]);
959 pm->SetMarkerStyle(8);
960 pm->SetMarkerSize(0.5);
961 pm->Draw("SAME");
962}
963
964////////////////////////////////////////////////////////////////////////////////
965
967
968////////////////////////////////////////////////////////////////////////////////
969/// Draw all volumes for a given path.
970
980
981////////////////////////////////////////////////////////////////////////////////
982/// Estimate camera movement between tmin and tmax for best track display
983
985{
986 if (!gPad)
987 return;
988 TIter next(gPad->GetListOfPrimitives());
990 TObject *obj;
991 Int_t ntracks = 0;
992 Double_t *point = nullptr;
993 AddTrackPoint(point, start, kTRUE);
994 while ((obj = next())) {
995 if (strcmp(obj->ClassName(), "TGeoTrack"))
996 continue;
997 track = (TVirtualGeoTrack *)obj;
998 ntracks++;
999 track->PaintCollect(tmin, start);
1000 }
1001
1002 if (!ntracks)
1003 return;
1004 next.Reset();
1005 AddTrackPoint(point, end, kTRUE);
1006 while ((obj = next())) {
1007 if (strcmp(obj->ClassName(), "TGeoTrack"))
1008 continue;
1009 track = (TVirtualGeoTrack *)obj;
1010 if (!track)
1011 continue;
1012 track->PaintCollect(tmax, end);
1013 }
1014}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// Execute mouse actions on a given volume.
1018
1019void TGeoPainter::ExecuteManagerEvent(TGeoManager * /*geom*/, Int_t event, Int_t /*px*/, Int_t /*py*/)
1020{
1021 if (!gPad)
1022 return;
1023 gPad->SetCursor(kPointer);
1024 switch (event) {
1025 case kButton1Down:
1026 if (!fIsEditable)
1027 CheckEdit();
1028 }
1029}
1030
1031////////////////////////////////////////////////////////////////////////////////
1032/// Execute mouse actions on a given shape.
1033
1034void TGeoPainter::ExecuteShapeEvent(TGeoShape * /*shape*/, Int_t event, Int_t /*px*/, Int_t /*py*/)
1035{
1036 if (!gPad)
1037 return;
1038 gPad->SetCursor(kHand);
1039 switch (event) {
1040 case kButton1Down:
1041 if (!fIsEditable)
1042 CheckEdit();
1043 }
1044}
1045
1046////////////////////////////////////////////////////////////////////////////////
1047/// Execute mouse actions on a given volume.
1048
1049void TGeoPainter::ExecuteVolumeEvent(TGeoVolume * /*volume*/, Int_t event, Int_t /*px*/, Int_t /*py*/)
1050{
1051 if (!gPad)
1052 return;
1053 if (!fIsEditable)
1054 CheckEdit();
1055 // if (fIsRaytracing) return;
1056 // Bool_t istop = (volume==fTopVolume)?kTRUE:kFALSE;
1057 // if (istop) gPad->SetCursor(kHand);
1058 // else gPad->SetCursor(kPointer);
1059 gPad->SetCursor(kHand);
1060 // static Int_t width, color;
1061 switch (event) {
1062 case kMouseEnter:
1063 // width = volume->GetLineWidth();
1064 // color = volume->GetLineColor();
1065 break;
1066
1067 case kMouseLeave:
1068 // volume->SetLineWidth(width);
1069 // volume->SetLineColor(color);
1070 break;
1071
1072 case kButton1Down:
1073 // volume->SetLineWidth(3);
1074 // volume->SetLineColor(2);
1075 // gPad->Modified();
1076 // gPad->Update();
1077 break;
1078
1079 case kButton1Up:
1080 // volume->SetLineWidth(width);
1081 // volume->SetLineColor(color);
1082 // gPad->Modified();
1083 // gPad->Update();
1084 break;
1085
1086 case kButton1Double:
1087 gPad->SetCursor(kWatch);
1088 GrabFocus();
1089 break;
1090 }
1091}
1092
1093////////////////////////////////////////////////////////////////////////////////
1094/// Get some info about the current selected volume.
1095
1096const char *TGeoPainter::GetVolumeInfo(const TGeoVolume *volume, Int_t /*px*/, Int_t /*py*/) const
1097{
1098 static TString info;
1099 info = "";
1100 if (!gPad)
1101 return info;
1102 if (fPaintingOverlaps) {
1103 if (!fOverlap) {
1104 info = "wrong overlapping flag";
1105 return info;
1106 }
1108 if (fOverlap->IsExtrusion())
1109 ovtype = "EXTRUSION";
1110 else
1111 ovtype = "OVERLAP";
1112 if (volume == fOverlap->GetFirstVolume())
1113 name = volume->GetName();
1114 else
1116 info = TString::Format("%s: %s of %g", name.Data(), ovtype.Data(), fOverlap->GetOverlap());
1117 return info;
1118 } else
1119 info = TString::Format("%s, shape=%s", fVolInfo.Data(), volume->GetShape()->ClassName());
1120 return info;
1121}
1122
1123////////////////////////////////////////////////////////////////////////////////
1124/// Get the current view angles.
1125
1127{
1128 if (!gPad)
1129 return;
1130 TView *view = gPad->GetView();
1131 if (!view)
1132 return;
1133 longitude = view->GetLongitude();
1134 latitude = view->GetLatitude();
1135 psi = view->GetPsi();
1136}
1137
1138////////////////////////////////////////////////////////////////////////////////
1139/// Move focus to current volume
1140
1142{
1143 if (!gPad)
1144 return;
1145 TView *view = gPad->GetView();
1146 if (!view)
1147 return;
1149 printf("Woops!!!\n");
1151 memcpy(&fCheckedBox[0], box->GetOrigin(), 3 * sizeof(Double_t));
1152 fCheckedBox[3] = box->GetDX();
1153 fCheckedBox[4] = box->GetDY();
1154 fCheckedBox[5] = box->GetDZ();
1155 }
1156 view->SetPerspective();
1158 Int_t nframes = nfr;
1159 if (nfr == 0) {
1160 nframes = 1;
1161 if (nvols < 1500)
1162 nframes = 10;
1163 if (nvols < 1000)
1164 nframes = 20;
1165 if (nvols < 200)
1166 nframes = 50;
1167 if (nvols < 100)
1168 nframes = 100;
1169 }
1171}
1172
1173////////////////////////////////////////////////////////////////////////////////
1174/// Convert a local vector according view rotation matrix
1175
1177{
1178 for (Int_t i = 0; i < 3; i++)
1179 master[i] = -local[0] * fMat[i] - local[1] * fMat[i + 3] - local[2] * fMat[i + 6];
1180}
1181
1182////////////////////////////////////////////////////////////////////////////////
1183/// Check if a pad and view are present and send signal "Modified" to pad.
1184
1186{
1187 if (!gPad)
1188 return;
1189 if (update) {
1190 gPad->Update();
1191 return;
1192 }
1193 TView *view = gPad->GetView();
1194 if (!view)
1195 return;
1196 view->SetViewChanged();
1197 gPad->Modified();
1198 if (gROOT->FromPopUp())
1199 gPad->Update();
1200}
1201
1202////////////////////////////////////////////////////////////////////////////////
1203/// Paint current geometry according to option.
1204
1206{
1207 if (!fGeoManager || !fTopVolume)
1208 return;
1210 if (gPad)
1211 is_padviewer = (!strcmp(gPad->GetViewer3D()->ClassName(), "TViewer3DPad")) ? kTRUE : kFALSE;
1212
1216 else if (fTopVolume->IsVisLeaves())
1218 else if (fTopVolume->IsVisOnly())
1220 else if (fTopVolume->IsVisBranch())
1222
1223 if (!fIsRaytracing || !is_padviewer) {
1224 if (fGeoManager->IsDrawingExtra()) {
1225 // loop the list of physical volumes
1226 fGeoManager->CdTop();
1228 Int_t nnodes = nodeList->GetEntriesFast();
1229 Int_t inode;
1230 TGeoPhysicalNode *node;
1231 for (inode = 0; inode < nnodes; inode++) {
1232 node = (TGeoPhysicalNode *)nodeList->UncheckedAt(inode);
1234 }
1235 } else {
1237 }
1238 fVisLock = kTRUE;
1239 }
1240 // Check if we have to raytrace (only in pad)
1242 Raytrace();
1243}
1244
1245////////////////////////////////////////////////////////////////////////////////
1246/// Paint an overlap.
1247
1249{
1250 if (!fGeoManager)
1251 return;
1253 if (!overlap)
1254 return;
1255 Int_t color, transparency;
1256 if (fOverlap != overlap)
1257 fOverlap = overlap;
1260 TGeoVolume *vol;
1261 TGeoVolume *vol1 = overlap->GetFirstVolume();
1262 TGeoVolume *vol2 = overlap->GetSecondVolume();
1263 TGeoHMatrix *matrix1 = overlap->GetFirstMatrix();
1264 TGeoHMatrix *matrix2 = overlap->GetSecondMatrix();
1265 //
1266 vol = vol1;
1267 *hmat = matrix1;
1268 fGeoManager->SetMatrixReflection(matrix1->IsReflection());
1269 if (!fVisLock)
1270 fVisVolumes->Add(vol);
1272 color = vol->GetLineColor();
1274 vol->SetLineColor(kGreen);
1275 vol->SetTransparency(40);
1276 if (!strstr(option, "range"))
1277 ((TAttLine *)vol)->Modify();
1278 PaintShape(*(vol->GetShape()), option);
1279 vol->SetLineColor(color);
1281 vol = vol2;
1282 *hmat = matrix2;
1283 fGeoManager->SetMatrixReflection(matrix2->IsReflection());
1284 if (!fVisLock)
1285 fVisVolumes->Add(vol);
1287 color = vol->GetLineColor();
1289 vol->SetLineColor(kBlue);
1290 vol->SetTransparency(40);
1291 if (!strstr(option, "range"))
1292 ((TAttLine *)vol)->Modify();
1293 PaintShape(*(vol->GetShape()), option);
1294 vol->SetLineColor(color);
1297 fVisLock = kTRUE;
1298}
1299
1300////////////////////////////////////////////////////////////////////////////////
1301/// Paint recursively a node and its content according to visualization options.
1302
1307
1308////////////////////////////////////////////////////////////////////////////////
1309/// Paint recursively a node and its content according to visualization options.
1310
1312{
1313 if (fTopVolume != top) {
1315 fVisLock = kFALSE;
1316 }
1317 fTopVolume = top;
1318 if (!fVisLevel)
1319 return;
1320 TGeoVolume *vol = top;
1321 if (global)
1322 *fGlobal = *global;
1323 else
1324 fGlobal->Clear();
1327 Bool_t vis = (top->IsVisible() && !top->IsAssembly());
1328 Int_t transparency = 0;
1329
1330 // Update pad attributes in case we need to paint VOL
1331 if (!strstr(option, "range"))
1332 ((TAttLine *)vol)->Modify();
1333
1334 // Do we need to draw a branch ?
1335 if (top->IsVisBranch()) {
1338 // while (fGeoManager->GetLevel()) {
1340 if (!fVisLock) {
1341 fVisVolumes->Add(vol);
1343 }
1346 vol->SetTransparency(40);
1347 if (!strstr(option, "range"))
1348 ((TAttLine *)vol)->Modify();
1349 if (global) {
1350 *fGlobal = *global;
1352 } else {
1354 }
1356 PaintShape(*(vol->GetShape()), option);
1358 fGeoManager->CdUp();
1359 // }
1360 fVisLock = kTRUE;
1363 return;
1364 }
1365
1366 // Do I need to draw the top volume ?
1367 if ((fTopVisible && vis) || !top->GetNdaughters() || !top->IsVisDaughters() || top->IsVisOnly()) {
1370 PaintShape(*(vol->GetShape()), option);
1371 if (!fVisLock && !vol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1372 fVisVolumes->Add(vol);
1374 }
1375 if (top->IsVisOnly() || !top->GetNdaughters() || !top->IsVisDaughters()) {
1376 fVisLock = kTRUE;
1377 return;
1378 }
1379 }
1380
1381 // Iterate the volume content
1382 TGeoIterator next(vol);
1383 if (fPlugin)
1384 next.SetUserPlugin(fPlugin);
1386 // TGeoMatrix *glmat;
1387 Int_t level, nd;
1388 Bool_t last;
1389 Int_t line_color = 0, line_width = 0, line_style = 0;
1390 while ((daughter = next())) {
1391 vol = daughter->GetVolume();
1393 level = next.GetLevel();
1394 nd = daughter->GetNdaughters();
1395 vis = daughter->IsVisible();
1397 if (top->IsVisContainers()) {
1398 if (vis && level <= fVisLevel) {
1399 if (fPlugin) {
1400 line_color = vol->GetLineColor();
1401 line_width = vol->GetLineWidth();
1402 line_style = vol->GetLineStyle();
1405 }
1406 if (!strstr(option, "range"))
1407 ((TAttLine *)vol)->Modify();
1408 if (global) {
1409 *fGlobal = *global;
1410 *fGlobal *= *next.GetCurrentMatrix();
1411 } else {
1412 *fGlobal = next.GetCurrentMatrix();
1413 }
1416 if (fPlugin) {
1421 }
1422 if (!fVisLock && !daughter->IsOnScreen()) {
1423 fVisVolumes->Add(vol);
1425 }
1426 }
1427 // Check if we have to skip this branch
1428 if (!drawDaughters || level == fVisLevel || !daughter->IsVisDaughters()) {
1429 next.Skip();
1430 continue;
1431 }
1432 } else if (top->IsVisLeaves()) {
1433 last = ((nd == 0) || (level == fVisLevel) || (!daughter->IsVisDaughters())) ? kTRUE : kFALSE;
1434 if (vis && last) {
1435 if (fPlugin) {
1436 line_color = vol->GetLineColor();
1437 line_width = vol->GetLineWidth();
1438 line_style = vol->GetLineStyle();
1441 }
1442 if (!strstr(option, "range"))
1443 ((TAttLine *)vol)->Modify();
1444 if (global) {
1445 *fGlobal = *global;
1446 *fGlobal *= *next.GetCurrentMatrix();
1447 } else {
1448 *fGlobal = next.GetCurrentMatrix();
1449 }
1452 if (fPlugin) {
1457 }
1458 if (!fVisLock && !daughter->IsOnScreen()) {
1459 fVisVolumes->Add(vol);
1461 }
1462 }
1463 // Check if we have to skip the branch
1464 if (!drawDaughters || last || !daughter->IsVisDaughters())
1465 next.Skip();
1466 }
1467 }
1468 if (fPlugin)
1469 fPlugin->SetIterator(nullptr);
1471 fVisLock = kTRUE;
1472}
1473
1474////////////////////////////////////////////////////////////////////////////////
1475/// Paint the supplied shape into the current 3D viewer
1476
1478{
1480
1481 TVirtualViewer3D *viewer = gPad->GetViewer3D();
1482
1483 if (!viewer || shape.IsA() == TGeoShapeAssembly::Class()) {
1484 return addDaughters;
1485 }
1486
1487 // For non-composite shapes we are the main paint method & perform the negotiation
1488 // with the viewer here
1489 if (!shape.IsComposite()) {
1490 // Does viewer prefer local frame positions?
1491 Bool_t localFrame = viewer->PreferLocalFrame();
1492 // Perform first fetch of buffer from the shape and try adding it
1493 // to the viewer
1494 const TBuffer3D &buffer =
1496 Int_t reqSections = viewer->AddObject(buffer, &addDaughters);
1497
1498 // If the viewer requires additional sections fetch from the shape (if possible)
1499 // and add again
1502 viewer->AddObject(buffer, &addDaughters);
1503 }
1504 }
1505 // Composite shapes have their own internal hierarchy of shapes, each
1506 // of which generate a filled TBuffer3D. Therefore we can't pass up a
1507 // single buffer to here. So as a special case the TGeoCompositeShape
1508 // performs it's own painting & negotiation with the viewer.
1509 else {
1510 const TGeoCompositeShape *composite = static_cast<const TGeoCompositeShape *>(&shape);
1511
1512 // We need the addDaughters flag returned from the viewer from paint
1513 // so can't use the normal TObject::Paint()
1514 // TGeoHMatrix *matrix = (TGeoHMatrix*)TGeoShape::GetTransform();
1515 // if (viewer->PreferLocalFrame()) matrix->Clear();
1516 addDaughters = composite->PaintComposite(option);
1517 }
1518
1519 return addDaughters;
1520}
1521
1522////////////////////////////////////////////////////////////////////////////////
1523/// Paint an overlap.
1524
1532
1533////////////////////////////////////////////////////////////////////////////////
1534/// Paints a physical node associated with a path.
1535
1537{
1538 if (!node->IsVisible())
1539 return;
1540 Int_t level = node->GetLevel();
1541 Int_t i, col, wid, sty;
1542 TGeoShape *shape;
1546 if (!node->IsVisibleFull()) {
1547 // Paint only last node in the branch
1548 vcrt = node->GetVolume();
1549 if (!strstr(option, "range"))
1550 ((TAttLine *)vcrt)->Modify();
1551 shape = vcrt->GetShape();
1552 *matrix = node->GetMatrix();
1553 fGeoManager->SetMatrixReflection(matrix->IsReflection());
1555 if (!node->IsVolAttributes() && !strstr(option, "range")) {
1556 col = vcrt->GetLineColor();
1557 wid = vcrt->GetLineWidth();
1558 sty = vcrt->GetLineStyle();
1559 vcrt->SetLineColor(node->GetLineColor());
1560 vcrt->SetLineWidth(node->GetLineWidth());
1561 vcrt->SetLineStyle(node->GetLineStyle());
1562 ((TAttLine *)vcrt)->Modify();
1563 PaintShape(*shape, option);
1564 vcrt->SetLineColor(col);
1565 vcrt->SetLineWidth(wid);
1566 vcrt->SetLineStyle(sty);
1567 } else {
1568 PaintShape(*shape, option);
1569 }
1570 } else {
1571 // Paint full branch, except top node
1572 for (i = 1; i <= level; i++) {
1573 vcrt = node->GetVolume(i);
1574 if (!strstr(option, "range"))
1575 ((TAttLine *)vcrt)->Modify();
1576 shape = vcrt->GetShape();
1577 *matrix = node->GetMatrix(i);
1578 fGeoManager->SetMatrixReflection(matrix->IsReflection());
1580 if (!node->IsVolAttributes() && !strstr(option, "range")) {
1581 col = vcrt->GetLineColor();
1582 wid = vcrt->GetLineWidth();
1583 sty = vcrt->GetLineStyle();
1584 vcrt->SetLineColor(node->GetLineColor());
1585 vcrt->SetLineWidth(node->GetLineWidth());
1586 vcrt->SetLineStyle(node->GetLineStyle());
1587 ((TAttLine *)vcrt)->Modify();
1588 PaintShape(*shape, option);
1589 vcrt->SetLineColor(col);
1590 vcrt->SetLineWidth(wid);
1591 vcrt->SetLineStyle(sty);
1592 } else {
1593 PaintShape(*shape, option);
1594 }
1595 }
1596 }
1598}
1599
1600////////////////////////////////////////////////////////////////////////////////
1601/// Raytrace current drawn geometry
1602
1604{
1605 if (!gPad || gPad->IsBatch())
1606 return;
1607 TView *view = gPad->GetView();
1608 if (!view)
1609 return;
1612 if (top != fTopVolume)
1614 if (!view->IsPerspective())
1615 view->SetPerspective();
1616 gVirtualX->SetMarkerSize(1);
1617 gVirtualX->SetMarkerStyle(1);
1620 Double_t lat = view->GetLatitude();
1621 Double_t longit = view->GetLongitude();
1622 Double_t psi = view->GetPsi();
1629 fMat[0] = c1 * c3 - s1 * c2 * s3;
1630 fMat[1] = c1 * s3 + s1 * c2 * c3;
1631 fMat[2] = s1 * s2;
1632
1633 fMat[3] = -s1 * c3 - c1 * c2 * s3;
1634 fMat[4] = -s1 * s3 + c1 * c2 * c3;
1635 fMat[5] = c1 * s2;
1636
1637 fMat[6] = s2 * s3;
1638 fMat[7] = -s2 * c3;
1639 fMat[8] = c2;
1640 Double_t u0, v0, du, dv;
1641 view->GetWindow(u0, v0, du, dv);
1642 Double_t dview = view->GetDview();
1643 Double_t dproj = view->GetDproj();
1644 Double_t local[3] = {0, 0, 1};
1645 Double_t dir[3], normal[3];
1647 Double_t min[3], max[3];
1648 view->GetRange(min, max);
1649 Double_t cov[3];
1650 for (Int_t i = 0; i < 3; i++)
1651 cov[i] = 0.5 * (min[i] + max[i]);
1652 Double_t cop[3];
1653 for (Int_t i = 0; i < 3; i++)
1654 cop[i] = cov[i] - dir[i] * dview;
1655 fGeoManager->InitTrack(cop, dir);
1658 if (fClippingShape)
1660 Int_t px, py;
1663 pxmin = gPad->UtoAbsPixel(0);
1664 pxmax = gPad->UtoAbsPixel(1);
1665 pymin = gPad->VtoAbsPixel(1);
1666 pymax = gPad->VtoAbsPixel(0);
1667 TGeoNode *next = nullptr;
1668 TGeoNode *nextnode = nullptr;
1669 Double_t step, steptot;
1670 Double_t *norm;
1671 const Double_t *point = fGeoManager->GetCurrentPoint();
1672 Double_t *ppoint = (Double_t *)point;
1673 Double_t tosource[3];
1674 Double_t calf;
1675 Double_t phi = 45. * krad;
1676 tosource[0] = -dir[0] * TMath::Cos(phi) + dir[1] * TMath::Sin(phi);
1677 tosource[1] = -dir[0] * TMath::Sin(phi) - dir[1] * TMath::Cos(phi);
1678 tosource[2] = -dir[2];
1679
1681
1682 Bool_t done;
1683 // Int_t istep;
1684 Int_t base_color, color;
1687 TPoint *pxy = new TPoint[1];
1689 Int_t up;
1690 Int_t ntotal = pxmax * pymax;
1691 Int_t nrays = 0;
1692 TStopwatch *timer = new TStopwatch();
1693 timer->Start();
1694 for (px = pxmin; px < pxmax; px++) {
1695 for (py = pymin; py < pymax; py++) {
1696 if ((nrays % 100) == 0)
1697 checker->OpProgress("Raytracing", nrays, ntotal, timer, kFALSE);
1698 nrays++;
1699 base_color = 1;
1700 steptot = 0;
1702 xloc = gPad->AbsPixeltoX(pxmin + pxmax - px);
1703 xloc = xloc * du - u0;
1704 yloc = gPad->AbsPixeltoY(pymin + pymax - py);
1705 yloc = yloc * dv - v0;
1707 local[0] = xloc / modloc;
1708 local[1] = yloc / modloc;
1709 local[2] = dproj / modloc;
1715 // fGeoManager->InitTrack(cop,dir);
1716 // current ray pointing to pixel (px,py)
1717 done = kFALSE;
1718 norm = nullptr;
1719 // propagate to the clipping shape if any
1720 if (fClippingShape) {
1721 if (inclip) {
1724 } else {
1726 stemin = 0;
1727 }
1728 }
1729
1730 while (!done) {
1731 if (fClippingShape) {
1732 if (stemin > 1E10)
1733 break;
1734 if (stemin > 0) {
1735 // we are inside clipping shape
1737 next = fGeoManager->Step();
1738 steptot = 0;
1739 stemin = 0;
1740 up = 0;
1741 while (next) {
1742 // we found something after clipping region
1743 nextvol = next->GetVolume();
1744 if (nextvol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1745 done = kTRUE;
1746 base_color = nextvol->GetLineColor();
1748 norm = normal;
1749 break;
1750 }
1751 up++;
1752 next = fGeoManager->GetMother(up);
1753 }
1754 if (done)
1755 break;
1757 fGeoManager->SetStep(1E-3);
1758 while (inclip) {
1759 fGeoManager->Step();
1761 }
1763 }
1764 }
1766 step = fGeoManager->GetStep();
1767 if (step > 1E10)
1768 break;
1769 steptot += step;
1770 next = nextnode;
1771 // Check the step
1772 if (fClippingShape) {
1773 if (steptot > stemax) {
1774 steptot = 0;
1776 if (inclip) {
1779 continue;
1780 } else {
1781 stemin = 0;
1783 }
1784 }
1785 }
1786 // Check if next node is visible
1787 if (!nextnode)
1788 continue;
1789 nextvol = nextnode->GetVolume();
1790 if (nextvol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1791 done = kTRUE;
1792 base_color = nextvol->GetLineColor();
1793 next = nextnode;
1794 break;
1795 }
1796 }
1797 if (!done)
1798 continue;
1799 // current ray intersect a visible volume having color=base_color
1800 if (rtMode > 0) {
1803 for (Int_t i = 0; i < 3; ++i)
1804 local[i] += 1.E-8 * dir[i];
1805 step = next->GetVolume()->GetShape()->DistFromInside(local, dir, 3);
1806 for (Int_t i = 0; i < 3; ++i)
1807 local[i] += step * dir[i];
1808 next->GetVolume()->GetShape()->ComputeNormal(local, dir, normal);
1809 norm = normal;
1810 } else {
1811 if (!norm)
1813 if (!norm)
1814 continue;
1815 }
1816 calf = norm[0] * tosource[0] + norm[1] * tosource[1] + norm[2] * tosource[2];
1818 color = GetColor(base_color, light);
1819 // Now we know the color of the pixel, just draw it
1820 gVirtualX->SetMarkerColor(color);
1821 pxy[0].fX = px;
1822 pxy[0].fY = py;
1823 gVirtualX->DrawPolyMarker(1, pxy);
1824 }
1825 }
1826 delete[] pxy;
1827 timer->Stop();
1828 checker->OpProgress("Raytracing", nrays, ntotal, timer, kTRUE);
1829 delete timer;
1830}
1831
1832////////////////////////////////////////////////////////////////////////////////
1833/// Set cartesian and radial bomb factors for translations.
1834
1844
1845////////////////////////////////////////////////////////////////////////////////
1846/// Set type of exploding view.
1847
1849{
1850 if ((ibomb < 0) || (ibomb > 3)) {
1851 Warning("SetExplodedView", "exploded view can be 0-3");
1852 return;
1853 }
1854 if ((Int_t)ibomb == fExplodedView)
1855 return;
1856 Bool_t change = (gPad == nullptr) ? kFALSE : kTRUE;
1857
1858 if (ibomb == kGeoNoBomb) {
1860 }
1861 if (ibomb == kGeoBombXYZ) {
1863 }
1864 if (ibomb == kGeoBombCyl) {
1866 }
1867 if (ibomb == kGeoBombSph) {
1869 }
1871 if (change)
1872 ModifiedPad();
1873}
1874
1875////////////////////////////////////////////////////////////////////////////////
1876/// Set number of segments to approximate circles.
1877
1879{
1880 if (nseg < 3) {
1881 Warning("SetNsegments", "number of segments should be > 2");
1882 return;
1883 }
1884 if (fNsegments == nseg)
1885 return;
1886 fNsegments = nseg;
1887 ModifiedPad();
1888}
1889
1890////////////////////////////////////////////////////////////////////////////////
1891/// Set default level down to which visualization is performed
1892
1894{
1895 if (level == fVisLevel && fLastVolume == fTopVolume)
1896 return;
1897 fVisLevel = level;
1898 if (!fTopVolume)
1899 return;
1900 if (fVisLock) {
1902 fVisLock = kFALSE;
1903 }
1904 if (!fLastVolume) {
1905 // printf("--- Drawing %6d nodes with %d visible levels\n",fNVisNodes,fVisLevel);
1906 return;
1907 }
1908 if (!gPad)
1909 return;
1910 if (gPad->GetView()) {
1911 // printf("--- Drawing %6d nodes with %d visible levels\n",fNVisNodes,fVisLevel);
1912 ModifiedPad();
1913 }
1914}
1915
1916////////////////////////////////////////////////////////////////////////////////
1917/// Set top geometry volume as visible.
1918
1920{
1921 if (fTopVisible == vis)
1922 return;
1923 fTopVisible = vis;
1924 ModifiedPad();
1925}
1926
1927////////////////////////////////////////////////////////////////////////////////
1928/// Set drawing mode :
1929/// - option=0 (default) all nodes drawn down to vislevel
1930/// - option=1 leaves and nodes at vislevel drawn
1931/// - option=2 path is drawn
1932
1934{
1935 if ((fVisOption < 0) || (fVisOption > 4)) {
1936 Warning("SetVisOption", "wrong visualization option");
1937 return;
1938 }
1939
1940 if (option == kGeoVisChanged) {
1941 if (fVisLock) {
1943 fVisLock = kFALSE;
1944 }
1945 ModifiedPad();
1946 return;
1947 }
1948
1949 if (fTopVolume) {
1951 att->SetAttBit(TGeoAtt::kVisBranch, kFALSE);
1952 att->SetAttBit(TGeoAtt::kVisContainers, kFALSE);
1953 att->SetAttBit(TGeoAtt::kVisOnly, kFALSE);
1954 switch (option) {
1955 case kGeoVisDefault: att->SetAttBit(TGeoAtt::kVisContainers, kTRUE); break;
1956 case kGeoVisLeaves: break;
1957 case kGeoVisOnly: att->SetAttBit(TGeoAtt::kVisOnly, kTRUE); break;
1958 }
1959 }
1960
1961 if (fVisOption == option)
1962 return;
1964 if (fVisLock) {
1966 fVisLock = kFALSE;
1967 }
1968 ModifiedPad();
1969}
1970
1971////////////////////////////////////////////////////////////////////////////////
1972/// Returns distance between point px,py on the pad an a shape.
1973
1975{
1976 const Int_t inaxis = 7;
1977 const Int_t maxdist = 5;
1978 const Int_t big = 9999;
1979 Int_t dist = big;
1980 if (!gPad)
1981 return dist;
1982 TView *view = gPad->GetView();
1983 if (!(numpoints && view))
1984 return dist;
1985 if (shape->IsA() == TGeoShapeAssembly::Class())
1986 return dist;
1987
1988 if (fIsPaintingShape) {
1989 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
1990 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
1991 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
1992 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
1993 // return if point not in user area
1994 if (px < puxmin - inaxis)
1995 return big;
1996 if (py > puymin + inaxis)
1997 return big;
1998 if (px > puxmax + inaxis)
1999 return big;
2000 if (py < puymax - inaxis)
2001 return big;
2002 if ((puxmax + inaxis - px) < 40) {
2003 // when the mouse points to the (40 pix) right edge of the pad, the manager class is selected
2004 gPad->SetSelected(fGeoManager);
2005 return 0;
2006 }
2007 }
2008
2009 fBuffer->SetRawSizes(numpoints, 3 * numpoints, 0, 0, 0, 0);
2011 shape->SetPoints(points);
2012 Double_t dpoint2, x1, y1, xndc[3];
2013 Double_t dmaster[3];
2014 Int_t j;
2015 for (Int_t i = 0; i < numpoints; i++) {
2016 j = 3 * i;
2017 TGeoShape::GetTransform()->LocalToMaster(&points[j], dmaster);
2018 points[j] = dmaster[0];
2019 points[j + 1] = dmaster[1];
2020 points[j + 2] = dmaster[2];
2021 view->WCtoNDC(&points[j], xndc);
2022 x1 = gPad->XtoAbsPixel(xndc[0]);
2023 y1 = gPad->YtoAbsPixel(xndc[1]);
2024 dpoint2 = (px - x1) * (px - x1) + (py - y1) * (py - y1);
2025 if (dpoint2 < dist)
2026 dist = (Int_t)dpoint2;
2027 }
2028 if (dist > 100)
2029 return dist;
2030 dist = Int_t(TMath::Sqrt(Double_t(dist)));
2031 if (dist < maxdist && fIsPaintingShape)
2032 gPad->SetSelected((TObject *)shape);
2033 return dist;
2034}
2035
2036////////////////////////////////////////////////////////////////////////////////
2037/// Get the new 'unbombed' translation vector according current exploded view mode.
2038
2040{
2041 memcpy(bombtr, tr, 3 * sizeof(Double_t));
2042 switch (fExplodedView) {
2043 case kGeoNoBomb: return;
2044 case kGeoBombXYZ:
2045 bombtr[0] /= fBombX;
2046 bombtr[1] /= fBombY;
2047 bombtr[2] /= fBombZ;
2048 return;
2049 case kGeoBombCyl:
2050 bombtr[0] /= fBombR;
2051 bombtr[1] /= fBombR;
2052 bombtr[2] /= fBombZ;
2053 return;
2054 case kGeoBombSph:
2055 bombtr[0] /= fBombR;
2056 bombtr[1] /= fBombR;
2057 bombtr[2] /= fBombR;
2058 return;
2059 default: return;
2060 }
2061}
@ 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:376
@ kHand
Definition GuiTypes.h:375
@ kPointer
Definition GuiTypes.h:376
#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:148
R__EXTERN TGeoManager * gGeoManager
float xmin
float xmax
#define gROOT
Definition TROOT.h:426
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
#define gPad
#define gVirtualX
Definition TVirtualX.h:368
const_iterator end() const
Line Attributes class.
Definition TAttLine.h:21
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:36
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:46
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:38
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:47
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:44
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:37
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:41
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:43
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:48
Generic 3D primitive description class.
Definition TBuffer3D.h:18
@ kBoundingBox
Definition TBuffer3D.h:60
@ kShapeSpecific
Definition TBuffer3D.h:61
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:122
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:2994
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:2743
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:3449
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:3613
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:3462
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:18
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:459
void Clear(Option_t *option="") override
clear the data for this matrix
void SetIterator(const TGeoIterator *iter)
Definition TGeoNode.h:238
virtual void ProcessNode()=0
A geometry iterator.
Definition TGeoNode.h:249
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:295
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:46
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:39
Bool_t IsReflection() const
Definition TGeoMatrix.h:67
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:100
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:95
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:139
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:181
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
Bool_t IsVisContainers() const
Definition TGeoVolume.h:158
void SetLineWidth(Width_t lwidth) override
Set the line width.
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:175
Int_t GetNdaughters() const
Definition TGeoVolume.h:363
void SetTransparency(Char_t transparency=0)
Definition TGeoVolume.h:378
void SetLineColor(Color_t lcolor) override
Set the line color.
TGeoShape * GetShape() const
Definition TGeoVolume.h:191
Bool_t IsRaytracing() const
Check if the painter is currently ray-tracing the content of this volume.
Bool_t IsVisLeaves() const
Definition TGeoVolume.h:159
void SetLineStyle(Style_t lstyle) override
Set the line style.
Bool_t IsVisOnly() const
Definition TGeoVolume.h:160
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
Char_t GetTransparency() const
Definition TGeoVolume.h:370
virtual Bool_t IsVisible() const
Definition TGeoVolume.h:156
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:42
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:224
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1081
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:201
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1095
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:641
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:249
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:122
TLine l
Definition textangle.C:4