Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoShape.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 31/01/02
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoShape
13\ingroup Shapes_classes
14Base abstract class for all shapes.
15
16 Shapes are geometrical objects that provide the basic modelling
17functionality. They provide the definition of the LOCAL frame of coordinates,
18with respect to which they are defined. Any implementation of a shape deriving
19from the base TGeoShape class has to provide methods for :
20
21 - finding out if a point defined in their local frame is or not contained
22 inside;
23 - computing the distance from a local point to getting outside/entering the
24 shape, given a known direction;
25 - computing the maximum distance in any direction from a local point that
26 does NOT result in a boundary crossing of the shape (safe distance);
27 - computing the cosines of the normal vector to the crossed shape surface,
28 given a starting local point and an ongoing direction.
29 All the features above are globally managed by the modeller in order to
30 provide navigation functionality. In addition to those, shapes have also to
31 implement additional specific abstract methods :
32 - computation of the minimal box bounding the shape, given that this box have
33 to be aligned with the local coordinates;
34 - algorithms for dividing the shape along a given axis and producing resulting
35 divisions volumes.
36
37 The modeler currently provides a set of 16 basic shapes, which we will call
38primitives. It also provides a special class allowing the creation of shapes
39made as a result of boolean operations between primitives. These are called
40composite shapes and the composition operation can be recursive (composition
41of composites). This allows the creation of a quite large number of different
42shape topologies and combinations.
43
44 Named shapes register themselves to the manager class at creation time. The
45manager is responsible for their final deletion. Shapes can be created using their
46default constructor if their retrieval by name is not needed, but in this case
47they are owned by the user. A shape may be referenced by several volumes,
48therefore its deletion is not possible once volumes were defined based on it.
49
50### Creating shapes
51
52 Shape objects embed only the minimum set of parameters that are fully
53describing a valid physical shape. For instance, a tube is represented by
54its half length, the minimum radius and the maximum radius. Shapes are used
55together with media in order to create volumes, which in their turn
56are the main components of the geometrical tree. A specific shape can be created
57stand-alone :
58
59~~~ {.cpp}
60 TGeoBBox *box = new TGeoBBox("s_box", halfX, halfY, halfZ); // named
61 TGeoTube *tub = new TGeoTube(rmin, rmax, halfZ); // no name
62 ... (see each specific shape constructors)
63~~~
64
65 Sometimes it is much easier to create a volume having a given shape in one
66step, since shapes are not directly linked in the geometrical tree but volumes
67are :
68
69~~~ {.cpp}
70 TGeoVolume *vol_box = gGeoManager->MakeBox("BOX_VOL", "mat1", halfX, halfY, halfZ);
71 TGeoVolume *vol_tub = gGeoManager->MakeTube("TUB_VOL", "mat2", rmin, rmax, halfZ);
72 ... (see MakeXXX() utilities in TGeoManager class)
73~~~
74
75### Shape queries
76
77Note that global queries related to a geometry are handled by the manager class.
78However, shape-related queries might be sometimes useful.
79
80#### `Bool_t TGeoShape::Contains(const Double_t *point[3])`
81
82this method returns true if POINT is actually inside the shape. The point
83has to be defined in the local shape reference. For instance, for a box having
84DX, DY and DZ half-lengths a point will be considered inside if :
85
86~~~ {.cpp}
87 | -DX <= point[0] <= DX
88 | -DY <= point[1] <= DY
89 | -DZ <= point[2] <= DZ
90~~~
91
92#### `Double_t TGeoShape::DistFromInside(Double_t *point[3], Double_t *dir[3], Int_t iact, Double_t step, Double_t
93*safe)`
94
95computes the distance to exiting a shape from a given point INSIDE, along
96a given direction. The direction is given by its director cosines with respect
97to the local shape coordinate system. This method provides additional
98information according the value of IACT input parameter :
99
100 - IACT = 0 => compute only safe distance and fill it at the location
101 given by SAFE
102 - IACT = 1 => a proposed STEP is supplied. The safe distance is computed
103 first. If this is bigger than STEP than the proposed step
104 is approved and returned by the method since it does not
105 cross the shape boundaries. Otherwise, the distance to
106 exiting the shape is computed and returned.
107 - IACT = 2 => compute both safe distance and distance to exiting, ignoring
108 the proposed step.
109 - IACT > 2 => compute only the distance to exiting, ignoring anything else.
110
111#### `Double_t TGeoShape::DistFromOutside(Double_t *point[3], Double_t *dir[3], Int_t iact, Double_t step, Double_t
112*safe)`
113
114computes the distance to entering a shape from a given point OUTSIDE. Acts
115in the same way as B).
116
117#### `Double_t Safety(const Double_t *point[3], Bool_t inside)`
118
119compute maximum shift of a point in any direction that does not change its
120INSIDE/OUTSIDE state (does not cross shape boundaries). The state of the point
121have to be properly supplied.
122
123#### `Double_t *Normal(Double_t *point[3], Double_t *dir[3], Bool_t inside)`
124
125returns director cosines of normal to the crossed shape surface from a
126given point towards a direction. One has to specify if the point is inside
127or outside shape. According to this, the normal will be outwards or inwards
128shape respectively. Normal components are statically stored by shape class,
129so it has to be copied after retrieval in a different array.
130
131### Dividing shapes
132
133Shapes can generally be divided along a given axis. Supported axis are
134X, Y, Z, Rxy, Phi, Rxyz. A given shape cannot be divided however on any axis.
135The general rule is that that divisions are possible on whatever axis that
136produces still known shapes as slices. The division of shapes should not be
137performed by TGeoShape::Divide() calls, but rather by TGeoVolume::Divide().
138The algorithm for dividing a specific shape is known by the shape object, but
139is always invoked in a generic way from the volume level. Details on how to
140do that can be found in TGeoVolume class. One can see how all division options
141are interpreted and which is their result inside specific shape classes.
142
143\image html geom_t_shape.png width=600px
144*/
145
146#include "TObjArray.h"
147#include "TEnv.h"
148#include "TError.h"
149
150#include "TGeoMatrix.h"
151#include "TGeoManager.h"
152#include "TGeoVolume.h"
153#include "TGeoShape.h"
154#include "TVirtualGeoPainter.h"
155#include "TBuffer3D.h"
156#include "TBuffer3DTypes.h"
157#include "TMath.h"
158
160
162Double_t TGeoShape::fgEpsMch = 2.220446049250313e-16;
163
164////////////////////////////////////////////////////////////////////////////////
165/// Default constructor
166
168{
169 fShapeBits = 0;
170 fShapeId = 0;
171 if (!gGeoManager) {
172 gGeoManager = new TGeoManager("Geometry", "default geometry");
173 // gROOT->AddGeoManager(gGeoManager);
174 }
175 // fShapeId = gGeoManager->GetListOfShapes()->GetSize();
176 // gGeoManager->AddShape(this);
177}
178
179////////////////////////////////////////////////////////////////////////////////
180/// Default constructor
181
183{
184 fShapeBits = 0;
185 fShapeId = 0;
186 if (!gGeoManager) {
187 gGeoManager = new TGeoManager("Geometry", "default geometry");
188 // gROOT->AddGeoManager(gGeoManager);
189 }
191 gGeoManager->AddShape(this);
192}
193
194////////////////////////////////////////////////////////////////////////////////
195/// Destructor
196
198{
201}
202
203////////////////////////////////////////////////////////////////////////////////
204/// Test for shape navigation methods. Summary for test numbers:
205///
206/// - 1: DistFromInside/Outside. Sample points inside the shape. Generate
207/// directions randomly in cos(theta). Compute DistFromInside and move the
208/// point with bigger distance. Compute DistFromOutside back from new point.
209/// Plot d-(d1+d2)
210
212{
213 if (!gGeoManager) {
214 Error("CheckShape", "No geometry manager");
215 return;
216 }
217 TGeoShape *shape = (TGeoShape *)this;
218 gGeoManager->CheckShape(shape, testNo, nsamples, option);
219}
220
221////////////////////////////////////////////////////////////////////////////////
222/// Compute machine round-off double precision error as the smallest number that
223/// if added to 1.0 is different than 1.0.
224
226{
227 Double_t temp1 = 1.0;
228 Double_t temp2 = 1.0 + temp1;
229 Double_t mchEps = 0.;
230 while (temp2 > 1.0) {
231 mchEps = temp1;
232 temp1 /= 2;
233 temp2 = 1.0 + temp1;
234 }
235 fgEpsMch = mchEps;
236 return fgEpsMch;
237}
238
239////////////////////////////////////////////////////////////////////////////////
240/// static function returning the machine round-off error
241
243{
244 return fgEpsMch;
245}
246
247////////////////////////////////////////////////////////////////////////////////
248/// Get the shape name.
249
250const char *TGeoShape::GetName() const
251{
252 if (!fName[0]) {
253 return ((TObject *)this)->ClassName();
254 }
255 return TNamed::GetName();
256}
257
258////////////////////////////////////////////////////////////////////////////////
259/// Returns distance to shape primitive mesh.
260
262{
264 if (!painter)
265 return 9999;
266 return painter->ShapeDistancetoPrimitive(this, numpoints, px, py);
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,s2
271
272Bool_t
274{
275 Double_t saf1 = TGeoShape::Big();
276 Double_t saf2 = TGeoShape::Big();
277 if (point[0] * c1 + point[1] * s1 >= 0)
278 saf1 = TMath::Abs(-point[0] * s1 + point[1] * c1);
279 if (point[0] * c2 + point[1] * s2 >= 0)
280 saf2 = TMath::Abs(point[0] * s2 - point[1] * c2);
281 Double_t saf = TMath::Min(saf1, saf2);
282 if (saf < epsil)
283 return kTRUE;
284 return kFALSE;
285}
286
287////////////////////////////////////////////////////////////////////////////////
288/// Static method to check if a point is in the phi range (phi1, phi2) [degrees]
289
291{
292 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
293 while (phi < phi1)
294 phi += 360.;
295 Double_t ddp = phi - phi1;
296 if (ddp > phi2 - phi1)
297 return kFALSE;
298 return kTRUE;
299}
300
301////////////////////////////////////////////////////////////////////////////////
302/// Compute distance from POINT to semiplane defined by PHI angle along DIR. Computes
303/// also radius at crossing point. This might be negative in case the crossing is
304/// on the other side of the semiplane.
305
307 Double_t &snext, Double_t &rxy)
308{
309 snext = rxy = TGeoShape::Big();
310 Double_t nx = -sphi;
311 Double_t ny = cphi;
312 Double_t rxy0 = point[0] * cphi + point[1] * sphi;
313 Double_t rdotn = point[0] * nx + point[1] * ny;
314 if (TMath::Abs(rdotn) < TGeoShape::Tolerance()) {
315 snext = 0.0;
316 rxy = rxy0;
317 return kTRUE;
318 }
319 if (rdotn < 0) {
320 rdotn = -rdotn;
321 } else {
322 nx = -nx;
323 ny = -ny;
324 }
325 Double_t ddotn = dir[0] * nx + dir[1] * ny;
326 if (ddotn <= 0)
327 return kFALSE;
328 snext = rdotn / ddotn;
329 rxy = rxy0 + snext * (dir[0] * cphi + dir[1] * sphi);
330 if (rxy < 0)
331 return kFALSE;
332 return kTRUE;
333}
334
335////////////////////////////////////////////////////////////////////////////////
336/// Check if two numbers differ with less than a tolerance.
337
339{
340 if (TMath::Abs(a - b) < 1.E-10)
341 return kTRUE;
342 return kFALSE;
343}
344
345////////////////////////////////////////////////////////////////////////////////
346/// Check if segments (A,B) and (C,D) are crossing,
347/// where: A(x1,y1), B(x2,y2), C(x3,y3), D(x4,y4)
348
350 Double_t x4, Double_t y4)
351{
353 Bool_t stand1 = kFALSE;
354 Double_t dx1 = x2 - x1;
355 Bool_t stand2 = kFALSE;
356 Double_t dx2 = x4 - x3;
357 Double_t xm = 0.;
358 Double_t ym = 0.;
359 Double_t a1 = 0.;
360 Double_t b1 = 0.;
361 Double_t a2 = 0.;
362 Double_t b2 = 0.;
363 if (TMath::Abs(dx1) < eps)
364 stand1 = kTRUE;
365 if (TMath::Abs(dx2) < eps)
366 stand2 = kTRUE;
367 if (!stand1) {
368 a1 = (x2 * y1 - x1 * y2) / dx1;
369 b1 = (y2 - y1) / dx1;
370 }
371 if (!stand2) {
372 a2 = (x4 * y3 - x3 * y4) / dx2;
373 b2 = (y4 - y3) / dx2;
374 }
375 if (stand1 && stand2) {
376 // Segments parallel and vertical
377 if (TMath::Abs(x1 - x3) < eps) {
378 // Check if segments are overlapping
379 if ((y3 - y1) * (y3 - y2) < -eps || (y4 - y1) * (y4 - y2) < -eps || (y1 - y3) * (y1 - y4) < -eps ||
380 (y2 - y3) * (y2 - y4) < -eps)
381 return kTRUE;
382 return kFALSE;
383 }
384 // Different x values
385 return kFALSE;
386 }
387
388 if (stand1) {
389 // First segment vertical
390 xm = x1;
391 ym = a2 + b2 * xm;
392 } else {
393 if (stand2) {
394 // Second segment vertical
395 xm = x3;
396 ym = a1 + b1 * xm;
397 } else {
398 // Normal crossing
399 if (TMath::Abs(b1 - b2) < eps) {
400 // Parallel segments, are they aligned
401 if (TMath::Abs(y3 - (a1 + b1 * x3)) > eps)
402 return kFALSE;
403 // Aligned segments, are they overlapping
404 if ((x3 - x1) * (x3 - x2) < -eps || (x4 - x1) * (x4 - x2) < -eps || (x1 - x3) * (x1 - x4) < -eps ||
405 (x2 - x3) * (x2 - x4) < -eps)
406 return kTRUE;
407 return kFALSE;
408 }
409 xm = (a1 - a2) / (b2 - b1);
410 ym = (a1 * b2 - a2 * b1) / (b2 - b1);
411 }
412 }
413 // Check if crossing point is both between A,B and C,D
414 Double_t check = (xm - x1) * (xm - x2) + (ym - y1) * (ym - y2);
415 if (check > -eps)
416 return kFALSE;
417 check = (xm - x3) * (xm - x4) + (ym - y3) * (ym - y4);
418 if (check > -eps)
419 return kFALSE;
420 return kTRUE;
421}
422
423////////////////////////////////////////////////////////////////////////////////
424/// compute distance from point (inside phi) to both phi planes. Return minimum.
425
427 Double_t c2, Double_t sm, Double_t cm, Bool_t in)
428{
429 Double_t sfi1 = TGeoShape::Big();
430 Double_t sfi2 = TGeoShape::Big();
431 Double_t s = 0;
432 Double_t un = dir[0] * s1 - dir[1] * c1;
433 if (!in)
434 un = -un;
435 if (un > 0) {
436 s = -point[0] * s1 + point[1] * c1;
437 if (!in)
438 s = -s;
439 if (s >= 0) {
440 s /= un;
441 if (((point[0] + s * dir[0]) * sm - (point[1] + s * dir[1]) * cm) >= 0)
442 sfi1 = s;
443 }
444 }
445 un = -dir[0] * s2 + dir[1] * c2;
446 if (!in)
447 un = -un;
448 if (un > 0) {
449 s = point[0] * s2 - point[1] * c2;
450 if (!in)
451 s = -s;
452 if (s >= 0) {
453 s /= un;
454 if ((-(point[0] + s * dir[0]) * sm + (point[1] + s * dir[1]) * cm) >= 0)
455 sfi2 = s;
456 }
457 }
458 return TMath::Min(sfi1, sfi2);
459}
460
461////////////////////////////////////////////////////////////////////////////////
462/// Static method to compute normal to phi planes.
463
464void TGeoShape::NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1,
465 Double_t c2, Double_t s2)
466{
467 Double_t saf1 = TGeoShape::Big();
468 Double_t saf2 = TGeoShape::Big();
469 if (point[0] * c1 + point[1] * s1 >= 0)
470 saf1 = TMath::Abs(-point[0] * s1 + point[1] * c1);
471 if (point[0] * c2 + point[1] * s2 >= 0)
472 saf2 = TMath::Abs(point[0] * s2 - point[1] * c2);
473 Double_t c, s;
474 if (saf1 < saf2) {
475 c = c1;
476 s = s1;
477 } else {
478 c = c2;
479 s = s2;
480 }
481 norm[2] = 0;
482 norm[0] = -s;
483 norm[1] = c;
484 if (dir[0] * norm[0] + dir[1] * norm[1] < 0) {
485 norm[0] = s;
486 norm[1] = -c;
487 }
488}
489
490////////////////////////////////////////////////////////////////////////////////
491/// Static method to compute safety w.r.t a phi corner defined by cosines/sines
492/// of the angles phi1, phi2.
493
495{
496 Bool_t inphi = TGeoShape::IsInPhiRange(point, phi1, phi2);
497 if (inphi && !in)
498 return -TGeoShape::Big();
499 phi1 *= TMath::DegToRad();
500 phi2 *= TMath::DegToRad();
501 Double_t c1 = TMath::Cos(phi1);
502 Double_t s1 = TMath::Sin(phi1);
503 Double_t c2 = TMath::Cos(phi2);
504 Double_t s2 = TMath::Sin(phi2);
505 Double_t rsq = point[0] * point[0] + point[1] * point[1];
506 Double_t rproj = point[0] * c1 + point[1] * s1;
507 Double_t safsq = rsq - rproj * rproj;
508 if (safsq < 0)
509 return 0.;
510 Double_t saf1 = (rproj < 0) ? TGeoShape::Big() : TMath::Sqrt(safsq);
511 rproj = point[0] * c2 + point[1] * s2;
512 safsq = rsq - rproj * rproj;
513 if (safsq < 0)
514 return 0.;
515 Double_t saf2 = (rproj < 0) ? TGeoShape::Big() : TMath::Sqrt(safsq);
516 Double_t safe = TMath::Min(saf1, saf2); // >0
517 if (safe > 1E10) {
518 if (in)
519 return TGeoShape::Big();
520 return -TGeoShape::Big();
521 }
522 return safe;
523}
524
525////////////////////////////////////////////////////////////////////////////////
526/// Compute distance from point of coordinates (r,z) to segment (r1,z1):(r2,z2)
527
529{
530 Double_t crossp = (z2 - z1) * (r - r1) - (z - z1) * (r2 - r1);
531 crossp *= (outer) ? 1. : -1.;
532 // Positive crossp means point on the requested side of the (1,2) segment
533 if (crossp < -TGeoShape::Tolerance()) {
534 // if (((z-z1)*(z2-z)) > -1.E-10) return 0;
535 if (outer)
536 return TGeoShape::Big();
537 else
538 return 0.;
539 }
540 // Compute (1,P) dot (1,2)
541 Double_t c1 = (z - z1) * (z2 - z1) + (r - r1) * (r2 - r1);
542 // Negative c1 means point (1) is closest
543 if (c1 < 1.E-10)
544 return TMath::Sqrt((r - r1) * (r - r1) + (z - z1) * (z - z1));
545 // Compute (2,P) dot (1,2)
546 Double_t c2 = (z - z2) * (z2 - z1) + (r - r2) * (r2 - r1);
547 // Positive c2 means point (2) is closest
548 if (c2 > -1.E-10)
549 return TMath::Sqrt((r - r2) * (r - r2) + (z - z2) * (z - z2));
550 // The closest point is between (1) and (2)
551 c2 = (z2 - z1) * (z2 - z1) + (r2 - r1) * (r2 - r1);
552 // projected length factor with respect to (1,2) length
553 Double_t alpha = c1 / c2;
554 Double_t rp = r1 + alpha * (r2 - r1);
555 Double_t zp = z1 + alpha * (z2 - z1);
556 return TMath::Sqrt((r - rp) * (r - rp) + (z - zp) * (z - zp));
557}
558
559////////////////////////////////////////////////////////////////////////////////
560/// Equivalent of TObject::SetBit.
561
563{
564 if (set) {
565 SetShapeBit(f);
566 } else {
568 }
569}
570
571////////////////////////////////////////////////////////////////////////////////
572/// Returns current transformation matrix that applies to shape.
573
575{
576 return fgTransform;
577}
578
579////////////////////////////////////////////////////////////////////////////////
580/// Set current transformation matrix that applies to shape.
581
583{
584 fgTransform = matrix;
585}
586
587////////////////////////////////////////////////////////////////////////////////
588/// Tranform a set of points (LocalToMaster)
589
591{
592 UInt_t i, j;
593 Double_t dmaster[3];
594 if (fgTransform) {
595 for (j = 0; j < NbPnts; j++) {
596 i = 3 * j;
597 fgTransform->LocalToMaster(&points[i], dmaster);
598 points[i] = dmaster[0];
599 points[i + 1] = dmaster[1];
600 points[i + 2] = dmaster[2];
601 }
602 return;
603 }
604 if (!gGeoManager)
605 return;
606 Bool_t bomb = (gGeoManager->GetBombMode() == 0) ? kFALSE : kTRUE;
607
608 for (j = 0; j < NbPnts; j++) {
609 i = 3 * j;
612 if (bomb)
613 glmat->LocalToMasterBomb(&points[i], dmaster);
614 else
615 glmat->LocalToMaster(&points[i], dmaster);
616 } else {
617 if (bomb)
618 gGeoManager->LocalToMasterBomb(&points[i], dmaster);
619 else
620 gGeoManager->LocalToMaster(&points[i], dmaster);
621 }
622 points[i] = dmaster[0];
623 points[i + 1] = dmaster[1];
624 points[i + 2] = dmaster[2];
625 }
626}
627
628////////////////////////////////////////////////////////////////////////////////
629/// Fill the supplied buffer, with sections in desired frame
630/// See TBuffer3D.h for explanation of sections, frame etc.
631
632void TGeoShape::FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
633{
634 // Catch this common potential error here
635 // We have to set kRawSize (unless already done) to allocate buffer space
636 // before kRaw can be filled
637 if (reqSections & TBuffer3D::kRaw) {
638 if (!(reqSections & TBuffer3D::kRawSizes) && !buffer.SectionsValid(TBuffer3D::kRawSizes)) {
640 }
641 }
642
643 if (reqSections & TBuffer3D::kCore) {
644 // If writing core section all others will be invalid
645 buffer.ClearSectionsValid();
646
647 // Check/grab some objects we need
648 if (!gGeoManager) {
650 return;
651 }
652 const TGeoVolume *paintVolume = gGeoManager->GetPaintVolume();
653 if (!paintVolume)
654 paintVolume = gGeoManager->GetTopVolume();
655 if (!paintVolume) {
656 buffer.fID = const_cast<TGeoShape *>(this);
657 buffer.fColor = 0;
658 buffer.fTransparency = 0;
659 // R__ASSERT(kFALSE);
660 // return;
661 } else {
662 buffer.fID = const_cast<TGeoVolume *>(paintVolume);
663 buffer.fColor = paintVolume->GetLineColor();
664
665 buffer.fTransparency = paintVolume->GetTransparency();
666 Double_t visdensity = gGeoManager->GetVisDensity();
667 if (visdensity > 0 && paintVolume->GetMedium()) {
668 if (paintVolume->GetMaterial()->GetDensity() < visdensity) {
669 buffer.fTransparency = 90;
670 }
671 }
672 }
673
674 buffer.fLocalFrame = localFrame;
675 Bool_t r1, r2 = kFALSE;
677 if (paintVolume && paintVolume->GetShape()) {
678 if (paintVolume->GetShape()->IsReflected()) {
679 // Temporary trick to deal with reflected shapes.
680 // Still lighting gets wrong...
681 if (buffer.Type() < TBuffer3DTypes::kTube)
682 r2 = kTRUE;
683 }
684 }
685 buffer.fReflection = ((r1 & (!r2)) | (r2 & !(r1)));
686
687 // Set up local -> master translation matrix
688 if (localFrame) {
689 TGeoMatrix *localMasterMat = nullptr;
691 localMasterMat = TGeoShape::GetTransform();
692 } else {
693 localMasterMat = gGeoManager->GetCurrentMatrix();
694
695 // For overlap drawing the correct matrix needs to obtained in
696 // from GetGLMatrix() - this should not be applied in the case
697 // of composite shapes
699 localMasterMat = gGeoManager->GetGLMatrix();
700 }
701 }
702 if (!localMasterMat) {
704 return;
705 }
706 localMasterMat->GetHomogenousMatrix(buffer.fLocalMaster);
707 } else {
708 buffer.SetLocalMasterIdentity();
709 }
710
712 }
713}
714
715////////////////////////////////////////////////////////////////////////////////
716/// Get the basic color (0-7).
717
719{
720 Int_t basicColor = 0; // TODO: Check on sensible fallback
721 if (gGeoManager) {
722 const TGeoVolume *volume = gGeoManager->GetPaintVolume();
723 if (volume) {
724 basicColor = ((volume->GetLineColor() % 8) - 1) * 4;
725 if (basicColor < 0)
726 basicColor = 0;
727 }
728 }
729 return basicColor;
730}
731
732////////////////////////////////////////////////////////////////////////////////
733/// Stub implementation to avoid forcing implementation at this stage
734
735const TBuffer3D &TGeoShape::GetBuffer3D(Int_t /*reqSections*/, Bool_t /*localFrame*/) const
736{
737 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
738 Warning("GetBuffer3D", "this must be implemented for shapes in a TGeoPainter hierarchy. This will be come a pure "
739 "virtual fn eventually.");
740 return buffer;
741}
742
743////////////////////////////////////////////////////////////////////////////////
744/// Provide a pointer name containing uid.
745
746const char *TGeoShape::GetPointerName() const
747{
748 static TString name;
749 Int_t uid = GetUniqueID();
750 if (uid)
751 name = TString::Format("p%s_%d", GetName(), uid);
752 else
753 name = TString::Format("p%s", GetName());
754 return name.Data();
755}
756
757////////////////////////////////////////////////////////////////////////////////
758/// Execute mouse actions on this shape.
759
761{
762 if (!gGeoManager)
763 return;
765 painter->ExecuteShapeEvent(this, event, px, py);
766}
767
768////////////////////////////////////////////////////////////////////////////////
769/// Draw this shape.
770
772{
774 if (option && option[0]) {
775 painter->DrawShape(this, option);
776 } else {
777 painter->DrawShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption", ""));
778 }
779}
780
781////////////////////////////////////////////////////////////////////////////////
782/// Paint this shape.
783
785{
787 if (option && option[0]) {
788 painter->PaintShape(this, option);
789 } else {
790 painter->PaintShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption", ""));
791 }
792}
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define s1(x)
Definition RSha256.hxx:91
bool Bool_t
Definition RtypesCore.h:63
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
#define R__ASSERT(e)
Definition TError.h:118
Option_t Option_t option
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 x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:33
Generic 3D primitive description class.
Definition TBuffer3D.h:18
void SetLocalMasterIdentity()
Set kRaw tessellation section of buffer with supplied sizes.
Int_t Type() const
Definition TBuffer3D.h:85
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
Double_t fLocalMaster[16]
Definition TBuffer3D.h:93
void ClearSectionsValid()
Clear any sections marked valid.
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Int_t fColor
Definition TBuffer3D.h:88
Short_t fTransparency
Definition TBuffer3D.h:89
Bool_t fReflection
Definition TBuffer3D.h:91
TObject * fID
Definition TBuffer3D.h:87
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:458
The manager class for any TGeo geometry.
Definition TGeoManager.h:44
TGeoHMatrix * GetGLMatrix() const
Bool_t IsMatrixTransform() const
void LocalToMaster(const Double_t *local, Double_t *master) const
TGeoVolume * GetPaintVolume() const
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
Double_t GetVisDensity() const
void LocalToMasterBomb(const Double_t *local, Double_t *master) const
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
Test for shape navigation methods.
TGeoHMatrix * GetCurrentMatrix() const
Bool_t IsMatrixReflection() const
TVirtualGeoPainter * GetPainter() const
Int_t GetBombMode() const
TObjArray * GetListOfShapes() const
TGeoVolume * GetTopVolume() const
Bool_t IsCleaning() const
Int_t AddShape(const TGeoShape *shape)
Add a shape to the list. Returns index of the shape in list.
virtual Double_t GetDensity() const
Geometrical transformation package.
Definition TGeoMatrix.h:38
virtual void LocalToMasterBomb(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
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
void GetHomogenousMatrix(Double_t *hmat) const
The homogenous matrix associated with the transformation is used for piling up's and visualization.
Base abstract class for all shapes.
Definition TGeoShape.h:25
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Stub implementation to avoid forcing implementation at this stage.
UInt_t fShapeBits
Definition TGeoShape.h:72
static Double_t Big()
Definition TGeoShape.h:87
Int_t GetBasicColor() const
Get the basic color (0-7).
static Bool_t IsSegCrossing(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t x3, Double_t y3, Double_t x4, Double_t y4)
Check if segments (A,B) and (C,D) are crossing, where: A(x1,y1), B(x2,y2), C(x3,y3),...
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
void ResetShapeBit(UInt_t f)
Definition TGeoShape.h:166
static Double_t DistToPhiMin(const Double_t *point, const Double_t *dir, Double_t s1, Double_t c1, Double_t s2, Double_t c2, Double_t sm, Double_t cm, Bool_t in=kTRUE)
compute distance from point (inside phi) to both phi planes. Return minimum.
TGeoShape()
Default constructor.
static Double_t SafetyPhi(const Double_t *point, Bool_t in, Double_t phi1, Double_t phi2)
Static method to compute safety w.r.t a phi corner defined by cosines/sines of the angles phi1,...
Int_t fShapeId
Definition TGeoShape.h:71
static void SetTransform(TGeoMatrix *matrix)
Set current transformation matrix that applies to shape.
virtual Bool_t IsComposite() const
Definition TGeoShape.h:130
void Draw(Option_t *option="") override
Draw this shape.
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
const char * GetPointerName() const
Provide a pointer name containing uid.
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
Fill the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections,...
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Test for shape navigation methods.
static Double_t EpsMch()
static function returning the machine round-off error
static Bool_t IsInPhiRange(const Double_t *point, Double_t phi1, Double_t phi2)
Static method to check if a point is in the phi range (phi1, phi2) [degrees].
static Double_t ComputeEpsMch()
Compute machine round-off double precision error as the smallest number that if added to 1....
void Paint(Option_t *option="") override
Paint this shape.
~TGeoShape() override
Destructor.
static TGeoMatrix * fgTransform
Definition TGeoShape.h:27
virtual Bool_t IsReflected() const
Definition TGeoShape.h:140
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
static Double_t fgEpsMch
Definition TGeoShape.h:28
static void NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Static method to compute normal to phi planes.
static Double_t SafetySeg(Double_t r, Double_t z, Double_t r1, Double_t z1, Double_t r2, Double_t z2, Bool_t outer)
Compute distance from point of coordinates (r,z) to segment (r1,z1):(r2,z2)
static Bool_t IsCrossingSemiplane(const Double_t *point, const Double_t *dir, Double_t cphi, Double_t sphi, Double_t &snext, Double_t &rxy)
Compute distance from POINT to semiplane defined by PHI angle along DIR.
const char * GetName() const override
Get the shape name.
static TGeoMatrix * GetTransform()
Returns current transformation matrix that applies to shape.
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute mouse actions on this shape.
static Double_t Tolerance()
Definition TGeoShape.h:90
static Bool_t IsCloseToPhi(Double_t epsil, const Double_t *point, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,...
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:175
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:174
TGeoShape * GetShape() const
Definition TGeoVolume.h:190
Char_t GetTransparency() const
Definition TGeoVolume.h:369
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
TString fName
Definition TNamed.h:32
TObject * Remove(TObject *obj) override
Remove object from array.
Mother of all ROOT objects.
Definition TObject.h:41
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition TObject.cxx:457
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:962
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:976
Basic string class.
Definition TString.h:139
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:2378
Abstract class for geometry painters.
virtual void DrawShape(TGeoShape *shape, Option_t *option="")=0
virtual void PaintShape(TGeoShape *shape, Option_t *option="")=0
virtual void ExecuteShapeEvent(TGeoShape *shape, Int_t event, Int_t px, Int_t py)=0
virtual Int_t ShapeDistancetoPrimitive(const TGeoShape *shape, Int_t numpoints, Int_t px, Int_t py) const =0
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:646
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:594
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:588
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:72
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123