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 "TCanvas.h"
28#include "TCanvasImp.h"
29#include "TH2F.h"
30#include "TF1.h"
31#include "TGraph.h"
32#include "TPluginManager.h"
33#include "TVirtualPadEditor.h"
34#include "TStopwatch.h"
35
36#include "TPolyMarker3D.h"
37
38#include "TGeoAtt.h"
39#include "TGeoVolume.h"
40#include "TGeoNode.h"
41#include "TGeoElement.h"
42#include "TGeoManager.h"
43#include "TGeoTrack.h"
44#include "TGeoOverlap.h"
45#include "TGeoChecker.h"
46#include "TGeoPhysicalNode.h"
47#include "TGeoPolygon.h"
48#include "TGeoCompositeShape.h"
49#include "TGeoShapeAssembly.h"
50#include "TGeoPainter.h"
51#include "TMath.h"
52
53#include "X3DBuffer.h"
54
55#include "TBuffer3D.h"
56#include "TBuffer3DTypes.h"
57#include "TVirtualViewer3D.h"
58#include "TVirtualX.h"
59
61
62////////////////////////////////////////////////////////////////////////////////
63/// Default constructor.
64
66{
68 if (manager) fGeoManager = manager;
69 else {
70 Error("ctor", "No geometry loaded");
71 return;
72 }
74 fNVisNodes = 0;
75 fBombX = 1.3;
76 fBombY = 1.3;
77 fBombZ = 1.3;
78 fBombR = 1.3;
82 fVisBranch = "";
83 fVolInfo = "";
88 fPlugin = nullptr;
89 fVisVolumes = new TObjArray();
90 fOverlap = nullptr;
91 fGlobal = new TGeoHMatrix();
92 fBuffer = new TBuffer3D(TBuffer3DTypes::kGeneric,20,3*20,0,0,0,0);
93 fClippingShape = nullptr;
94 fLastVolume = nullptr;
95 fTopVolume = nullptr;
97 memset(&fCheckedBox[0], 0, 6*sizeof(Double_t));
98
102 DefineColors();
103}
104////////////////////////////////////////////////////////////////////////////////
105/// Default destructor.
106
108{
109 if (fChecker) delete fChecker;
110 delete fVisVolumes;
111 delete fGlobal;
112 delete fBuffer;
113 if (fPlugin) delete fPlugin;
114}
115////////////////////////////////////////////////////////////////////////////////
116/// Add numpoints, numsegs, numpolys to the global 3D size.
117
118void TGeoPainter::AddSize3D(Int_t numpoints, Int_t numsegs, Int_t numpolys)
119{
120 gSize3D.numPoints += numpoints;
121 gSize3D.numSegs += numsegs;
122 gSize3D.numPolys += numpolys;
123}
124////////////////////////////////////////////////////////////////////////////////
125/// Create a primary TGeoTrack.
126
128{
129 return (TVirtualGeoTrack*)(new TGeoTrack(id,pdgcode,0,particle));
130}
131
132////////////////////////////////////////////////////////////////////////////////
133/// Average center of view of all painted tracklets and compute view box.
134
136{
137 static Int_t npoints = 0;
138 static Double_t xmin[3] = {0,0,0};
139 static Double_t xmax[3] = {0,0,0};
140 Int_t i;
141 if (reset) {
142 memset(box, 0, 6*sizeof(Double_t));
143 memset(xmin, 0, 3*sizeof(Double_t));
144 memset(xmax, 0, 3*sizeof(Double_t));
145 npoints = 0;
146 return;
147 }
148 if (npoints==0) {
149 for (i=0; i<3; i++) xmin[i]=xmax[i]=0;
150 npoints++;
151 }
152 npoints++;
153 Double_t ninv = 1./Double_t(npoints);
154 for (i=0; i<3; i++) {
155 box[i] += ninv*(point[i]-box[i]);
156 if (point[i]<xmin[i]) xmin[i]=point[i];
157 if (point[i]>xmax[i]) xmax[i]=point[i];
158 box[i+3] = 0.5*(xmax[i]-xmin[i]);
159 }
160}
161
162////////////////////////////////////////////////////////////////////////////////
163/// Get the new 'bombed' translation vector according current exploded view mode.
164
166{
167 memcpy(bombtr, tr, 3*sizeof(Double_t));
168 switch (fExplodedView) {
169 case kGeoNoBomb:
170 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:
187 return;
188 }
189}
190
191////////////////////////////////////////////////////////////////////////////////
192/// Check pushes and pulls needed to cross the next boundary with respect to the
193/// position given by FindNextBoundary. If radius is not mentioned the full bounding
194/// box will be sampled.
195
197{
198 fChecker->CheckBoundaryErrors(ntracks, radius);
199}
200
201////////////////////////////////////////////////////////////////////////////////
202/// Check the boundary errors reference file created by CheckBoundaryErrors method.
203/// The shape for which the crossing failed is drawn with the starting point in red
204/// and the extrapolated point to boundary (+/- failing push/pull) in yellow.
205
207{
209}
210
211////////////////////////////////////////////////////////////////////////////////
212/// Geometry checking method (see: TGeoManager::CheckGeometry())
213
214void TGeoPainter::CheckGeometryFull(Bool_t checkoverlaps, Bool_t checkcrossings, Int_t ntracks, const Double_t *vertex)
215{
216 fChecker->CheckGeometryFull(checkoverlaps,checkcrossings,ntracks,vertex);
217}
218
219////////////////////////////////////////////////////////////////////////////////
220/// Geometry checking method (see TGeoChecker).
221
222void TGeoPainter::CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
223{
224 fChecker->CheckGeometry(nrays, startx, starty, startz);
225}
226
227////////////////////////////////////////////////////////////////////////////////
228/// Check overlaps for the top volume of the geometry, within a limit OVLP.
229
230void TGeoPainter::CheckOverlaps(const TGeoVolume *vol, Double_t ovlp, Option_t *option) const
231{
232 fChecker->CheckOverlaps(vol, ovlp, option);
233}
234
235////////////////////////////////////////////////////////////////////////////////
236/// Check current point in the geometry.
237
239{
240 fChecker->CheckPoint(x,y,z,option);
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Test for shape navigation methods. Summary for test numbers:
245/// - 1: DistFromInside/Outside. Sample points inside the shape. Generate
246/// directions randomly in cos(theta). Compute DistFromInside and move the
247/// point with bigger distance. Compute DistFromOutside back from new point.
248/// Plot d-(d1+d2)
249
250void TGeoPainter::CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
251{
252 fChecker->CheckShape(shape, testNo, nsamples, option);
253}
254
255////////////////////////////////////////////////////////////////////////////////
256///Clear the list of visible volumes
257///reset the kVisOnScreen bit for volumes previously in the list
258
260{
261 if (!fVisVolumes) return;
262 TIter next(fVisVolumes);
263 TGeoVolume *vol;
264 while ((vol = (TGeoVolume*)next())) {
266 }
268}
269
270
271////////////////////////////////////////////////////////////////////////////////
272/// Define 100 colors with increasing light intensities for each basic color (1-7)
273/// Register these colors at indexes starting with 1000.
274
276{
277 static Int_t color = 0;
278 if (!color) {
280 for (auto icol=1; icol<10; ++icol)
281 color = GetColor(icol, 0.5);
282 }
283}
284
285////////////////////////////////////////////////////////////////////////////////
286/// Get index of a base color with given light intensity (0,1)
287
289{
290 using IntMap_t = std::map<Int_t, Int_t>;
291 constexpr Int_t ncolors = 100;
292 constexpr Float_t lmin = 0.25;
293 constexpr Float_t lmax = 0.75;
294 static IntMap_t colmap;
295 Int_t color = base;
296 // Search color in the map
297 auto it = colmap.find(base);
298 if (it != colmap.end()) return (it->second + light*(ncolors-1));
299 // Get color pointer if stored
300 TColor* col_base = gROOT->GetColor(base);
301 if (!col_base) {
302 // If color not defined, use gray palette
303 it = colmap.find(kBlack);
304 if (it != colmap.end()) return (it->second + light*(ncolors-1));
305 col_base = gROOT->GetColor(kBlack);
306 color = 1;
307 }
308 // Create a color palette for col_base
309 Float_t r=0., g=0., b=0., h=0., l=0., s=0.;
310 Double_t red[2], green[2], blue[2];
311 Double_t stop[] = {0., 1.0};
312
313 if (col_base) col_base->GetRGB(r,g,b);
314 TColor::RGB2HLS(r,g,b,h,l,s);
315 TColor::HLS2RGB(h,lmin,s,r,g,b);
316 red[0] = r;
317 green[0] = g;
318 blue[0] = b;
319 TColor::HLS2RGB(h,lmax,s,r,g,b);
320 red[1] = r;
321 green[1] = g;
322 blue[1] = b;
323 Int_t color_map_idx = TColor::CreateGradientColorTable(2, stop, red, green, blue, ncolors);
324 colmap[color] = color_map_idx;
325 return (color_map_idx + light*(ncolors-1));
326}
327
328////////////////////////////////////////////////////////////////////////////////
329/// Get currently drawn volume.
330
332{
333 if (!gPad) return nullptr;
334 return fTopVolume;
335}
336
337////////////////////////////////////////////////////////////////////////////////
338/// Compute the closest distance of approach from point px,py to a volume.
339
341{
342 const Int_t big = 9999;
343 const Int_t inaxis = 7;
344 const Int_t maxdist = 5;
345
346 if (fTopVolume != volume) fTopVolume = volume;
347 TView *view = gPad->GetView();
348 if (!view) return big;
349 TGeoBBox *box;
350 fGlobal->Clear();
352
353 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
354 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
355 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
356 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
357 // return if point not in user area
358 if (px < puxmin - inaxis) return big;
359 if (py > puymin + inaxis) return big;
360 if (px > puxmax + inaxis) return big;
361 if (py < puymax - inaxis) return big;
362
364 gPad->SetSelected(view);
365 Int_t dist = big;
366// Int_t id;
367
368 if (fPaintingOverlaps) {
369 TGeoVolume *crt;
370 crt = fOverlap->GetFirstVolume();
372 dist = crt->GetShape()->DistancetoPrimitive(px,py);
373 if (dist<maxdist) {
374 gPad->SetSelected(crt);
375 box = (TGeoBBox*)crt->GetShape();
376 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
377 fCheckedBox[3] = box->GetDX();
378 fCheckedBox[4] = box->GetDY();
379 fCheckedBox[5] = box->GetDZ();
380 return 0;
381 }
382 crt = fOverlap->GetSecondVolume();
384 dist = crt->GetShape()->DistancetoPrimitive(px,py);
385 if (dist<maxdist) {
386 gPad->SetSelected(crt);
387 box = (TGeoBBox*)crt->GetShape();
388 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
389 fCheckedBox[3] = box->GetDX();
390 fCheckedBox[4] = box->GetDY();
391 fCheckedBox[5] = box->GetDZ();
392 return 0;
393 }
394 return big;
395 }
396 // Compute distance to the right edge
397 if ((puxmax+inaxis-px) < 40) {
398 if ((py-puymax+inaxis) < 40) {
399 // when the mouse points to the (40x40) right corner of the pad, the manager class is selected
400 gPad->SetSelected(fGeoManager);
402 box = (TGeoBBox*)volume->GetShape();
403 memcpy(fCheckedBox, box->GetOrigin(), 3*sizeof(Double_t));
404 fCheckedBox[3] = box->GetDX();
405 fCheckedBox[4] = box->GetDY();
406 fCheckedBox[5] = box->GetDZ();
407 return 0;
408 }
409 // when the mouse points to the (40 pix) right edge of the pad, the top volume is selected
410 gPad->SetSelected(volume);
411 fVolInfo = volume->GetName();
412 box = (TGeoBBox*)volume->GetShape();
413 memcpy(fCheckedBox, box->GetOrigin(), 3*sizeof(Double_t));
414 fCheckedBox[3] = box->GetDX();
415 fCheckedBox[4] = box->GetDY();
416 fCheckedBox[5] = box->GetDZ();
417 return 0;
418 }
419
420 TGeoVolume *vol = volume;
421 Bool_t vis = vol->IsVisible();
422// Bool_t drawDaughters = kTRUE;
423 // Do we need to check a branch only?
424 if (volume->IsVisBranch()) {
425 if (!fGeoManager->IsClosed()) return big;
428 while (fGeoManager->GetLevel()) {
431 dist = vol->GetShape()->DistancetoPrimitive(px,py);
432 if (dist<maxdist) {
434 box = (TGeoBBox*)vol->GetShape();
435 fGeoManager->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
437 if (fGeoManager->IsNodeSelectable()) gPad->SetSelected(fCheckedNode);
438 else gPad->SetSelected(vol);
439 fCheckedBox[3] = box->GetDX();
440 fCheckedBox[4] = box->GetDY();
441 fCheckedBox[5] = box->GetDZ();
443 return 0;
444 }
445 fGeoManager->CdUp();
446 }
448 return dist;
449 }
450
451 // Do I need to look for the top volume ?
452 if ((fTopVisible && vis) || !vol->GetNdaughters() || !vol->IsVisDaughters() || vol->IsVisOnly()) {
453 dist = vol->GetShape()->DistancetoPrimitive(px,py);
454 if (dist<maxdist) {
455 fVolInfo = vol->GetName();
456 gPad->SetSelected(vol);
457 box = (TGeoBBox*)vol->GetShape();
458 memcpy(fCheckedBox, box->GetOrigin(), 3*sizeof(Double_t));
459 fCheckedBox[3] = box->GetDX();
460 fCheckedBox[4] = box->GetDY();
461 fCheckedBox[5] = box->GetDZ();
462 return 0;
463 }
464 if (vol->IsVisOnly() || !vol->GetNdaughters() || !vol->IsVisDaughters())
465 return dist;
466 }
467
468 // Iterate the volume content
469 TGeoIterator next(vol);
470 next.SetTopName(TString::Format("%s_1",vol->GetName()));
471 TGeoNode *daughter;
472
473 Int_t level, nd;
474 Bool_t last;
475
476 while ((daughter=next())) {
477 vol = daughter->GetVolume();
478 level = next.GetLevel();
479 nd = daughter->GetNdaughters();
480 vis = daughter->IsVisible();
481 if (volume->IsVisContainers()) {
482 if (vis && level<=fVisLevel) {
483 *fGlobal = next.GetCurrentMatrix();
484 dist = vol->GetShape()->DistancetoPrimitive(px,py);
485 if (dist<maxdist) {
486 next.GetPath(fVolInfo);
487 box = (TGeoBBox*)vol->GetShape();
488 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
489 fCheckedNode = daughter;
490 if (fGeoManager->IsNodeSelectable()) gPad->SetSelected(fCheckedNode);
491 else gPad->SetSelected(vol);
492 fCheckedBox[3] = box->GetDX();
493 fCheckedBox[4] = box->GetDY();
494 fCheckedBox[5] = box->GetDZ();
496 return 0;
497 }
498 }
499 // Check if we have to skip this branch
500 if (level==fVisLevel || !daughter->IsVisDaughters()) {
501 next.Skip();
502 continue;
503 }
504 } else if (volume->IsVisLeaves()) {
505 last = ((nd==0) || (level==fVisLevel) || (!daughter->IsVisDaughters()))?kTRUE:kFALSE;
506 if (vis && last) {
507 *fGlobal = next.GetCurrentMatrix();
508 dist = vol->GetShape()->DistancetoPrimitive(px,py);
509 if (dist<maxdist) {
510 next.GetPath(fVolInfo);
511 box = (TGeoBBox*)vol->GetShape();
512 fGlobal->LocalToMaster(box->GetOrigin(), &fCheckedBox[0]);
513 fCheckedNode = daughter;
514 if (fGeoManager->IsNodeSelectable()) gPad->SetSelected(fCheckedNode);
515 else gPad->SetSelected(vol);
516 fCheckedBox[3] = box->GetDX();
517 fCheckedBox[4] = box->GetDY();
518 fCheckedBox[5] = box->GetDZ();
520 return 0;
521 }
522 }
523 // Check if we have to skip the branch
524 if (last || !daughter->IsVisDaughters()) next.Skip();
525 }
526 }
527 return dist;
528}
529
530////////////////////////////////////////////////////////////////////////////////
531/// Set default angles for the current view.
532
534{
535 if (gPad) {
536 Int_t irep;
537 TView *view = gPad->GetView();
538 if (!view) return;
539 view->SetView(-206,126,75,irep);
540 ModifiedPad();
541 }
542}
543
544////////////////////////////////////////////////////////////////////////////////
545/// Set default volume colors according to tracking media
546
548{
550 TGeoVolume *vol;
551 while ((vol=(TGeoVolume*)next()))
553 ModifiedPad();
554}
555
556////////////////////////////////////////////////////////////////////////////////
557/// Count number of visible nodes down to a given level.
558
560{
561 TGeoVolume *vol = volume;
562 Int_t count = 0;
563 Bool_t vis = vol->IsVisible();
564 // Do I need to look for the top volume ?
565 if ((fTopVisible && vis) || !vol->GetNdaughters() || !vol->IsVisDaughters() || vol->IsVisOnly())
566 count++;
567 // Is this the only volume?
568 if (volume->IsVisOnly()) return count;
569
570 // Do we need to check a branch only?
571 if (volume->IsVisBranch()) {
574 count = fGeoManager->GetLevel() + 1;
576 return count;
577 }
578 // Iterate the volume content
579 TGeoIterator next(vol);
580 TGeoNode *daughter;
581 Int_t level, nd;
582 Bool_t last;
583
584 while ((daughter=next())) {
585// vol = daughter->GetVolume();
586 level = next.GetLevel();
587 nd = daughter->GetNdaughters();
588 vis = daughter->IsVisible();
589 if (volume->IsVisContainers()) {
590 if (vis && level<=rlevel) count++;
591 // Check if we have to skip this branch
592 if (level==rlevel || !daughter->IsVisDaughters()) {
593 next.Skip();
594 continue;
595 }
596 } else if (volume->IsVisLeaves()) {
597 last = ((nd==0) || (level==rlevel) || (!daughter->IsVisDaughters()))?kTRUE:kFALSE;
598 if (vis && last) count++;
599 // Check if we have to skip the branch
600 if (last) next.Skip();
601 }
602 }
603 return count;
604}
605
606////////////////////////////////////////////////////////////////////////////////
607/// Count total number of visible nodes.
608
610{
611 Int_t maxnodes = fGeoManager->GetMaxVisNodes();
612 Int_t vislevel = fGeoManager->GetVisLevel();
613// TGeoVolume *top = fGeoManager->GetTopVolume();
614 TGeoVolume *top = fTopVolume;
615 if (maxnodes <= 0 && top) {
616 fNVisNodes = CountNodes(top, vislevel);
617 SetVisLevel(vislevel);
618 return fNVisNodes;
619 }
620 //if (the total number of nodes of the top volume is less than maxnodes
621 // we can visualize everything.
622 //recompute the best visibility level
623 if (!top) {
624 SetVisLevel(vislevel);
625 return 0;
626 }
627 fNVisNodes = -1;
628 Bool_t again = kFALSE;
629 for (Int_t level = 1;level<20;level++) {
630 vislevel = level;
631 Int_t nnodes = CountNodes(top, level);
632 if (top->IsVisOnly() || top->IsVisBranch()) {
633 vislevel = fVisLevel;
634 fNVisNodes = nnodes;
635 break;
636 }
637 if (nnodes > maxnodes) {
638 vislevel--;
639 break;
640 }
641 if (nnodes == fNVisNodes) {
642 if (again) break;
643 again = kTRUE;
644 }
645 fNVisNodes = nnodes;
646 }
647 SetVisLevel(vislevel);
648 return fNVisNodes;
649}
650
651////////////////////////////////////////////////////////////////////////////////
652/// Check if Ged library is loaded and load geometry editor classe.
653
655{
656 if (fIsEditable) return;
657 if (!TClass::GetClass("TGedEditor")) return;
659 if ((h = gROOT->GetPluginManager()->FindHandler("TGeoManagerEditor"))) {
660 if (h->LoadPlugin() == -1) return;
661 h->ExecPlugin(0);
662 }
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// Start the geometry editor.
668
670{
671 if (!gPad) return;
672 if (!fIsEditable) {
673 if (!option[0]) gPad->GetCanvas()->GetCanvasImp()->ShowEditor();
675 CheckEdit();
676 }
677 gPad->SetSelected(fGeoManager);
678 gPad->GetCanvas()->Selected(gPad,fGeoManager,kButton1Down);
679}
680
681////////////////////////////////////////////////////////////////////////////////
682/// Draw method.
683
685{
687}
688
689////////////////////////////////////////////////////////////////////////////////
690/// Draw the time evolution of a radionuclide.
691
693{
694 Int_t ncoeff = sol->GetNcoeff();
695 if (!ncoeff) return;
696 Double_t tlo=0., thi=0.;
697 Double_t cn=0., lambda=0.;
698 Int_t i;
699 sol->GetRange(tlo, thi);
700 Bool_t autorange = (thi==0.)?kTRUE:kFALSE;
701
702 // Try to find the optimum range in time.
703 if (autorange) tlo = 0.;
704 sol->GetCoeff(0, cn, lambda);
705 Double_t lambdamin = lambda;
706 TString formula = "";
707 for (i=0; i<ncoeff; i++) {
708 sol->GetCoeff(i, cn, lambda);
709 formula += TString::Format("%g*exp(-%g*x)",cn, lambda);
710 if (i < ncoeff-1) formula += "+";
711 if (lambda < lambdamin &&
712 lambda > 0.) lambdamin = lambda;
713 }
714 if (autorange) thi = 10./lambdamin;
715 // Create a function
716 TF1 *func = new TF1(TString::Format("conc%s",sol->GetElement()->GetName()), formula.Data(), tlo,thi);
717 func->SetTitle(formula + ";time[s]" + TString::Format(";Concentration_of_%s",sol->GetElement()->GetName()));
718 func->SetMinimum(1.e-3);
719 func->SetMaximum(1.25*TMath::Max(sol->Concentration(tlo), sol->Concentration(thi)));
720 func->SetLineColor(sol->GetLineColor());
721 func->SetLineStyle(sol->GetLineStyle());
722 func->SetLineWidth(sol->GetLineWidth());
723 func->SetMarkerColor(sol->GetMarkerColor());
724 func->SetMarkerStyle(sol->GetMarkerStyle());
725 func->SetMarkerSize(sol->GetMarkerSize());
726 func->Draw(option);
727}
728
729////////////////////////////////////////////////////////////////////////////////
730/// Draw a polygon in 3D.
731
733{
734 Int_t nvert = poly->GetNvert();
735 if (!nvert) {
736 Error("DrawPolygon", "No vertices defined");
737 return;
738 }
739 Int_t nconv = poly->GetNconvex();
740 Double_t *x = new Double_t[nvert+1];
741 Double_t *y = new Double_t[nvert+1];
742 poly->GetVertices(x,y);
743 x[nvert] = x[0];
744 y[nvert] = y[0];
745 TGraph *g1 = new TGraph(nvert+1, x,y);
746 g1->SetTitle(Form("Polygon with %d vertices (outscribed %d)",nvert, nconv));
747 g1->SetLineColor(kRed);
748 g1->SetMarkerColor(kRed);
749 g1->SetMarkerStyle(4);
750 g1->SetMarkerSize(0.8);
751 delete [] x;
752 delete [] y;
753 Double_t *xc = 0;
754 Double_t *yc = 0;
755 TGraph *g2 = 0;
756 if (nconv && !poly->IsConvex()) {
757 xc = new Double_t[nconv+1];
758 yc = new Double_t[nconv+1];
759 poly->GetConvexVertices(xc,yc);
760 xc[nconv] = xc[0];
761 yc[nconv] = yc[0];
762 g2 = new TGraph(nconv+1, xc,yc);
763 g2->SetLineColor(kBlue);
764 g2->SetLineColor(kBlue);
766 g2->SetMarkerStyle(21);
767 g2->SetMarkerSize(0.4);
768 delete [] xc;
769 delete [] yc;
770 }
771 if (!gPad) {
772 gROOT->MakeDefCanvas();
773 }
774 g1->Draw("ALP");
775 if (g2) g2->Draw("LP");
776}
777
778////////////////////////////////////////////////////////////////////////////////
779/// Draw method.
780
782{
783 fTopVolume = vol;
784 fLastVolume = nullptr;
786// if (fVisOption==kGeoVisOnly ||
787// fVisOption==kGeoVisBranch) fGeoManager->SetVisOption(kGeoVisLeaves);
789 TString opt = option;
790 opt.ToLower();
792 fOverlap = nullptr;
793
794 if (fVisLock) {
797 }
798 Bool_t has_pad = (gPad==0)?kFALSE:kTRUE;
799 // Clear pad if option "same" not given
800 if (!gPad) {
801 gROOT->MakeDefCanvas();
802 }
803 if (!opt.Contains("same")) gPad->Clear();
804 // append this volume to pad
805 fTopVolume->AppendPad(option);
806
807 // Create a 3-D view
808 TView *view = gPad->GetView();
809 if (!view) {
810 view = TView::CreateView(11,0,0);
811 // Set the view to perform a first autorange (frame) draw.
812 // TViewer3DPad will revert view to normal painting after this
813 view->SetAutoRange(kTRUE);
814 if (has_pad) gPad->Update();
815 }
816 if (!opt.Contains("same")) Paint("range");
817 else Paint(opt);
818 view->SetAutoRange(kFALSE);
819 // If we are drawing into the pad, then the view needs to be
820 // set to perspective
821// if (!view->IsPerspective()) view->SetPerspective();
822
824
825 // Create a 3D viewer to paint us
826 gPad->GetViewer3D(option);
827}
828
829////////////////////////////////////////////////////////////////////////////////
830/// Draw a shape.
831
833{
834 TString opt = option;
835 opt.ToLower();
837 fOverlap = nullptr;
839
840 Bool_t has_pad = (gPad==0)?kFALSE:kTRUE;
841 // Clear pad if option "same" not given
842 if (!gPad) {
843 gROOT->MakeDefCanvas();
844 }
845 if (!opt.Contains("same")) gPad->Clear();
846 // append this shape to pad
847 shape->AppendPad(option);
848
849 // Create a 3-D view
850 TView *view = gPad->GetView();
851 if (!view) {
852 view = TView::CreateView(11,0,0);
853 // Set the view to perform a first autorange (frame) draw.
854 // TViewer3DPad will revert view to normal painting after this
855 view->SetAutoRange(kTRUE);
856 if (has_pad) gPad->Update();
857 }
858 PaintShape(shape,"range");
859 view->SetAutoRange(kFALSE);
860 view->SetPerspective();
861 // Create a 3D viewer to paint us
862 gPad->GetViewer3D(option);
863}
864
865////////////////////////////////////////////////////////////////////////////////
866/// Draw an overlap.
867
868void TGeoPainter::DrawOverlap(void *ovlp, Option_t *option)
869{
870 TString opt = option;
872 TGeoOverlap *overlap = (TGeoOverlap*)ovlp;
873 if (!overlap) return;
874
876 fOverlap = overlap;
877 opt.ToLower();
878 if (fVisLock) {
881 }
882 Bool_t has_pad = (gPad==0)?kFALSE:kTRUE;
883 // Clear pad if option "same" not given
884 if (!gPad) {
885 gROOT->MakeDefCanvas();
886 }
887 if (!opt.Contains("same")) gPad->Clear();
888 // append this volume to pad
889 overlap->AppendPad(option);
890
891 // Create a 3-D view
892 // Create a 3D viewer to paint us
893 gPad->GetViewer3D(option);
894 TView *view = gPad->GetView();
895 if (!view) {
896 view = TView::CreateView(11,0,0);
897 // Set the view to perform a first autorange (frame) draw.
898 // TViewer3DPad will revert view to normal painting after this
899 view->SetAutoRange(kTRUE);
900 PaintOverlap(ovlp, "range");
901 overlap->GetPolyMarker()->Draw("SAME");
902 if (has_pad) gPad->Update();
903 }
904
905 // If we are drawing into the pad, then the view needs to be
906 // set to perspective
907// if (!view->IsPerspective()) view->SetPerspective();
908 fVisLock = kTRUE;
909}
910
911
912////////////////////////////////////////////////////////////////////////////////
913/// Draw only one volume.
914
916{
917 TString opt = option;
918 opt.ToLower();
919 if (fVisLock) {
922 }
925 Bool_t has_pad = (gPad==0)?kFALSE:kTRUE;
926 // Clear pad if option "same" not given
927 if (!gPad) {
928 gROOT->MakeDefCanvas();
929 }
930 if (!opt.Contains("same")) gPad->Clear();
931 // append this volume to pad
933 fTopVolume->AppendPad(option);
934
935 // Create a 3-D view
936 TView *view = gPad->GetView();
937 if (!view) {
938 view = TView::CreateView(11,0,0);
939 // Set the view to perform a first autorange (frame) draw.
940 // TViewer3DPad will revert view to normal painting after this
941 view->SetAutoRange(kTRUE);
943 if (has_pad) gPad->Update();
944 }
945
946 // If we are drawing into the pad, then the view needs to be
947 // set to perspective
948// if (!view->IsPerspective()) view->SetPerspective();
949 fVisLock = kTRUE;
950}
951
952////////////////////////////////////////////////////////////////////////////////
953/// Draw current point in the same view.
954
956{
957 if (!gPad) return;
958 if (!gPad->GetView()) return;
959 TPolyMarker3D *pm = new TPolyMarker3D();
960 pm->SetMarkerColor(color);
961 const Double_t *point = fGeoManager->GetCurrentPoint();
962 pm->SetNextPoint(point[0], point[1], point[2]);
963 pm->SetMarkerStyle(8);
964 pm->SetMarkerSize(0.5);
965 pm->Draw("SAME");
966}
967
968////////////////////////////////////////////////////////////////////////////////
969
971{
972}
973
974////////////////////////////////////////////////////////////////////////////////
975/// Draw all volumes for a given path.
976
977void TGeoPainter::DrawPath(const char *path, Option_t *option)
978{
980 fVisBranch=path;
984 DrawVolume(fTopVolume,option);
985}
986
987////////////////////////////////////////////////////////////////////////////////
988/// Estimate camera movement between tmin and tmax for best track display
989
991{
992 if (!gPad) return;
993 TIter next(gPad->GetListOfPrimitives());
994 TVirtualGeoTrack *track;
995 TObject *obj;
996 Int_t ntracks = 0;
997 Double_t *point = 0;
998 AddTrackPoint(point, start, kTRUE);
999 while ((obj=next())) {
1000 if (strcmp(obj->ClassName(), "TGeoTrack")) continue;
1001 track = (TVirtualGeoTrack*)obj;
1002 ntracks++;
1003 track->PaintCollect(tmin, start);
1004 }
1005
1006 if (!ntracks) return;
1007 next.Reset();
1008 AddTrackPoint(point, end, kTRUE);
1009 while ((obj=next())) {
1010 if (strcmp(obj->ClassName(), "TGeoTrack")) continue;
1011 track = (TVirtualGeoTrack*)obj;
1012 if (!track) continue;
1013 track->PaintCollect(tmax, end);
1014 }
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018/// Execute mouse actions on a given volume.
1019
1021{
1022 if (!gPad) return;
1023 gPad->SetCursor(kPointer);
1024 switch (event) {
1025 case kButton1Down:
1026 if (!fIsEditable) CheckEdit();
1027 }
1028}
1029
1030////////////////////////////////////////////////////////////////////////////////
1031/// Execute mouse actions on a given shape.
1032
1034{
1035 if (!gPad) return;
1036 gPad->SetCursor(kHand);
1037 switch (event) {
1038 case kButton1Down:
1039 if (!fIsEditable) CheckEdit();
1040 }
1041}
1042
1043////////////////////////////////////////////////////////////////////////////////
1044/// Execute mouse actions on a given volume.
1045
1047{
1048 if (!gPad) return;
1049 if (!fIsEditable) CheckEdit();
1050// if (fIsRaytracing) return;
1051// Bool_t istop = (volume==fTopVolume)?kTRUE:kFALSE;
1052// if (istop) gPad->SetCursor(kHand);
1053// else gPad->SetCursor(kPointer);
1054 gPad->SetCursor(kHand);
1055// static Int_t width, color;
1056 switch (event) {
1057 case kMouseEnter:
1058// width = volume->GetLineWidth();
1059// color = volume->GetLineColor();
1060 break;
1061
1062 case kMouseLeave:
1063// volume->SetLineWidth(width);
1064// volume->SetLineColor(color);
1065 break;
1066
1067 case kButton1Down:
1068// volume->SetLineWidth(3);
1069// volume->SetLineColor(2);
1070// gPad->Modified();
1071// gPad->Update();
1072 break;
1073
1074 case kButton1Up:
1075// volume->SetLineWidth(width);
1076// volume->SetLineColor(color);
1077// gPad->Modified();
1078// gPad->Update();
1079 break;
1080
1081 case kButton1Double:
1082 gPad->SetCursor(kWatch);
1083 GrabFocus();
1084 break;
1085 }
1086}
1087
1088////////////////////////////////////////////////////////////////////////////////
1089/// Get some info about the current selected volume.
1090
1091const char *TGeoPainter::GetVolumeInfo(const TGeoVolume *volume, Int_t /*px*/, Int_t /*py*/) const
1092{
1093 static TString info;
1094 info = "";
1095 if (!gPad) return info;
1096 if (fPaintingOverlaps) {
1097 if (!fOverlap) {
1098 info = "wrong overlapping flag";
1099 return info;
1100 }
1101 TString ovtype, name;
1102 if (fOverlap->IsExtrusion()) ovtype="EXTRUSION";
1103 else ovtype = "OVERLAP";
1104 if (volume==fOverlap->GetFirstVolume()) name=volume->GetName();
1106 info = TString::Format("%s: %s of %g", name.Data(), ovtype.Data(), fOverlap->GetOverlap());
1107 return info;
1108 }
1109 else info = TString::Format("%s, shape=%s", fVolInfo.Data(), volume->GetShape()->ClassName());
1110 return info;
1111}
1112
1113////////////////////////////////////////////////////////////////////////////////
1114/// Create/return geometry checker.
1115
1117{
1119 return fChecker;
1120}
1121
1122////////////////////////////////////////////////////////////////////////////////
1123/// Get the current view angles.
1124
1125void TGeoPainter::GetViewAngles(Double_t &longitude, Double_t &latitude, Double_t &psi)
1126{
1127 if (!gPad) return;
1128 TView *view = gPad->GetView();
1129 if (!view) return;
1130 longitude = view->GetLongitude();
1131 latitude = view->GetLatitude();
1132 psi = view->GetPsi();
1133}
1134
1135////////////////////////////////////////////////////////////////////////////////
1136/// Move focus to current volume
1137
1139{
1140 if (!gPad) return;
1141 TView *view = gPad->GetView();
1142 if (!view) return;
1144 printf("Woops!!!\n");
1146 memcpy(&fCheckedBox[0], box->GetOrigin(), 3*sizeof(Double_t));
1147 fCheckedBox[3] = box->GetDX();
1148 fCheckedBox[4] = box->GetDY();
1149 fCheckedBox[5] = box->GetDZ();
1150 }
1151 view->SetPerspective();
1152 Int_t nvols = fVisVolumes->GetEntriesFast();
1153 Int_t nframes = nfr;
1154 if (nfr==0) {
1155 nframes = 1;
1156 if (nvols<1500) nframes=10;
1157 if (nvols<1000) nframes=20;
1158 if (nvols<200) nframes = 50;
1159 if (nvols<100) nframes = 100;
1160 }
1161 view->MoveFocus(&fCheckedBox[0], fCheckedBox[3], fCheckedBox[4], fCheckedBox[5], nframes, dlong, dlat, dpsi);
1162}
1163
1164////////////////////////////////////////////////////////////////////////////////
1165/// Generate a lego plot fot the top volume, according to option.
1166
1168 Int_t nphi, Double_t phimin, Double_t phimax,
1169 Double_t rmin, Double_t rmax, Option_t *option)
1170{
1171 return fChecker->LegoPlot(ntheta, themin, themax, nphi, phimin, phimax, rmin, rmax, option);
1172}
1173////////////////////////////////////////////////////////////////////////////////
1174/// Convert a local vector according view rotation matrix
1175
1176void TGeoPainter::LocalToMasterVect(const Double_t *local, Double_t *master) const
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) return;
1188 if (update) {
1189 gPad->Update();
1190 return;
1191 }
1192 TView *view = gPad->GetView();
1193 if (!view) return;
1194 view->SetViewChanged();
1195 gPad->Modified();
1196 if (gROOT->FromPopUp()) gPad->Update();
1197}
1198
1199////////////////////////////////////////////////////////////////////////////////
1200/// Paint current geometry according to option.
1201
1203{
1204 if (!fGeoManager || !fTopVolume) return;
1205 Bool_t is_padviewer = kTRUE;
1206 if (gPad) is_padviewer = (!strcmp(gPad->GetViewer3D()->ClassName(),"TViewer3DPad"))?kTRUE:kFALSE;
1207
1213
1214
1215 if (!fIsRaytracing || !is_padviewer) {
1216 if (fGeoManager->IsDrawingExtra()) {
1217 // loop the list of physical volumes
1218 fGeoManager->CdTop();
1220 Int_t nnodes = nodeList->GetEntriesFast();
1221 Int_t inode;
1222 TGeoPhysicalNode *node;
1223 for (inode=0; inode<nnodes; inode++) {
1224 node = (TGeoPhysicalNode*)nodeList->UncheckedAt(inode);
1225 PaintPhysicalNode(node, option);
1226 }
1227 } else {
1228 PaintVolume(fTopVolume,option);
1229 }
1230 fVisLock = kTRUE;
1231 }
1232 // Check if we have to raytrace (only in pad)
1233 if (fIsRaytracing && is_padviewer) Raytrace();
1234}
1235
1236////////////////////////////////////////////////////////////////////////////////
1237/// Paint an overlap.
1238
1239void TGeoPainter::PaintOverlap(void *ovlp, Option_t *option)
1240{
1241 if (!fGeoManager) return;
1242 TGeoOverlap *overlap = (TGeoOverlap *)ovlp;
1243 if (!overlap) return;
1244 Int_t color, transparency;
1245 if (fOverlap != overlap) fOverlap = overlap;
1247 TGeoHMatrix *hmat = fGlobal;
1248 TGeoVolume *vol;
1249 TGeoVolume *vol1 = overlap->GetFirstVolume();
1250 TGeoVolume *vol2 = overlap->GetSecondVolume();
1251 TGeoHMatrix *matrix1 = overlap->GetFirstMatrix();
1252 TGeoHMatrix *matrix2 = overlap->GetSecondMatrix();
1253 //
1254 vol = vol1;
1255 *hmat = matrix1;
1257 if (!fVisLock) fVisVolumes->Add(vol);
1259 color = vol->GetLineColor();
1260 transparency = vol->GetTransparency();
1261 vol->SetLineColor(kGreen);
1262 vol->SetTransparency(40);
1263 if (!strstr(option,"range")) ((TAttLine*)vol)->Modify();
1264 PaintShape(*(vol->GetShape()),option);
1265 vol->SetLineColor(color);
1266 vol->SetTransparency(transparency);
1267 vol = vol2;
1268 *hmat = matrix2;
1270 if (!fVisLock) fVisVolumes->Add(vol);
1272 color = vol->GetLineColor();
1273 transparency = vol->GetTransparency();
1274 vol->SetLineColor(kBlue);
1275 vol->SetTransparency(40);
1276 if (!strstr(option,"range")) ((TAttLine*)vol)->Modify();
1277 PaintShape(*(vol->GetShape()),option);
1278 vol->SetLineColor(color);
1279 vol->SetTransparency(transparency);
1281 fVisLock = kTRUE;
1282}
1283
1284////////////////////////////////////////////////////////////////////////////////
1285/// Paint recursively a node and its content according to visualization options.
1286
1288{
1289 PaintVolume(node->GetVolume(), option, global);
1290}
1291
1292////////////////////////////////////////////////////////////////////////////////
1293/// Paint recursively a node and its content according to visualization options.
1294
1296{
1297 if (fTopVolume != top) {
1299 fVisLock = kFALSE;
1300 }
1301 fTopVolume = top;
1302 if (!fVisLevel) return;
1303 TGeoVolume *vol = top;
1304 if(global)
1305 *fGlobal = *global;
1306 else
1307 fGlobal->Clear();
1309 Bool_t drawDaughters = kTRUE;
1310 Bool_t vis = (top->IsVisible() && !top->IsAssembly());
1311 Int_t transparency = 0;
1312
1313 // Update pad attributes in case we need to paint VOL
1314 if (!strstr(option,"range")) ((TAttLine*)vol)->Modify();
1315
1316 // Do we need to draw a branch ?
1317 if (top->IsVisBranch()) {
1320// while (fGeoManager->GetLevel()) {
1322 if (!fVisLock) {
1323 fVisVolumes->Add(vol);
1325 }
1327 transparency = vol->GetTransparency();
1328 vol->SetTransparency(40);
1329 if (!strstr(option,"range")) ((TAttLine*)vol)->Modify();
1330 if (global) {
1331 *fGlobal = *global;
1333 } else {
1335 }
1337 PaintShape(*(vol->GetShape()),option);
1338 vol->SetTransparency(transparency);
1339 fGeoManager->CdUp();
1340// }
1341 fVisLock = kTRUE;
1344 return;
1345 }
1346
1347 // Do I need to draw the top volume ?
1348 if ((fTopVisible && vis) || !top->GetNdaughters() || !top->IsVisDaughters() || top->IsVisOnly()) {
1351 PaintShape(*(vol->GetShape()),option);
1352 if (!fVisLock && !vol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1353 fVisVolumes->Add(vol);
1355 }
1356 if (top->IsVisOnly() || !top->GetNdaughters() || !top->IsVisDaughters()) {
1357 fVisLock = kTRUE;
1358 return;
1359 }
1360 }
1361
1362 // Iterate the volume content
1363 TGeoIterator next(vol);
1364 if (fPlugin) next.SetUserPlugin(fPlugin);
1365 TGeoNode *daughter;
1366// TGeoMatrix *glmat;
1367 Int_t level, nd;
1368 Bool_t last;
1369 Int_t line_color=0, line_width=0, line_style=0;
1370 while ((daughter=next())) {
1371 vol = daughter->GetVolume();
1373 level = next.GetLevel();
1374 nd = daughter->GetNdaughters();
1375 vis = daughter->IsVisible();
1376 drawDaughters = kTRUE;
1377 if (top->IsVisContainers()) {
1378 if (vis && level<=fVisLevel) {
1379 if (fPlugin) {
1380 line_color = vol->GetLineColor();
1381 line_width = vol->GetLineWidth();
1382 line_style = vol->GetLineStyle();
1383 transparency = vol->GetTransparency();
1385 }
1386 if (!strstr(option,"range")) ((TAttLine*)vol)->Modify();
1387 if (global) {
1388 *fGlobal = *global;
1389 *fGlobal *= *next.GetCurrentMatrix();
1390 } else {
1391 *fGlobal = next.GetCurrentMatrix();
1392 }
1394 drawDaughters = PaintShape(*(vol->GetShape()),option);
1395 if (fPlugin) {
1396 vol->SetLineColor(line_color);
1397 vol->SetLineWidth(line_width);
1398 vol->SetLineStyle(line_style);
1399 vol->SetTransparency(transparency);
1400 }
1401 if (!fVisLock && !daughter->IsOnScreen()) {
1402 fVisVolumes->Add(vol);
1404 }
1405 }
1406 // Check if we have to skip this branch
1407 if (!drawDaughters || level==fVisLevel || !daughter->IsVisDaughters()) {
1408 next.Skip();
1409 continue;
1410 }
1411 } else if (top->IsVisLeaves()) {
1412 last = ((nd==0) || (level==fVisLevel) || (!daughter->IsVisDaughters()))?kTRUE:kFALSE;
1413 if (vis && last) {
1414 if (fPlugin) {
1415 line_color = vol->GetLineColor();
1416 line_width = vol->GetLineWidth();
1417 line_style = vol->GetLineStyle();
1418 transparency = vol->GetTransparency();
1420 }
1421 if (!strstr(option,"range")) ((TAttLine*)vol)->Modify();
1422 if (global) {
1423 *fGlobal = *global;
1424 *fGlobal *= *next.GetCurrentMatrix();
1425 } else {
1426 *fGlobal = next.GetCurrentMatrix();
1427 }
1429 drawDaughters = PaintShape(*(vol->GetShape()),option);
1430 if (fPlugin) {
1431 vol->SetLineColor(line_color);
1432 vol->SetLineWidth(line_width);
1433 vol->SetLineStyle(line_style);
1434 vol->SetTransparency(transparency);
1435 }
1436 if (!fVisLock && !daughter->IsOnScreen()) {
1437 fVisVolumes->Add(vol);
1439 }
1440 }
1441 // Check if we have to skip the branch
1442 if (!drawDaughters || last || !daughter->IsVisDaughters()) next.Skip();
1443 }
1444 }
1445 if (fPlugin) fPlugin->SetIterator(0);
1447 fVisLock = kTRUE;
1448}
1449
1450////////////////////////////////////////////////////////////////////////////////
1451/// Paint the supplied shape into the current 3D viewer
1452
1454{
1455 Bool_t addDaughters = kTRUE;
1456
1457 TVirtualViewer3D *viewer = gPad->GetViewer3D();
1458
1459 if (!viewer || shape.IsA()==TGeoShapeAssembly::Class()) {
1460 return addDaughters;
1461 }
1462
1463 // For non-composite shapes we are the main paint method & perform the negotiation
1464 // with the viewer here
1465 if (!shape.IsComposite()) {
1466 // Does viewer prefer local frame positions?
1467 Bool_t localFrame = viewer->PreferLocalFrame();
1468 // Perform first fetch of buffer from the shape and try adding it
1469 // to the viewer
1470 const TBuffer3D & buffer =
1472 Int_t reqSections = viewer->AddObject(buffer, &addDaughters);
1473
1474 // If the viewer requires additional sections fetch from the shape (if possible)
1475 // and add again
1476 if (reqSections != TBuffer3D::kNone) {
1477 shape.GetBuffer3D(reqSections, localFrame);
1478 viewer->AddObject(buffer, &addDaughters);
1479 }
1480 }
1481 // Composite shapes have their own internal hierarchy of shapes, each
1482 // of which generate a filled TBuffer3D. Therefore we can't pass up a
1483 // single buffer to here. So as a special case the TGeoCompositeShape
1484 // performs it's own painting & negotiation with the viewer.
1485 else {
1486 const TGeoCompositeShape * composite = static_cast<const TGeoCompositeShape *>(&shape);
1487
1488 // We need the addDaughters flag returned from the viewer from paint
1489 // so can't use the normal TObject::Paint()
1490// TGeoHMatrix *matrix = (TGeoHMatrix*)TGeoShape::GetTransform();
1491// if (viewer->PreferLocalFrame()) matrix->Clear();
1492 addDaughters = composite->PaintComposite(option);
1493 }
1494
1495 return addDaughters;
1496}
1497
1498////////////////////////////////////////////////////////////////////////////////
1499/// Paint an overlap.
1500
1502{
1504 fGlobal->Clear();
1506 PaintShape(*shape,option);
1507}
1508
1509////////////////////////////////////////////////////////////////////////////////
1510/// Paints a physical node associated with a path.
1511
1513{
1514 if (!node->IsVisible()) return;
1515 Int_t level = node->GetLevel();
1516 Int_t i, col, wid, sty;
1517 TGeoShape *shape;
1519 TGeoHMatrix *matrix = fGlobal;
1520 TGeoVolume *vcrt;
1521 if (!node->IsVisibleFull()) {
1522 // Paint only last node in the branch
1523 vcrt = node->GetVolume();
1524 if (!strstr(option,"range")) ((TAttLine*)vcrt)->Modify();
1525 shape = vcrt->GetShape();
1526 *matrix = node->GetMatrix();
1529 if (!node->IsVolAttributes() && !strstr(option,"range")) {
1530 col = vcrt->GetLineColor();
1531 wid = vcrt->GetLineWidth();
1532 sty = vcrt->GetLineStyle();
1533 vcrt->SetLineColor(node->GetLineColor());
1534 vcrt->SetLineWidth(node->GetLineWidth());
1535 vcrt->SetLineStyle(node->GetLineStyle());
1536 ((TAttLine*)vcrt)->Modify();
1537 PaintShape(*shape,option);
1538 vcrt->SetLineColor(col);
1539 vcrt->SetLineWidth(wid);
1540 vcrt->SetLineStyle(sty);
1541 } else {
1542 PaintShape(*shape,option);
1543 }
1544 } else {
1545 // Paint full branch, except top node
1546 for (i=1;i<=level; i++) {
1547 vcrt = node->GetVolume(i);
1548 if (!strstr(option,"range")) ((TAttLine*)vcrt)->Modify();
1549 shape = vcrt->GetShape();
1550 *matrix = node->GetMatrix(i);
1553 if (!node->IsVolAttributes() && !strstr(option,"range")) {
1554 col = vcrt->GetLineColor();
1555 wid = vcrt->GetLineWidth();
1556 sty = vcrt->GetLineStyle();
1557 vcrt->SetLineColor(node->GetLineColor());
1558 vcrt->SetLineWidth(node->GetLineWidth());
1559 vcrt->SetLineStyle(node->GetLineStyle());
1560 ((TAttLine*)vcrt)->Modify();
1561 PaintShape(*shape,option);
1562 vcrt->SetLineColor(col);
1563 vcrt->SetLineWidth(wid);
1564 vcrt->SetLineStyle(sty);
1565 } else {
1566 PaintShape(*shape,option);
1567 }
1568 }
1569 }
1571}
1572
1573////////////////////////////////////////////////////////////////////////////////
1574/// Print overlaps (see TGeoChecker::PrintOverlaps())
1575
1577{
1579}
1580
1581////////////////////////////////////////////////////////////////////////////////
1582/// Text progress bar.
1583
1584void TGeoPainter::OpProgress(const char *opname, Long64_t current, Long64_t size, TStopwatch *watch, Bool_t last, Bool_t refresh, const char *msg)
1585{
1586 fChecker->OpProgress(opname,current,size,watch,last,refresh, msg);
1587}
1588
1589////////////////////////////////////////////////////////////////////////////////
1590/// Draw random points in the bounding box of a volume.
1591
1592void TGeoPainter::RandomPoints(const TGeoVolume *vol, Int_t npoints, Option_t *option)
1593{
1594 fChecker->RandomPoints((TGeoVolume*)vol, npoints, option);
1595}
1596
1597////////////////////////////////////////////////////////////////////////////////
1598/// Shoot nrays in the current drawn geometry
1599
1600void TGeoPainter::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol, Bool_t check_norm)
1601{
1602 fChecker->RandomRays(nrays, startx, starty, startz, target_vol, check_norm);
1603}
1604
1605////////////////////////////////////////////////////////////////////////////////
1606/// Raytrace current drawn geometry
1607
1609{
1610 if (!gPad || gPad->IsBatch()) return;
1611 TView *view = gPad->GetView();
1612 if (!view) return;
1613 Int_t rtMode = fGeoManager->GetRTmode();
1616 if (!view->IsPerspective()) view->SetPerspective();
1617 gVirtualX->SetMarkerSize(1);
1618 gVirtualX->SetMarkerStyle(1);
1619 Bool_t inclipst=kFALSE;
1620 Double_t krad = TMath::DegToRad();
1621 Double_t lat = view->GetLatitude();
1622 Double_t longit = view->GetLongitude();
1623 Double_t psi = view->GetPsi();
1624 Double_t c1 = TMath::Cos(psi*krad);
1625 Double_t s1 = TMath::Sin(psi*krad);
1626 Double_t c2 = TMath::Cos(lat*krad);
1627 Double_t s2 = TMath::Sin(lat*krad);
1628 Double_t s3 = TMath::Cos(longit*krad);
1629 Double_t c3 = -TMath::Sin(longit*krad);
1630 fMat[0] = c1*c3 - s1*c2*s3;
1631 fMat[1] = c1*s3 + s1*c2*c3;
1632 fMat[2] = s1*s2;
1633
1634 fMat[3] = -s1*c3 - c1*c2*s3;
1635 fMat[4] = -s1*s3 + c1*c2*c3;
1636 fMat[5] = c1*s2;
1637
1638 fMat[6] = s2*s3;
1639 fMat[7] = -s2*c3;
1640 fMat[8] = c2;
1641 Double_t u0, v0, du, dv;
1642 view->GetWindow(u0,v0,du,dv);
1643 Double_t dview = view->GetDview();
1644 Double_t dproj = view->GetDproj();
1645 Double_t local[3] = {0,0,1};
1646 Double_t dir[3], normal[3];
1647 LocalToMasterVect(local,dir);
1648 Double_t min[3], max[3];
1649 view->GetRange(min, max);
1650 Double_t cov[3];
1651 for (Int_t i=0; i<3; i++) cov[i] = 0.5*(min[i]+max[i]);
1652 Double_t cop[3];
1653 for (Int_t i=0; i<3; i++) cop[i] = cov[i] - dir[i]*dview;
1654 fGeoManager->InitTrack(cop, dir);
1655 Bool_t outside = fGeoManager->IsOutside();
1657 if (fClippingShape) inclipst = fClippingShape->Contains(cop);
1658 Int_t px, py;
1659 Double_t xloc, yloc, modloc;
1660 Int_t pxmin,pxmax, pymin,pymax;
1661 pxmin = gPad->UtoAbsPixel(0);
1662 pxmax = gPad->UtoAbsPixel(1);
1663 pymin = gPad->VtoAbsPixel(1);
1664 pymax = gPad->VtoAbsPixel(0);
1665 TGeoNode *next = nullptr;
1666 TGeoNode *nextnode = nullptr;
1667 Double_t step,steptot;
1668 Double_t *norm;
1669 const Double_t *point = fGeoManager->GetCurrentPoint();
1670 Double_t *ppoint = (Double_t*)point;
1671 Double_t tosource[3];
1672 Double_t calf;
1673 Double_t phi = 45.*krad;
1674 tosource[0] = -dir[0]*TMath::Cos(phi)+dir[1]*TMath::Sin(phi);
1675 tosource[1] = -dir[0]*TMath::Sin(phi)-dir[1]*TMath::Cos(phi);
1676 tosource[2] = -dir[2];
1677
1678 Bool_t done;
1679// Int_t istep;
1680 Int_t base_color, color;
1681 Double_t light;
1682 Double_t stemin=0, stemax=TGeoShape::Big();
1683 TPoint *pxy = new TPoint[1];
1684 TGeoVolume *nextvol;
1685 Int_t up;
1686 Int_t ntotal = pxmax*pymax;
1687 Int_t nrays = 0;
1688 TStopwatch *timer = new TStopwatch();
1689 timer->Start();
1690 for (px=pxmin; px<pxmax; px++) {
1691 for (py=pymin; py<pymax; py++) {
1692 if ((nrays%100)==0) OpProgress("Raytracing",nrays,ntotal,timer,kFALSE);
1693 nrays++;
1694 base_color = 1;
1695 steptot = 0;
1696 Bool_t inclip = inclipst;
1697 xloc = gPad->AbsPixeltoX(pxmin+pxmax-px);
1698 xloc = xloc*du-u0;
1699 yloc = gPad->AbsPixeltoY(pymin+pymax-py);
1700 yloc = yloc*dv-v0;
1701 modloc = TMath::Sqrt(xloc*xloc+yloc*yloc+dproj*dproj);
1702 local[0] = xloc/modloc;
1703 local[1] = yloc/modloc;
1704 local[2] = dproj/modloc;
1705 LocalToMasterVect(local,dir);
1707 fGeoManager->SetOutside(outside);
1710// fGeoManager->InitTrack(cop,dir);
1711 // current ray pointing to pixel (px,py)
1712 done = kFALSE;
1713 norm = 0;
1714 // propagate to the clipping shape if any
1715 if (fClippingShape) {
1716 if (inclip) {
1717 stemin = fClippingShape->DistFromInside(cop,dir,3);
1718 stemax = TGeoShape::Big();
1719 } else {
1720 stemax = fClippingShape->DistFromOutside(cop,dir,3);
1721 stemin = 0;
1722 }
1723 }
1724
1725 while (!done) {
1726 if (fClippingShape) {
1727 if (stemin>1E10) break;
1728 if (stemin>0) {
1729 // we are inside clipping shape
1730 fGeoManager->SetStep(stemin);
1731 next = fGeoManager->Step();
1732 steptot = 0;
1733 stemin = 0;
1734 up = 0;
1735 while (next) {
1736 // we found something after clipping region
1737 nextvol = next->GetVolume();
1738 if (nextvol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1739 done = kTRUE;
1740 base_color = nextvol->GetLineColor();
1741 fClippingShape->ComputeNormal(ppoint, dir, normal);
1742 norm = normal;
1743 break;
1744 }
1745 up++;
1746 next = fGeoManager->GetMother(up);
1747 }
1748 if (done) break;
1749 inclip = fClippingShape->Contains(ppoint);
1750 fGeoManager->SetStep(1E-3);
1751 while (inclip) {
1752 fGeoManager->Step();
1753 inclip = fClippingShape->Contains(ppoint);
1754 }
1755 stemax = fClippingShape->DistFromOutside(ppoint,dir,3);
1756 }
1757 }
1759 step = fGeoManager->GetStep();
1760 if (step>1E10) break;
1761 steptot += step;
1762 next = nextnode;
1763 // Check the step
1764 if (fClippingShape) {
1765 if (steptot>stemax) {
1766 steptot = 0;
1767 inclip = fClippingShape->Contains(ppoint);
1768 if (inclip) {
1769 stemin = fClippingShape->DistFromInside(ppoint,dir,3);
1770 stemax = TGeoShape::Big();
1771 continue;
1772 } else {
1773 stemin = 0;
1774 stemax = fClippingShape->DistFromOutside(ppoint,dir,3);
1775 }
1776 }
1777 }
1778 // Check if next node is visible
1779 if (!nextnode) continue;
1780 nextvol = nextnode->GetVolume();
1781 if (nextvol->TestAttBit(TGeoAtt::kVisOnScreen)) {
1782 done = kTRUE;
1783 base_color = nextvol->GetLineColor();
1784 next = nextnode;
1785 break;
1786 }
1787 }
1788 if (!done) continue;
1789 // current ray intersect a visible volume having color=base_color
1790 if (rtMode > 0) {
1793 for (Int_t i=0; i<3; ++i) local[i] += 1.E-8*dir[i];
1794 step = next->GetVolume()->GetShape()->DistFromInside(local,dir,3);
1795 for (Int_t i=0; i<3; ++i) local[i] += step*dir[i];
1796 next->GetVolume()->GetShape()->ComputeNormal(local, dir, normal);
1797 norm = normal;
1798 } else {
1799 if (!norm) norm = fGeoManager->FindNormalFast();
1800 if (!norm) continue;
1801 }
1802 calf = norm[0]*tosource[0]+norm[1]*tosource[1]+norm[2]*tosource[2];
1803 light = TMath::Abs(calf);
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 fChecker->OpProgress("Raytracing",nrays,ntotal,timer,kTRUE);
1815 delete timer;
1816}
1817
1818////////////////////////////////////////////////////////////////////////////////
1819/// Shoot npoints randomly in a box of 1E-5 around current point.
1820/// Return minimum distance to points outside.
1821
1823 const char* g3path)
1824{
1825 return fChecker->SamplePoints(npoints, dist, epsil, g3path);
1826}
1827
1828////////////////////////////////////////////////////////////////////////////////
1829/// Set cartesian and radial bomb factors for translations.
1830
1832{
1833 fBombX = bombx;
1834 fBombY = bomby;
1835 fBombZ = bombz;
1836 fBombR = bombr;
1837 if (IsExplodedView()) ModifiedPad();
1838}
1839
1840////////////////////////////////////////////////////////////////////////////////
1841/// Set type of exploding view.
1842
1844{
1845 if ((ibomb<0) || (ibomb>3)) {
1846 Warning("SetExplodedView", "exploded view can be 0-3");
1847 return;
1848 }
1849 if ((Int_t)ibomb==fExplodedView) return;
1850 Bool_t change = (gPad==0)?kFALSE:kTRUE;
1851
1852 if (ibomb==kGeoNoBomb) {
1853 change &= ((fExplodedView==kGeoNoBomb)?kFALSE:kTRUE);
1854 }
1855 if (ibomb==kGeoBombXYZ) {
1856 change &= ((fExplodedView==kGeoBombXYZ)?kFALSE:kTRUE);
1857 }
1858 if (ibomb==kGeoBombCyl) {
1859 change &= ((fExplodedView==kGeoBombCyl)?kFALSE:kTRUE);
1860 }
1861 if (ibomb==kGeoBombSph) {
1862 change &= ((fExplodedView==kGeoBombSph)?kFALSE:kTRUE);
1863 }
1864 fExplodedView = ibomb;
1865 if (change) ModifiedPad();
1866}
1867
1868////////////////////////////////////////////////////////////////////////////////
1869/// Set number of segments to approximate circles.
1870
1872{
1873 if (nseg<3) {
1874 Warning("SetNsegments", "number of segments should be > 2");
1875 return;
1876 }
1877 if (fNsegments==nseg) return;
1878 fNsegments = nseg;
1879 ModifiedPad();
1880}
1881
1882////////////////////////////////////////////////////////////////////////////////
1883/// Set number of points to be generated on the shape outline when checking for overlaps.
1884
1886 fChecker->SetNmeshPoints(npoints);
1887}
1888
1889////////////////////////////////////////////////////////////////////////////////
1890/// Select a node to be checked for overlaps. All overlaps not involving it will
1891/// be ignored.
1892
1895}
1896
1897////////////////////////////////////////////////////////////////////////////////
1898/// Set default level down to which visualization is performed
1899
1901 if (level==fVisLevel && fLastVolume==fTopVolume) return;
1902 fVisLevel=level;
1903 if (!fTopVolume) return;
1904 if (fVisLock) {
1906 fVisLock = kFALSE;
1907 }
1908 if (!fLastVolume) {
1909// printf("--- Drawing %6d nodes with %d visible levels\n",fNVisNodes,fVisLevel);
1910 return;
1911 }
1912 if (!gPad) return;
1913 if (gPad->GetView()) {
1914// printf("--- Drawing %6d nodes with %d visible levels\n",fNVisNodes,fVisLevel);
1915 ModifiedPad();
1916 }
1917}
1918
1919////////////////////////////////////////////////////////////////////////////////
1920/// Set top geometry volume as visible.
1921
1923{
1924 if (fTopVisible==vis) return;
1925 fTopVisible = vis;
1926 ModifiedPad();
1927}
1928
1929////////////////////////////////////////////////////////////////////////////////
1930/// Set drawing mode :
1931/// - option=0 (default) all nodes drawn down to vislevel
1932/// - option=1 leaves and nodes at vislevel drawn
1933/// - option=2 path is drawn
1934
1936 if ((fVisOption<0) || (fVisOption>4)) {
1937 Warning("SetVisOption", "wrong visualization option");
1938 return;
1939 }
1940
1941 if (option == kGeoVisChanged) {
1942 if (fVisLock) {
1944 fVisLock = kFALSE;
1945 }
1946 ModifiedPad();
1947 return;
1948 }
1949
1950 if (fTopVolume) {
1951 TGeoAtt *att = (TGeoAtt*)fTopVolume;
1955 switch (option) {
1956 case kGeoVisDefault:
1958 break;
1959 case kGeoVisLeaves:
1960 break;
1961 case kGeoVisOnly:
1963 break;
1964 }
1965 }
1966
1967 if (fVisOption==option) return;
1968 fVisOption=option;
1969 if (fVisLock) {
1971 fVisLock = kFALSE;
1972 }
1973 ModifiedPad();
1974}
1975
1976////////////////////////////////////////////////////////////////////////////////
1977/// Returns distance between point px,py on the pad an a shape.
1978
1980{
1981 const Int_t inaxis = 7;
1982 const Int_t maxdist = 5;
1983 const Int_t big = 9999;
1984 Int_t dist = big;
1985 if (!gPad) return dist;
1986 TView *view = gPad->GetView();
1987 if (!(numpoints && view)) return dist;
1988 if (shape->IsA()==TGeoShapeAssembly::Class()) return dist;
1989
1990 if (fIsPaintingShape) {
1991 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
1992 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
1993 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
1994 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
1995 // return if point not in user area
1996 if (px < puxmin - inaxis) return big;
1997 if (py > puymin + inaxis) return big;
1998 if (px > puxmax + inaxis) return big;
1999 if (py < puymax - inaxis) return big;
2000 if ((puxmax+inaxis-px) < 40) {
2001 // when the mouse points to the (40 pix) right edge of the pad, the manager class is selected
2002 gPad->SetSelected(fGeoManager);
2003 return 0;
2004 }
2005 }
2006
2007 fBuffer->SetRawSizes(numpoints, 3*numpoints, 0, 0, 0, 0);
2009 shape->SetPoints(points);
2010 Double_t dpoint2, x1, y1, xndc[3];
2011 Double_t dmaster[3];
2012 Int_t j;
2013 for (Int_t i=0; i<numpoints; i++) {
2014 j = 3*i;
2016 points[j]=dmaster[0]; points[j+1]=dmaster[1]; points[j+2]=dmaster[2];
2017 view->WCtoNDC(&points[j], xndc);
2018 x1 = gPad->XtoAbsPixel(xndc[0]);
2019 y1 = gPad->YtoAbsPixel(xndc[1]);
2020 dpoint2 = (px-x1)*(px-x1) + (py-y1)*(py-y1);
2021 if (dpoint2 < dist) dist=(Int_t)dpoint2;
2022 }
2023 if (dist > 100) return dist;
2024 dist = Int_t(TMath::Sqrt(Double_t(dist)));
2025 if (dist<maxdist && fIsPaintingShape) gPad->SetSelected((TObject*)shape);
2026 return dist;
2027}
2028
2029////////////////////////////////////////////////////////////////////////////////
2030/// Check time of finding "Where am I" for n points.
2031
2032void TGeoPainter::Test(Int_t npoints, Option_t *option)
2033{
2034 fChecker->Test(npoints, option);
2035}
2036
2037////////////////////////////////////////////////////////////////////////////////
2038/// Geometry overlap checker based on sampling.
2039
2040void TGeoPainter::TestOverlaps(const char* path)
2041{
2042 fChecker->TestOverlaps(path);
2043}
2044
2045////////////////////////////////////////////////////////////////////////////////
2046/// Check voxels efficiency per volume.
2047
2049{
2050 return fChecker->TestVoxels(vol);
2051}
2052
2053////////////////////////////////////////////////////////////////////////////////
2054/// Get the new 'unbombed' translation vector according current exploded view mode.
2055
2057{
2058 memcpy(bombtr, tr, 3*sizeof(Double_t));
2059 switch (fExplodedView) {
2060 case kGeoNoBomb:
2061 return;
2062 case kGeoBombXYZ:
2063 bombtr[0] /= fBombX;
2064 bombtr[1] /= fBombY;
2065 bombtr[2] /= fBombZ;
2066 return;
2067 case kGeoBombCyl:
2068 bombtr[0] /= fBombR;
2069 bombtr[1] /= fBombR;
2070 bombtr[2] /= fBombZ;
2071 return;
2072 case kGeoBombSph:
2073 bombtr[0] /= fBombR;
2074 bombtr[1] /= fBombR;
2075 bombtr[2] /= fBombR;
2076 return;
2077 default:
2078 return;
2079 }
2080}
2081
2082////////////////////////////////////////////////////////////////////////////////
2083/// Compute weight [kg] of the current volume.
2084
2086{
2087 return fChecker->Weight(precision, option);
2088}
2089
2090
@ 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
ROOT::R::TRInterface & r
Definition Object.C:4
#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)
static const double x1[5]
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
long long Long64_t
Definition RtypesCore.h:80
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
@ kRed
Definition Rtypes.h:66
@ kBlack
Definition Rtypes.h:65
@ kGreen
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float xmax
#define gROOT
Definition TROOT.h:404
char * Form(const char *fmt,...)
#define gPad
#define gVirtualX
Definition TVirtualX.h:338
point * points
Definition X3DBuffer.c:22
#define gSize3D
Definition X3DBuffer.h:40
Line Attributes class.
Definition TAttLine.h:18
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:35
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:34
virtual void Modify()
Change current line attributes if necessary.
Definition TAttLine.cxx:245
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition TAttMarker.h:32
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:31
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition TAttMarker.h:33
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:41
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:112
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:2966
The color creation and management class.
Definition TColor.h:19
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:1510
virtual void GetRGB(Float_t &r, Float_t &g, Float_t &b) const
Definition TColor.h:52
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1822
static void InitializeColors()
Initialize colors used by the TCanvas based graphics (via TColor objects).
Definition TColor.cxx:1142
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:1650
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.)
Static function creating a color table with several connected linear gradients.
Definition TColor.cxx:2352
1-Dim function class
Definition TF1.h:213
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:3414
virtual void SetTitle(const char *title="")
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition TF1.cxx:3569
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:3427
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition TF1.cxx:1328
Visualization and tracking attributes for volumes and nodes.
Definition TGeoAtt.h:18
Bool_t TestAttBit(UInt_t f) const
Definition TGeoAtt.h:68
Bool_t IsVisBranch() const
Definition TGeoAtt.h:90
@ kVisContainers
Definition TGeoAtt.h:33
@ kVisOnly
Definition TGeoAtt.h:34
@ kVisOnScreen
Definition TGeoAtt.h:32
@ kVisBranch
Definition TGeoAtt.h:35
void ResetAttBit(UInt_t f)
Definition TGeoAtt.h:67
void SetVisRaytrace(Bool_t flag=kTRUE)
Definition TGeoAtt.h:70
Bool_t IsVisDaughters() const
Definition TGeoAtt.h:89
void SetAttBit(UInt_t f)
Definition TGeoAtt.h:65
Box class.
Definition TGeoBBox.h:18
Double_t Concentration(Double_t time) const
Find concentration of the element at a given time.
Int_t GetNcoeff() const
void GetCoeff(Int_t i, Double_t &cn, Double_t &lambda) const
void GetRange(Double_t &tmin, Double_t &tmax) const
TGeoElementRN * GetElement() const
Geometry checking package.
Definition TGeoChecker.h:38
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
Test for shape navigation methods.
void CheckOverlaps(const TGeoVolume *vol, Double_t ovlp=0.1, Option_t *option="")
Check illegal overlaps for volume VOL within a limit OVLP.
void TestOverlaps(const char *path)
Geometry overlap checker based on sampling.
void SetSelectedNode(TGeoNode *node)
Definition TGeoChecker.h:88
TGeoNode * SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil, const char *g3path)
shoot npoints randomly in a box of 1E-5 around current point.
Double_t Weight(Double_t precision=0.01, Option_t *option="v")
Estimate weight of top level volume with a precision SIGMA(W)/W better than PRECISION.
void CheckGeometryFull(Bool_t checkoverlaps=kTRUE, Bool_t checkcrossings=kTRUE, Int_t nrays=10000, const Double_t *vertex=NULL)
Geometry checking.
void CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
Shoot nrays with random directions from starting point (startx, starty, startz) in the reference fram...
virtual void CheckBoundaryReference(Int_t icheck=-1)
Check the boundary errors reference file created by CheckBoundaryErrors method.
void CheckPoint(Double_t x=0, Double_t y=0, Double_t z=0, Option_t *option="")
Draw point (x,y,z) over the picture of the daughters of the volume containing this point.
void SetNmeshPoints(Int_t npoints=1000)
Set number of points to be generated on the shape outline when checking for overlaps.
void PrintOverlaps() const
Print the current list of overlaps held by the manager class.
void Test(Int_t npoints, Option_t *option)
Check time of finding "Where am I" for n points.
Bool_t TestVoxels(TGeoVolume *vol, Int_t npoints=1000000)
Returns optimal voxelization type for volume vol.
void OpProgress(const char *opname, Long64_t current, Long64_t size, TStopwatch *watch=0, Bool_t last=kFALSE, Bool_t refresh=kFALSE, const char *msg="")
Print current operation progress.
void RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol=0, Bool_t check_norm=kFALSE)
Randomly shoot nrays from point (startx,starty,startz) and plot intersections with surfaces for curre...
TH2F * LegoPlot(Int_t ntheta=60, Double_t themin=0., Double_t themax=180., Int_t nphi=90, Double_t phimin=0., Double_t phimax=360., Double_t rmin=0., Double_t rmax=9999999, Option_t *option="")
Generate a lego plot fot the top volume, according to option.
void RandomPoints(TGeoVolume *vol, Int_t npoints, Option_t *option)
Draw random points in the bounding box of a volume.
virtual void CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.)
Check pushes and pulls needed to cross the next boundary with respect to the position given by FindNe...
Composite shapes are Boolean combinations of two or more shape components.
virtual Bool_t PaintComposite(Option_t *option="") const
Paint this composite shape into the current 3D viewer Returns bool flag indicating if the caller shou...
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:421
void Clear(Option_t *option="")
clear the data for this matrix
void SetIterator(const TGeoIterator *iter)
Definition TGeoNode.h:233
virtual void ProcessNode()=0
A geometry iterator.
Definition TGeoNode.h:245
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:276
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
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:41
Bool_t IsReflection() const
Definition TGeoMatrix.h:69
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:41
Bool_t IsVisDaughters() const
Definition TGeoNode.h:108
Bool_t IsOnScreen() const
check if this node is drawn. Assumes that this node is current
Definition TGeoNode.cxx:273
TGeoVolume * GetVolume() const
Definition TGeoNode.h:97
Int_t GetNdaughters() const
Definition TGeoNode.h:93
Bool_t IsVisible() const
Definition TGeoNode.h:107
Base class describing geometry overlaps.
Definition TGeoOverlap.h:41
TGeoVolume * GetSecondVolume() const
Definition TGeoOverlap.h:74
TPolyMarker3D * GetPolyMarker() const
Definition TGeoOverlap.h:72
TGeoHMatrix * GetFirstMatrix() const
Definition TGeoOverlap.h:75
Bool_t IsExtrusion() const
Definition TGeoOverlap.h:78
Double_t GetOverlap() const
Definition TGeoOverlap.h:77
TGeoHMatrix * GetSecondMatrix() const
Definition TGeoOverlap.h:76
TGeoVolume * GetFirstVolume() const
Definition TGeoOverlap.h:73
Class implementing all draw interfaces for a generic 3D viewer using TBuffer3D mechanism.
Definition TGeoPainter.h:40
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:44
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:70
TGeoIteratorPlugin * fPlugin
Definition TGeoPainter.h:69
Double_t fMat[9]
Definition TGeoPainter.h:47
TGeoNode * SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil, const char *g3path) override
Shoot npoints randomly in a box of 1E-5 around current point.
void PrintOverlaps() const override
Print overlaps (see TGeoChecker::PrintOverlaps())
TBuffer3D * fBuffer
Definition TGeoPainter.h:63
void RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol=nullptr, Bool_t check_norm=kFALSE) override
Shoot nrays in the current drawn geometry.
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:64
TH2F * LegoPlot(Int_t ntheta=60, Double_t themin=0., Double_t themax=180., Int_t nphi=90, Double_t phimin=0., Double_t phimax=360., Double_t rmin=0., Double_t rmax=9999999, Option_t *option="") override
Generate a lego plot fot the top volume, according to option.
void Raytrace(Option_t *option="") override
Raytrace current drawn geometry.
TString fVisBranch
Definition TGeoPainter.h:58
TGeoVolume * fTopVolume
Definition TGeoPainter.h:67
TGeoChecker * GetChecker()
Create/return geometry checker.
virtual ~TGeoPainter()
Default destructor.
void SetVisOption(Int_t option=0) override
Set drawing mode :
Bool_t IsExplodedView() const override
Bool_t fIsPaintingShape
Definition TGeoPainter.h:57
Double_t fBombR
Definition TGeoPainter.h:45
void TestOverlaps(const char *path) override
Geometry overlap checker based on sampling.
void CheckGeometryFull(Bool_t checkoverlaps=kTRUE, Bool_t checkcrossings=kTRUE, Int_t nrays=10000, const Double_t *vertex=nullptr) override
Geometry checking method (see: TGeoManager::CheckGeometry())
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:68
void SetNmeshPoints(Int_t npoints) override
Set number of points to be generated on the shape outline when checking for overlaps.
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:53
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:66
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:71
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:46
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:61
void CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const override
Geometry checking method (see TGeoChecker).
Double_t fBombY
Definition TGeoPainter.h:43
Bool_t fPaintingOverlaps
Definition TGeoPainter.h:55
void CheckOverlaps(const TGeoVolume *vol, Double_t ovlp=0.1, Option_t *option="") const override
Check overlaps for the top volume of the geometry, within a limit OVLP.
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 CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.) override
Check pushes and pulls needed to cross the next boundary with respect to the position given by FindNe...
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:48
void Test(Int_t npoints, Option_t *option) override
Check time of finding "Where am I" for n points.
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:59
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option) override
Test for shape navigation methods.
void RandomPoints(const TGeoVolume *vol, Int_t npoints, Option_t *option="") override
Draw random points in the bounding box of a volume.
Bool_t fTopVisible
Definition TGeoPainter.h:54
TGeoHMatrix * fGlobal
Definition TGeoPainter.h:62
Int_t fVisLevel
Definition TGeoPainter.h:50
void SetNsegments(Int_t nseg=20) override
Set number of segments to approximate circles.
void CheckBoundaryReference(Int_t icheck=-1) override
Check the boundary errors reference file created by CheckBoundaryErrors method.
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:60
void PaintPhysicalNode(TGeoPhysicalNode *node, Option_t *option="")
Paints a physical node associated with a path.
Int_t fNVisNodes
Definition TGeoPainter.h:49
Int_t fExplodedView
Definition TGeoPainter.h:52
const char * GetVolumeInfo(const TGeoVolume *volume, Int_t px, Int_t py) const override
Get some info about the current selected volume.
void CheckPoint(Double_t x=0, Double_t y=0, Double_t z=0, Option_t *option="") override
Check current point in the geometry.
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 OpProgress(const char *opname, Long64_t current, Long64_t size, TStopwatch *watch=nullptr, Bool_t last=kFALSE, Bool_t refresh=kFALSE, const char *msg="") override
Text progress bar.
void SetCheckedNode(TGeoNode *node) override
Select a node to be checked for overlaps.
void DrawBatemanSol(TGeoBatemanSol *sol, Option_t *option="") override
Draw the time evolution of a radionuclide.
TGeoChecker * fChecker
Definition TGeoPainter.h:65
Bool_t fIsRaytracing
Definition TGeoPainter.h:56
void ExecuteShapeEvent(TGeoShape *shape, Int_t event, Int_t px, Int_t py) override
Execute mouse actions on a given shape.
Double_t Weight(Double_t precision, Option_t *option="v") override
Compute weight [kg] of the current volume.
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:42
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:51
Bool_t TestVoxels(TGeoVolume *vol) override
Check voxels efficiency per volume.
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:20
Bool_t IsConvex() const
Definition TGeoPolygon.h:62
Int_t GetNconvex() const
Definition TGeoPolygon.h:56
void GetVertices(Double_t *x, Double_t *y) const
Fill list of vertices into provided arrays.
Int_t GetNvert() const
Definition TGeoPolygon.h:55
void GetConvexVertices(Double_t *x, Double_t *y) const
Fill list of vertices of the convex outscribed polygon into provided arrays.
Base abstract class for all shapes.
Definition TGeoShape.h:26
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)=0
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)=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:88
static void SetTransform(TGeoMatrix *matrix)
Set current transformation matrix that applies to shape.
virtual Bool_t IsComposite() const
Definition TGeoShape.h:130
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const =0
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) 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
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
Bool_t IsVisContainers() const
Definition TGeoVolume.h:156
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:173
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
void SetTransparency(Char_t transparency=0)
Definition TGeoVolume.h:218
TGeoShape * GetShape() const
Definition TGeoVolume.h:189
Bool_t IsRaytracing() const
Check if the painter is currently ray-tracing the content of this volume.
Bool_t IsVisLeaves() const
Definition TGeoVolume.h:157
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Bool_t IsVisOnly() const
Definition TGeoVolume.h:158
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
Char_t GetTransparency() const
Definition TGeoVolume.h:188
virtual Bool_t IsVisible() const
Definition TGeoVolume.h:154
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual void SetTitle(const char *title="")
Change (i.e.
Definition TGraph.cxx:2353
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition TGraph.cxx:769
2-D histogram with a float per channel (see TH1 documentation)}
Definition TH2.h:251
void Reset()
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
void Add(TObject *obj)
Definition TObjArray.h:68
virtual void Clear(Option_t *option="")
Remove all objects from the array.
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
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:200
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:949
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:177
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
Longptr_t ExecPlugin(int nargs, const T &... params)
SCoord_t fY
Definition TPoint.h:36
SCoord_t fX
Definition TPoint.h:35
A 3D polymarker.
virtual Int_t SetNextPoint(Double_t x, Double_t y, Double_t z)
Set point following LastPoint to x, y, z.
virtual void Draw(Option_t *option="")
Draws 3-D polymarker with its current attributes.
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Stop()
Stop the stopwatch.
Basic string class.
Definition TString.h:136
void ToLower()
Change string to lower-case.
Definition TString.cxx:1150
const char * Data() const
Definition TString.h:369
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:2336
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
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
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
static TView * CreateView(Int_t system=1, const Double_t *rmin=0, const Double_t *rmax=0)
Create a concrete default 3-d view via the plug-in manager.
Definition TView.cxx:27
virtual 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.
virtual void PaintCollect(Double_t, Double_t *)
static void ShowEditor()
Show the global pad editor. Static method.
Abstract 3D shapes viewer.
virtual Bool_t PreferLocalFrame() const =0
virtual Int_t AddObject(const TBuffer3D &buffer, Bool_t *addChildren=0)=0
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
return c1
Definition legend1.C:41
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
return c2
Definition legend2.C:14
return c3
Definition legend3.C:15
Short_t Max(Short_t a, Short_t b)
Definition TMathBase.h:208
constexpr Double_t DegToRad()
Conversion from degree to radian:
Definition TMath.h:81
Double_t Sqrt(Double_t x)
Definition TMath.h:641
Double_t Cos(Double_t)
Definition TMath.h:593
Double_t Sin(Double_t)
Definition TMath.h:589
Short_t Abs(Short_t d)
Definition TMathBase.h:120
auto * l
Definition textangle.C:4
REAL * vertex
Definition triangle.c:513