Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoBBox.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$// Author: Andrei Gheata 24/10/01
2
3// Contains() and DistFromOutside/Out() implemented by Mihaela Gheata
4
5/*************************************************************************
6 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13/** \class TGeoBBox
14\ingroup Shapes_classes
15\brief Box class.
16
17 - [Building boxes](\ref GEOB00)
18 - [Creation of boxes](\ref GEOB01)
19 - [Divisions of boxes](\ref GEOB02)
20
21All shape primitives inherit from this, their
22 constructor filling automatically the parameters of the box that bounds
23 the given shape. Defined by 6 parameters :
24
25```
26TGeoBBox(Double_t dx,Double_t dy,Double_t dz,Double_t *origin=0);
27```
28
29 - `fDX`, `fDY`, `fDZ` : half lengths on X, Y and Z axis
30 - `fOrigin[3]` : position of box origin
31
32\anchor GEOB00
33### Building boxes
34
35Normally a box has to be built only with 3 parameters: `DX,DY,DZ`
36representing the half-lengths on X, Y and Z-axes. In this case, the
37origin of the box will match the one of its reference frame and the box
38will range from: `-DX` to `DX` on X-axis, from `-DY` to `DY` on Y and
39from `-DZ` to `DZ` on Z. On the other hand, any other shape needs to
40compute and store the parameters of their minimal bounding box. The
41bounding boxes are essential to optimize navigation algorithms.
42Therefore all other primitives derive from **`TGeoBBox`**. Since the
43minimal bounding box is not necessary centered in the origin, any box
44allows an origin translation `(Ox`,`Oy`,`Oz)`. All primitive
45constructors automatically compute the bounding box parameters. Users
46should be aware that building a translated box that will represent a
47primitive shape by itself would affect any further positioning of other
48shapes inside. Therefore it is highly recommendable to build
49non-translated boxes as primitives and translate/rotate their
50corresponding volumes only during positioning stage.
51
52\anchor GEOB01
53#### Creation of boxes
54
55```
56 TGeoBBox *box = new TGeoBBox("BOX", 20, 30, 40);
57```
58
59Begin_Macro
60{
61 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
62 new TGeoManager("box", "poza1");
63 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
64 TGeoMedium *med = new TGeoMedium("MED",1,mat);
65 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
66 gGeoManager->SetTopVolume(top);
67 TGeoVolume *vol = gGeoManager->MakeBox("BOX",med, 20,30,40);
68 vol->SetLineWidth(2);
69 top->AddNode(vol,1);
70 gGeoManager->CloseGeometry();
71 gGeoManager->SetNsegments(80);
72 top->Draw();
73 TView *view = gPad->GetView();
74 if (view) view->ShowAxis();
75}
76End_Macro
77
78A volume having a box shape can be built in one step:
79
80```
81 TGeoVolume *vbox = gGeoManager->MakeBox("vbox", ptrMed, 20,30,40);
82```
83
84\anchor GEOB02
85#### Divisions of boxes
86
87 Volumes having box shape can be divided with equal-length slices on
88X, Y or Z axis. The following options are supported:
89
90 - Dividing the full range of one axis in N slices
91```
92 TGeoVolume *divx = vbox->Divide("SLICEX", 1, N);
93```
94here `1` stands for the division axis (1-X, 2-Y, 3-Z)
95
96Begin_Macro
97{
98 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
99 new TGeoManager("box", "poza1");
100 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
101 TGeoMedium *med = new TGeoMedium("MED",1,mat);
102 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
103 gGeoManager->SetTopVolume(top);
104 TGeoVolume *vol = gGeoManager->MakeBox("BOX",med, 20,30,40);
105 vol->SetLineWidth(2);
106 top->AddNode(vol,1);
107 TGeoVolume *divx = vol->Divide("SLICE",1,8,0,0);
108 gGeoManager->CloseGeometry();
109 gGeoManager->SetNsegments(80);
110 top->Draw();
111 TView *view = gPad->GetView();
112 if (view) view->ShowAxis();
113}
114End_Macro
115
116 - Dividing in a limited range - general case.
117```
118 TGeoVolume *divy = vbox->Divide("SLICEY",2,N,start,step);
119```
120 - start = starting offset within (-fDY, fDY)
121 - step = slicing step
122
123Begin_Macro
124{
125 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
126 new TGeoManager("box", "poza1");
127 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
128 TGeoMedium *med = new TGeoMedium("MED",1,mat);
129 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
130 gGeoManager->SetTopVolume(top);
131 TGeoVolume *vol = gGeoManager->MakeBox("BOX",med, 20,30,40);
132 vol->SetLineWidth(2);
133 top->AddNode(vol,1);
134 TGeoVolume *divx = vol->Divide("SLICE",2,8,2,3);
135 gGeoManager->CloseGeometry();
136 gGeoManager->SetNsegments(80);
137 top->Draw();
138 TView *view = gPad->GetView();
139 if (view) view->ShowAxis();
140}
141End_Macro
142
143Both cases are supported by all shapes.
144See also class TGeoShape for utility methods provided by any particular shape.
145*/
146
147#include <iostream>
148
149#include "TGeoManager.h"
150#include "TGeoMatrix.h"
151#include "TGeoVolume.h"
152#include "TVirtualGeoPainter.h"
153#include "TGeoBBox.h"
154#include "TBuffer3D.h"
155#include "TBuffer3DTypes.h"
156#include "TMath.h"
157#include "TRandom.h"
158
159////////////////////////////////////////////////////////////////////////////////
160/// Default constructor
161
163{
165 fDX = fDY = fDZ = 0;
166 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
167}
168
169////////////////////////////////////////////////////////////////////////////////
170/// Constructor where half-lengths are provided.
171
178
179////////////////////////////////////////////////////////////////////////////////
180/// Constructor with shape name.
181
188
189////////////////////////////////////////////////////////////////////////////////
190/// Constructor based on the array of parameters.
191/// - param[0] - half-length in x
192/// - param[1] - half-length in y
193/// - param[2] - half-length in z
194
196{
198 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
199 SetDimensions(param);
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Destructor
204
206
207////////////////////////////////////////////////////////////////////////////////
208/// Computes capacity of the shape in [length^3].
209
211{
212 return (8. * fDX * fDY * fDZ);
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Computes normal to closest surface from POINT.
217
218void TGeoBBox::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const
219{
220 memset(norm, 0, 3 * sizeof(Double_t));
221 Double_t saf[3];
222 Int_t i;
223 saf[0] = TMath::Abs(TMath::Abs(point[0] - fOrigin[0]) - fDX);
224 saf[1] = TMath::Abs(TMath::Abs(point[1] - fOrigin[1]) - fDY);
225 saf[2] = TMath::Abs(TMath::Abs(point[2] - fOrigin[2]) - fDZ);
226 i = TMath::LocMin(3, saf);
227 norm[i] = (dir[i] > 0) ? 1 : (-1);
228}
229
230////////////////////////////////////////////////////////////////////////////////
231/// Decides fast if the bounding box could be crossed by a vector.
232
233Bool_t TGeoBBox::CouldBeCrossed(const Double_t *point, const Double_t *dir) const
234{
235 Double_t mind = fDX;
236 if (fDY < mind)
237 mind = fDY;
238 if (fDZ < mind)
239 mind = fDZ;
240 Double_t dx = fOrigin[0] - point[0];
241 Double_t dy = fOrigin[1] - point[1];
242 Double_t dz = fOrigin[2] - point[2];
243 Double_t do2 = dx * dx + dy * dy + dz * dz;
244 if (do2 <= (mind * mind))
245 return kTRUE;
246 Double_t rmax2 = fDX * fDX + fDY * fDY + fDZ * fDZ;
247 if (do2 <= rmax2)
248 return kTRUE;
249 // inside bounding sphere
250 Double_t doct = dx * dir[0] + dy * dir[1] + dz * dir[2];
251 // leaving ray
252 if (doct <= 0)
253 return kFALSE;
254 Double_t dirnorm = dir[0] * dir[0] + dir[1] * dir[1] + dir[2] * dir[2];
255 if ((doct * doct) >= (do2 - rmax2) * dirnorm)
256 return kTRUE;
257 return kFALSE;
258}
259
260////////////////////////////////////////////////////////////////////////////////
261/// Compute closest distance from point px,py to each corner.
262
264{
265 const Int_t numPoints = 8;
266 return ShapeDistancetoPrimitive(numPoints, px, py);
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Divide this box shape belonging to volume "voldiv" into ndiv equal volumes
271/// called divname, from start position with the given step. Returns pointer
272/// to created division cell volume. In case a wrong division axis is supplied,
273/// returns pointer to volume to be divided.
274
277{
278 TGeoShape *shape; //--- shape to be created
279 TGeoVolume *vol; //--- division volume to be created
280 TGeoVolumeMulti *vmulti; //--- generic divided volume
281 TGeoPatternFinder *finder; //--- finder to be attached
282 TString opt = ""; //--- option to be attached
283 Double_t end = start + ndiv * step;
284 switch (iaxis) {
285 case 1: //--- divide on X
286 shape = new TGeoBBox(step / 2., fDY, fDZ);
287 finder = new TGeoPatternX(voldiv, ndiv, start, end);
288 opt = "X";
289 break;
290 case 2: //--- divide on Y
291 shape = new TGeoBBox(fDX, step / 2., fDZ);
292 finder = new TGeoPatternY(voldiv, ndiv, start, end);
293 opt = "Y";
294 break;
295 case 3: //--- divide on Z
296 shape = new TGeoBBox(fDX, fDY, step / 2.);
297 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
298 opt = "Z";
299 break;
300 default: Error("Divide", "Wrong axis type for division"); return nullptr;
301 }
302 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
304 vmulti->AddVolume(vol);
305 voldiv->SetFinder(finder);
306 finder->SetDivIndex(voldiv->GetNdaughters());
307 for (Int_t ic = 0; ic < ndiv; ic++) {
308 voldiv->AddNodeOffset(vol, ic, start + step / 2. + ic * step, opt.Data());
309 ((TGeoNodeOffset *)voldiv->GetNodes()->At(voldiv->GetNdaughters() - 1))->SetFinder(finder);
310 }
311 return vmulti;
312}
313
314////////////////////////////////////////////////////////////////////////////////
315/// Compute bounding box - nothing to do in this case.
316
318
319////////////////////////////////////////////////////////////////////////////////
320/// Test if point is inside this shape.
321
323{
324 if (TMath::Abs(point[2] - fOrigin[2]) > fDZ)
325 return kFALSE;
326 if (TMath::Abs(point[0] - fOrigin[0]) > fDX)
327 return kFALSE;
328 if (TMath::Abs(point[1] - fOrigin[1]) > fDY)
329 return kFALSE;
330 return kTRUE;
331}
332
333////////////////////////////////////////////////////////////////////////////////
334/// Static method to check if point[3] is located inside a box of having dx, dy, dz
335/// as half-lengths.
336
338{
339 if (TMath::Abs(point[2] - origin[2]) > dz)
340 return kFALSE;
341 if (TMath::Abs(point[0] - origin[0]) > dx)
342 return kFALSE;
343 if (TMath::Abs(point[1] - origin[1]) > dy)
344 return kFALSE;
345 return kTRUE;
346}
347
348////////////////////////////////////////////////////////////////////////////////
349/// Compute distance from inside point to surface of the box.
350/// Boundary safe algorithm.
351
354{
355 Double_t s, smin, saf[6];
356 Double_t newpt[3];
357 Int_t i;
358 for (i = 0; i < 3; i++)
359 newpt[i] = point[i] - fOrigin[i];
360 saf[0] = fDX + newpt[0];
361 saf[1] = fDX - newpt[0];
362 saf[2] = fDY + newpt[1];
363 saf[3] = fDY - newpt[1];
364 saf[4] = fDZ + newpt[2];
365 saf[5] = fDZ - newpt[2];
366 if (iact < 3 && safe) {
367 smin = saf[0];
368 // compute safe distance
369 for (i = 1; i < 6; i++)
370 if (saf[i] < smin)
371 smin = saf[i];
372 *safe = smin;
373 if (smin < 0)
374 *safe = 0.0;
375 if (iact == 0)
376 return TGeoShape::Big();
377 if (iact == 1 && step < *safe)
378 return TGeoShape::Big();
379 }
380 // compute distance to surface
382 for (i = 0; i < 3; i++) {
383 if (dir[i] != 0) {
384 s = (dir[i] > 0) ? (saf[(i << 1) + 1] / dir[i]) : (-saf[i << 1] / dir[i]);
385 if (s < 0)
386 return 0.0;
387 if (s < smin)
388 smin = s;
389 }
390 }
391 return smin;
392}
393
394////////////////////////////////////////////////////////////////////////////////
395/// Compute distance from inside point to surface of the box.
396/// Boundary safe algorithm.
397
399 const Double_t *origin, Double_t /*stepmax*/)
400{
401 Double_t s, smin, saf[6];
402 Double_t newpt[3];
403 Int_t i;
404 for (i = 0; i < 3; i++)
405 newpt[i] = point[i] - origin[i];
406 saf[0] = dx + newpt[0];
407 saf[1] = dx - newpt[0];
408 saf[2] = dy + newpt[1];
409 saf[3] = dy - newpt[1];
410 saf[4] = dz + newpt[2];
411 saf[5] = dz - newpt[2];
412 // compute distance to surface
414 for (i = 0; i < 3; i++) {
415 if (dir[i] != 0) {
416 s = (dir[i] > 0) ? (saf[(i << 1) + 1] / dir[i]) : (-saf[i << 1] / dir[i]);
417 if (s < 0)
418 return 0.0;
419 if (s < smin)
420 smin = s;
421 }
422 }
423 return smin;
424}
425
426////////////////////////////////////////////////////////////////////////////////
427/// Compute distance from outside point to surface of the box.
428/// Boundary safe algorithm.
429
432{
433 Bool_t in = kTRUE;
434 Double_t saf[3];
435 Double_t par[3];
436 Double_t newpt[3];
437 Int_t i, j;
438 for (i = 0; i < 3; i++)
439 newpt[i] = point[i] - fOrigin[i];
440 par[0] = fDX;
441 par[1] = fDY;
442 par[2] = fDZ;
443 for (i = 0; i < 3; i++) {
444 saf[i] = TMath::Abs(newpt[i]) - par[i];
445 if (saf[i] >= step)
446 return TGeoShape::Big();
447 if (in && saf[i] > 0)
448 in = kFALSE;
449 }
450 if (iact < 3 && safe) {
451 // compute safe distance
452 if (in) {
453 *safe = 0.0;
454 } else {
455 *safe = saf[0];
456 if (saf[1] > *safe)
457 *safe = saf[1];
458 if (saf[2] > *safe)
459 *safe = saf[2];
460 }
461 if (iact == 0)
462 return TGeoShape::Big();
463 if (iact == 1 && step < *safe)
464 return TGeoShape::Big();
465 }
466 // compute distance from point to box
467 // protection in case point is actually inside box
468 if (in) {
469 j = 0;
470 Double_t ss = saf[0];
471 if (saf[1] > ss) {
472 ss = saf[1];
473 j = 1;
474 }
475 if (saf[2] > ss)
476 j = 2;
477 if (newpt[j] * dir[j] > 0)
478 return TGeoShape::Big(); // in fact exiting
479 return 0.0;
480 }
481 for (i = 0; i < 3; i++) {
482 if (saf[i] < 0)
483 continue;
484 if (newpt[i] * dir[i] >= 0)
485 continue;
486 Double_t snxt = saf[i] / TMath::Abs(dir[i]);
487 Int_t ibreak = 0;
488 for (j = 0; j < 3; j++) {
489 if (j == i)
490 continue;
491 Double_t coord = newpt[j] + snxt * dir[j];
492 if (TMath::Abs(coord) > par[j]) {
493 ibreak = 1;
494 break;
495 }
496 }
497 if (!ibreak)
498 return snxt;
499 }
500 return TGeoShape::Big();
501}
502
503////////////////////////////////////////////////////////////////////////////////
504/// Compute distance from outside point to surface of the box.
505/// Boundary safe algorithm.
506
509{
510 Bool_t in = kTRUE;
511 Double_t saf[3];
512 Double_t par[3];
513 Double_t newpt[3];
514 Int_t i, j;
515 for (i = 0; i < 3; i++)
516 newpt[i] = point[i] - origin[i];
517 par[0] = dx;
518 par[1] = dy;
519 par[2] = dz;
520 for (i = 0; i < 3; i++) {
521 saf[i] = TMath::Abs(newpt[i]) - par[i];
522 if (saf[i] >= stepmax)
523 return TGeoShape::Big();
524 if (in && saf[i] > 0)
525 in = kFALSE;
526 }
527 // In case point is inside return ZERO
528 if (in)
529 return 0.0;
531 Int_t ibreak = 0;
532 for (i = 0; i < 3; i++) {
533 if (saf[i] < 0)
534 continue;
535 if (newpt[i] * dir[i] >= 0)
536 continue;
537 snxt = saf[i] / TMath::Abs(dir[i]);
538 ibreak = 0;
539 for (j = 0; j < 3; j++) {
540 if (j == i)
541 continue;
542 coord = newpt[j] + snxt * dir[j];
543 if (TMath::Abs(coord) > par[j]) {
544 ibreak = 1;
545 break;
546 }
547 }
548 if (!ibreak)
549 return snxt;
550 }
551 return TGeoShape::Big();
552}
553
554////////////////////////////////////////////////////////////////////////////////
555/// Returns name of axis IAXIS.
556
558{
559 switch (iaxis) {
560 case 1: return "X";
561 case 2: return "Y";
562 case 3: return "Z";
563 default: return "UNDEFINED";
564 }
565}
566
567////////////////////////////////////////////////////////////////////////////////
568/// Get range of shape for a given axis.
569
571{
572 xlo = 0;
573 xhi = 0;
574 Double_t dx = 0;
575 switch (iaxis) {
576 case 1:
577 xlo = fOrigin[0] - fDX;
578 xhi = fOrigin[0] + fDX;
579 dx = 2 * fDX;
580 return dx;
581 case 2:
582 xlo = fOrigin[1] - fDY;
583 xhi = fOrigin[1] + fDY;
584 dx = 2 * fDY;
585 return dx;
586 case 3:
587 xlo = fOrigin[2] - fDZ;
588 xhi = fOrigin[2] + fDZ;
589 dx = 2 * fDZ;
590 return dx;
591 }
592 return dx;
593}
594
595////////////////////////////////////////////////////////////////////////////////
596/// Fill vector param[4] with the bounding cylinder parameters. The order
597/// is the following : Rmin, Rmax, Phi1, Phi2
598
600{
601 param[0] = 0.; // Rmin
602 param[1] = fDX * fDX + fDY * fDY; // Rmax
603 param[2] = 0.; // Phi1
604 param[3] = 360.; // Phi2
605}
606
607////////////////////////////////////////////////////////////////////////////////
608/// Get area in internal units of the facet with a given index.
609/// Possible index values:
610/// - 0 - all facets together
611/// - 1 to 6 - facet index from bottom to top Z
612
614{
615 Double_t area = 0.;
616 switch (index) {
617 case 0: area = 8. * (fDX * fDY + fDX * fDZ + fDY * fDZ); return area;
618 case 1:
619 case 6: area = 4. * fDX * fDY; return area;
620 case 2:
621 case 4: area = 4. * fDX * fDZ; return area;
622 case 3:
623 case 5: area = 4. * fDY * fDZ; return area;
624 }
625 return area;
626}
627
628////////////////////////////////////////////////////////////////////////////////
629/// Fills array with n random points located on the surface of indexed facet.
630/// The output array must be provided with a length of minimum 3*npoints. Returns
631/// true if operation succeeded.
632/// Possible index values:
633/// - 0 - all facets together
634/// - 1 to 6 - facet index from bottom to top Z
635
637{
639 return kFALSE;
640 Double_t surf[6];
641 Double_t area = 0.;
642 if (index == 0) {
643 for (Int_t isurf = 0; isurf < 6; isurf++) {
645 if (isurf > 0)
646 surf[isurf] += surf[isurf - 1];
647 }
648 area = surf[5];
649 }
650
651 for (Int_t i = 0; i < npoints; i++) {
652 // Generate randomly a surface index if needed.
653 Double_t *point = &array[3 * i];
655 if (surfindex == 0) {
656 Double_t val = area * gRandom->Rndm();
657 surfindex = 2 + TMath::BinarySearch(6, surf, val);
658 if (surfindex > 6)
659 surfindex = 6;
660 }
661 switch (surfindex) {
662 case 1:
663 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
664 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
665 point[2] = -fDZ;
666 break;
667 case 2:
668 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
669 point[1] = -fDY;
670 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
671 break;
672 case 3:
673 point[0] = -fDX;
674 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
675 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
676 break;
677 case 4:
678 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
679 point[1] = fDY;
680 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
681 break;
682 case 5:
683 point[0] = fDX;
684 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
685 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
686 break;
687 case 6:
688 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
689 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
690 point[2] = fDZ;
691 break;
692 }
693 }
694 return kTRUE;
695}
696
697////////////////////////////////////////////////////////////////////////////////
698/// Fills array with n random points located on the line segments of the shape mesh.
699/// The output array must be provided with a length of minimum 3*npoints. Returns
700/// true if operation is implemented.
701
703{
704 if (!array || npoints <= 0)
705 return kFALSE;
706
707 auto buff = this->MakeBuffer3D();
708 if (!buff)
709 return kFALSE;
710
711 const Int_t npnts = (Int_t)buff->NbPnts();
712 const Int_t nsegs = (Int_t)buff->NbSegs();
713 if (npnts <= 0 || nsegs <= 0 || !buff->fPnts || !buff->fSegs) {
714 delete buff;
715 return kFALSE;
716 }
717
718 // Fill only points on segments into `array`
720 Int_t icrt = 0; // <--- start at 0, output-only buffer
721
723
724 for (Int_t i = 0; i < nsegs && remaining > 0; i++) {
725 const Int_t i0 = buff->fSegs[3 * i + 1];
726 const Int_t i1 = buff->fSegs[3 * i + 2];
727
729 delete buff;
730 Error("GetPointsOnSegments", "Segment index out of range: (%d,%d) npnts=%d", i0, i1, npnts);
731 return kFALSE;
732 }
733 const Double_t *p0 = &buff->fPnts[3 * i0];
734 const Double_t *p1 = &buff->fPnts[3 * i1];
735
736 if (i == (nsegs - 1))
738
739 const Double_t dx = (p1[0] - p0[0]) / (nperseg + 1);
740 const Double_t dy = (p1[1] - p0[1]) / (nperseg + 1);
741 const Double_t dz = (p1[2] - p0[2]) / (nperseg + 1);
742
743 for (Int_t j = 0; j < nperseg && remaining > 0; ++j) {
744 array[icrt++] = p0[0] + (j + 1) * dx;
745 array[icrt++] = p0[1] + (j + 1) * dy;
746 array[icrt++] = p0[2] + (j + 1) * dz;
747 --remaining;
748 }
749 }
750 delete buff;
751
752 if (remaining > 0) {
753 Error("GetPointsOnSegments", "Could not fill %d points (%d missing)", npoints, remaining);
754 return kFALSE;
755 }
756 return kTRUE;
757}
758
759////////////////////////////////////////////////////////////////////////////////
760/// Fills real parameters of a positioned box inside this one. Returns 0 if successful.
761
763{
764 dx = dy = dz = 0;
765 if (mat->IsRotation()) {
766 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
767 return 1; // ### rotation not accepted ###
768 }
769 //--> translate the origin of the parametrized box to the frame of this box.
770 Double_t origin[3];
771 mat->LocalToMaster(parambox->GetOrigin(), origin);
772 if (!Contains(origin)) {
773 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
774 return 1; // ### wrong matrix ###
775 }
776 //--> now we have to get the valid range for all parametrized axis
777 Double_t xlo = 0, xhi = 0;
778 Double_t dd[3];
779 dd[0] = parambox->GetDX();
780 dd[1] = parambox->GetDY();
781 dd[2] = parambox->GetDZ();
782 for (Int_t iaxis = 0; iaxis < 3; iaxis++) {
783 if (dd[iaxis] >= 0)
784 continue;
785 TGeoBBox::GetAxisRange(iaxis + 1, xlo, xhi);
786 //-> compute best fitting parameter
787 dd[iaxis] = TMath::Min(origin[iaxis] - xlo, xhi - origin[iaxis]);
788 if (dd[iaxis] < 0) {
789 Error("GetFittingBox", "wrong matrix");
790 return 1;
791 }
792 }
793 dx = dd[0];
794 dy = dd[1];
795 dz = dd[2];
796 return 0;
797}
798
799////////////////////////////////////////////////////////////////////////////////
800/// In case shape has some negative parameters, these has to be computed
801/// in order to fit the mother
802
804{
806 return nullptr;
807 Double_t dx, dy, dz;
808 Int_t ierr = mother->GetFittingBox(this, mat, dx, dy, dz);
809 if (ierr) {
810 Error("GetMakeRuntimeShape", "cannot fit this to mother");
811 return nullptr;
812 }
813 return (new TGeoBBox(dx, dy, dz));
814}
815
816////////////////////////////////////////////////////////////////////////////////
817/// Returns numbers of vertices, segments and polygons composing the shape mesh.
818
820{
821 nvert = 8;
822 nsegs = 12;
823 npols = 6;
824}
825
826////////////////////////////////////////////////////////////////////////////////
827/// @brief Compute the world-space center of a placed TGeoBBox.
828///
829/// A TGeoBBox does not necessarily have its center at local coordinates (0,0,0).
830/// Instead, its local center is given by TGeoBBox::GetOrigin().
831///
832/// @param[in] m Pointer to the placement TGeoMatrix.
833/// @return World-space center of the box.
834
836{
837 Double_t w[3];
838 m->LocalToMaster(fOrigin, w);
839 return ROOT::Math::XYZVector(w[0], w[1], w[2]);
840}
841
842////////////////////////////////////////////////////////////////////////////////
843/// Prints shape parameters
844
846{
847 printf("*** Shape %s: TGeoBBox ***\n", GetName());
848 printf(" dX = %11.5f\n", fDX);
849 printf(" dY = %11.5f\n", fDY);
850 printf(" dZ = %11.5f\n", fDZ);
851 printf(" origin: x=%11.5f y=%11.5f z=%11.5f\n", fOrigin[0], fOrigin[1], fOrigin[2]);
852}
853
854////////////////////////////////////////////////////////////////////////////////
855/// @brief Test whether a given axis is a separating axis for two oriented boxes.
856///
857/// The test is based on the Separating Axis Theorem (SAT). The boxes are projected
858/// onto the axis L using a center-plus-radius formulation, avoiding explicit
859/// vertex enumeration.
860///
861/// Tolerance handling:
862/// Two boxes that only touch (face, edge, or corner contact) within
863/// TGeoShape::Tolerance() are considered NON-intersecting.
864///
865/// Formally, the axis separates the boxes if:
866///
867/// |(oB - oA) · L| >= rA(L) + rB(L) - tol
868///
869/// where rA and rB are the projection radii of each box onto L.
870///
871/// @param[in] L Candidate separating axis in world coordinates.
872/// @param[in] D Vector from center of box A to center of box B (oB - oA).
873/// @param[in] Ax World-space X axis of box A.
874/// @param[in] Ay World-space Y axis of box A.
875/// @param[in] Az World-space Z axis of box A.
876/// @param[in] Bx World-space X axis of box B.
877/// @param[in] By World-space Y axis of box B.
878/// @param[in] Bz World-space Z axis of box B.
879/// @param[in] dA Half-length vector of box A.
880/// @param[in] dB Half-length vector of box B.
881/// @param[in] tol Geometric tolerance (typically TGeoShape::Tolerance()).
882///
883/// @return true If L is a separating axis (boxes do NOT intersect).
884/// @return false If L does not separate the boxes or is degenerate.
885
891{
892 // Degenerate axes can arise from cross products of nearly parallel axes.
893 const Double_t eps = 1e-18;
894 if (L.Mag2() < eps)
895 return false;
896
897 const Double_t dist = std::fabs(D.Dot(L));
898 const Double_t rA = TMath::Abs(dA.x() * Ax.Dot(L)) + TMath::Abs(dA.y() * Ay.Dot(L)) + TMath::Abs(dA.z() * Az.Dot(L));
899 const Double_t rB = TMath::Abs(dB.x() * Bx.Dot(L)) + TMath::Abs(dB.y() * By.Dot(L)) + TMath::Abs(dB.z() * Bz.Dot(L));
900
901 // Touching within tolerance is treated as non-intersecting
902 return dist >= (rA + rB - tol);
903}
904
905////////////////////////////////////////////////////////////////////////////////
906/// Creates a TBuffer3D describing *this* shape.
907/// Coordinates are in local reference frame.
908
910{
911 Int_t npnts = 0, nsegs = 0, npols = 0;
912 this->GetMeshNumbers(npnts, nsegs, npols);
913 if (npnts <= 0 || nsegs <= 0)
914 return nullptr;
916 this->SetPoints(buff->fPnts);
917 this->SetSegsAndPols(*buff);
918 return buff;
919}
920
921////////////////////////////////////////////////////////////////////////////////
922/// @brief Fast "may-intersect" test for two placed TGeoBBox objects.
923///
924/// This method implements a broad-phase oriented bounding box (OBB) intersection
925/// test using the Separating Axis Theorem (SAT).
926///
927/// The test is optimized for fast rejection:
928/// - bounding-sphere early exit,
929/// - face-normal axes tested first,
930/// - cross-product axes tested only if needed.
931///
932/// Geometry and tolerance semantics:
933/// - Boxes are treated as oriented bounding boxes defined by their half-lengths
934/// and placement matrices.
935/// - TGeoBBox::GetOrigin() is fully taken into account.
936/// - Boxes that only touch within TGeoShape::Tolerance() are considered
937/// NON-intersecting.
938///
939/// Intended usage:
940/// - Broad-phase culling in navigation or overlap checks.
941/// - Conservative test: returning true means the boxes MAY intersect.
942///
943/// @param[in] boxA Pointer to the first TGeoBBox shape.
944/// @param[in] mA Placement matrix of the first box.
945/// @param[in] boxB Pointer to the second TGeoBBox shape.
946/// @param[in] mB Placement matrix of the second box.
947///
948/// @return true If the boxes MAY intersect.
949/// @return false If the boxes are guaranteed NOT to intersect
950/// (or only touch within tolerance).
951
953{
955
956 // Half-lengths vectors
957 const ROOT::Math::XYZVector dA = boxA->GetDimensions();
958 const ROOT::Math::XYZVector dB = boxB->GetDimensions();
959
960 // Correct world-space centers (including TGeoBBox origin)
961 const ROOT::Math::XYZVector oA = boxA->GetWorldCenter(mA);
962 const ROOT::Math::XYZVector oB = boxB->GetWorldCenter(mB);
963 const ROOT::Math::XYZVector D = oB - oA;
964
965 // Very cheap early rejection using bounding spheres
966 const Double_t rA = TMath::Sqrt(dA.Mag2());
967 const Double_t rB = TMath::Sqrt(dB.Mag2());
968 const Double_t dmax = rA + rB - tol;
969
970 // Touching within tolerance is treated as non-intersecting
971 if (D.Mag2() >= dmax * dmax)
972 return kFALSE;
973
974 // World-space box axes
977 mA->GetWorldAxes(Ax, Ay, Az);
978 mB->GetWorldAxes(Bx, By, Bz);
979
980 // SAT: face-normal axes (6 tests)
981 if (IsSeparatingAxis(Ax, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
982 return kFALSE;
983 if (IsSeparatingAxis(Ay, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
984 return kFALSE;
985 if (IsSeparatingAxis(Az, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
986 return kFALSE;
987
988 if (IsSeparatingAxis(Bx, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
989 return kFALSE;
990 if (IsSeparatingAxis(By, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
991 return kFALSE;
992 if (IsSeparatingAxis(Bz, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
993 return kFALSE;
994
995 // SAT: cross-product axes (9 tests)
996 const ROOT::Math::XYZVector Aaxes[3] = {Ax, Ay, Az};
997 const ROOT::Math::XYZVector Baxes[3] = {Bx, By, Bz};
998
999 for (int i = 0; i < 3; ++i) {
1000 for (int j = 0; j < 3; ++j) {
1001 const ROOT::Math::XYZVector L = Aaxes[i].Cross(Baxes[j]);
1002 if (IsSeparatingAxis(L, D, Ax, Ay, Az, Bx, By, Bz, dA, dB, tol))
1003 return kFALSE;
1004 }
1005 }
1006
1007 // No separating axis found: boxes MAY intersect
1008 return kTRUE;
1009}
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// Fills TBuffer3D structure for segments and polygons.
1013
1015{
1016 Int_t c = GetBasicColor();
1017
1018 buff.fSegs[0] = c;
1019 buff.fSegs[1] = 0;
1020 buff.fSegs[2] = 1;
1021 buff.fSegs[3] = c + 1;
1022 buff.fSegs[4] = 1;
1023 buff.fSegs[5] = 2;
1024 buff.fSegs[6] = c + 1;
1025 buff.fSegs[7] = 2;
1026 buff.fSegs[8] = 3;
1027 buff.fSegs[9] = c;
1028 buff.fSegs[10] = 3;
1029 buff.fSegs[11] = 0;
1030 buff.fSegs[12] = c + 2;
1031 buff.fSegs[13] = 4;
1032 buff.fSegs[14] = 5;
1033 buff.fSegs[15] = c + 2;
1034 buff.fSegs[16] = 5;
1035 buff.fSegs[17] = 6;
1036 buff.fSegs[18] = c + 3;
1037 buff.fSegs[19] = 6;
1038 buff.fSegs[20] = 7;
1039 buff.fSegs[21] = c + 3;
1040 buff.fSegs[22] = 7;
1041 buff.fSegs[23] = 4;
1042 buff.fSegs[24] = c;
1043 buff.fSegs[25] = 0;
1044 buff.fSegs[26] = 4;
1045 buff.fSegs[27] = c + 2;
1046 buff.fSegs[28] = 1;
1047 buff.fSegs[29] = 5;
1048 buff.fSegs[30] = c + 1;
1049 buff.fSegs[31] = 2;
1050 buff.fSegs[32] = 6;
1051 buff.fSegs[33] = c + 3;
1052 buff.fSegs[34] = 3;
1053 buff.fSegs[35] = 7;
1054
1055 buff.fPols[0] = c;
1056 buff.fPols[1] = 4;
1057 buff.fPols[2] = 0;
1058 buff.fPols[3] = 9;
1059 buff.fPols[4] = 4;
1060 buff.fPols[5] = 8;
1061 buff.fPols[6] = c + 1;
1062 buff.fPols[7] = 4;
1063 buff.fPols[8] = 1;
1064 buff.fPols[9] = 10;
1065 buff.fPols[10] = 5;
1066 buff.fPols[11] = 9;
1067 buff.fPols[12] = c;
1068 buff.fPols[13] = 4;
1069 buff.fPols[14] = 2;
1070 buff.fPols[15] = 11;
1071 buff.fPols[16] = 6;
1072 buff.fPols[17] = 10;
1073 buff.fPols[18] = c + 1;
1074 buff.fPols[19] = 4;
1075 buff.fPols[20] = 3;
1076 buff.fPols[21] = 8;
1077 buff.fPols[22] = 7;
1078 buff.fPols[23] = 11;
1079 buff.fPols[24] = c + 2;
1080 buff.fPols[25] = 4;
1081 buff.fPols[26] = 0;
1082 buff.fPols[27] = 3;
1083 buff.fPols[28] = 2;
1084 buff.fPols[29] = 1;
1085 buff.fPols[30] = c + 3;
1086 buff.fPols[31] = 4;
1087 buff.fPols[32] = 4;
1088 buff.fPols[33] = 5;
1089 buff.fPols[34] = 6;
1090 buff.fPols[35] = 7;
1091}
1092
1093////////////////////////////////////////////////////////////////////////////////
1094/// Computes the closest distance from given point to this shape.
1095
1097{
1099 if (in) {
1100 safe = fDX - TMath::Abs(point[0] - fOrigin[0]);
1101 safy = fDY - TMath::Abs(point[1] - fOrigin[1]);
1102 safz = fDZ - TMath::Abs(point[2] - fOrigin[2]);
1103 if (safy < safe)
1104 safe = safy;
1105 if (safz < safe)
1106 safe = safz;
1107 } else {
1108 safe = -fDX + TMath::Abs(point[0] - fOrigin[0]);
1109 safy = -fDY + TMath::Abs(point[1] - fOrigin[1]);
1110 safz = -fDZ + TMath::Abs(point[2] - fOrigin[2]);
1111 if (safy > safe)
1112 safe = safy;
1113 if (safz > safe)
1114 safe = safz;
1115 }
1116 return safe;
1117}
1118
1119////////////////////////////////////////////////////////////////////////////////
1120/// Save a primitive as a C++ statement(s) on output stream "out".
1121
1122void TGeoBBox::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
1123{
1125 return;
1126 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
1127 out << " dx = " << fDX << ";" << std::endl;
1128 out << " dy = " << fDY << ";" << std::endl;
1129 out << " dz = " << fDZ << ";" << std::endl;
1132 out << " origin[0] = " << fOrigin[0] << ";" << std::endl;
1133 out << " origin[1] = " << fOrigin[1] << ";" << std::endl;
1134 out << " origin[2] = " << fOrigin[2] << ";" << std::endl;
1135 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz,origin);"
1136 << std::endl;
1137 } else {
1138 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz);" << std::endl;
1139 }
1141}
1142
1143////////////////////////////////////////////////////////////////////////////////
1144/// Set parameters of the box.
1145
1147{
1148 fDX = dx;
1149 fDY = dy;
1150 fDZ = dz;
1151 if (origin) {
1152 fOrigin[0] = origin[0];
1153 fOrigin[1] = origin[1];
1154 fOrigin[2] = origin[2];
1155 }
1158 return;
1159 if ((fDX < 0) || (fDY < 0) || (fDZ < 0))
1161}
1162
1163////////////////////////////////////////////////////////////////////////////////
1164/// Set dimensions based on the array of parameters
1165/// param[0] - half-length in x
1166/// param[1] - half-length in y
1167/// param[2] - half-length in z
1168
1170{
1171 if (!param) {
1172 Error("SetDimensions", "null parameters");
1173 return;
1174 }
1175 fDX = param[0];
1176 fDY = param[1];
1177 fDZ = param[2];
1180 return;
1181 if ((fDX < 0) || (fDY < 0) || (fDZ < 0))
1183}
1184
1185////////////////////////////////////////////////////////////////////////////////
1186/// Fill box vertices to an array.
1187
1192
1193////////////////////////////////////////////////////////////////////////////////
1194/// Fill box points.
1195
1197{
1198 if (!points)
1199 return;
1200 Double_t xmin, xmax, ymin, ymax, zmin, zmax;
1201 xmin = -fDX + fOrigin[0];
1202 xmax = fDX + fOrigin[0];
1203 ymin = -fDY + fOrigin[1];
1204 ymax = fDY + fOrigin[1];
1205 zmin = -fDZ + fOrigin[2];
1206 zmax = fDZ + fOrigin[2];
1207 points[0] = xmin;
1208 points[1] = ymin;
1209 points[2] = zmin;
1210 points[3] = xmin;
1211 points[4] = ymax;
1212 points[5] = zmin;
1213 points[6] = xmax;
1214 points[7] = ymax;
1215 points[8] = zmin;
1216 points[9] = xmax;
1217 points[10] = ymin;
1218 points[11] = zmin;
1219 points[12] = xmin;
1220 points[13] = ymin;
1221 points[14] = zmax;
1222 points[15] = xmin;
1223 points[16] = ymax;
1224 points[17] = zmax;
1225 points[18] = xmax;
1226 points[19] = ymax;
1227 points[20] = zmax;
1228 points[21] = xmax;
1229 points[22] = ymin;
1230 points[23] = zmax;
1231}
1232
1233////////////////////////////////////////////////////////////////////////////////
1234/// Fill box points.
1235
1237{
1238 if (!points)
1239 return;
1240 Double_t xmin, xmax, ymin, ymax, zmin, zmax;
1241 xmin = -fDX + fOrigin[0];
1242 xmax = fDX + fOrigin[0];
1243 ymin = -fDY + fOrigin[1];
1244 ymax = fDY + fOrigin[1];
1245 zmin = -fDZ + fOrigin[2];
1246 zmax = fDZ + fOrigin[2];
1247 points[0] = xmin;
1248 points[1] = ymin;
1249 points[2] = zmin;
1250 points[3] = xmin;
1251 points[4] = ymax;
1252 points[5] = zmin;
1253 points[6] = xmax;
1254 points[7] = ymax;
1255 points[8] = zmin;
1256 points[9] = xmax;
1257 points[10] = ymin;
1258 points[11] = zmin;
1259 points[12] = xmin;
1260 points[13] = ymin;
1261 points[14] = zmax;
1262 points[15] = xmin;
1263 points[16] = ymax;
1264 points[17] = zmax;
1265 points[18] = xmax;
1266 points[19] = ymax;
1267 points[20] = zmax;
1268 points[21] = xmax;
1269 points[22] = ymin;
1270 points[23] = zmax;
1271}
1272
1273////////////////////////////////////////////////////////////////////////////////
1274////// fill size of this 3-D object
1275//// TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
1276//// if (painter) painter->AddSize3D(8, 12, 6);
1277
1278void TGeoBBox::Sizeof3D() const {}
1279
1280////////////////////////////////////////////////////////////////////////////////
1281/// Fills a static 3D buffer and returns a reference.
1282
1284{
1285 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1286
1288
1289 // TODO: A box itself has has nothing more as already described
1290 // by bounding box. How will viewer interpret?
1292 if (buffer.SetRawSizes(8, 3 * 8, 12, 3 * 12, 6, 6 * 6)) {
1294 }
1295 }
1297 SetPoints(buffer.fPnts);
1298 if (!buffer.fLocalFrame) {
1299 TransformPoints(buffer.fPnts, buffer.NbPnts());
1300 }
1301
1302 SetSegsAndPols(buffer);
1304 }
1305
1306 return buffer;
1307}
1308
1309////////////////////////////////////////////////////////////////////////////////
1310/// Fills the supplied buffer, with sections in desired frame
1311/// See TBuffer3D.h for explanation of sections, frame etc.
1312
1314{
1316
1318 Double_t halfLengths[3] = {fDX, fDY, fDZ};
1320
1321 if (!buffer.fLocalFrame) {
1322 TransformPoints(buffer.fBBVertex[0], 8);
1323 }
1325 }
1326}
1327
1328////////////////////////////////////////////////////////////////////////////////
1329/// Check the inside status for each of the points in the array.
1330/// Input: Array of point coordinates + vector size
1331/// Output: Array of Booleans for the inside of each point
1332
1334{
1335 for (Int_t i = 0; i < vecsize; i++)
1336 inside[i] = Contains(&points[3 * i]);
1337}
1338
1339////////////////////////////////////////////////////////////////////////////////
1340/// Compute the normal for an array o points so that norm.dot.dir is positive
1341/// Input: Arrays of point coordinates and directions + vector size
1342/// Output: Array of normal directions
1343
1345{
1346 for (Int_t i = 0; i < vecsize; i++)
1347 ComputeNormal(&points[3 * i], &dirs[3 * i], &norms[3 * i]);
1348}
1349
1350////////////////////////////////////////////////////////////////////////////////
1351/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1352
1354 Double_t *step) const
1355{
1356 for (Int_t i = 0; i < vecsize; i++)
1357 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1358}
1359
1360////////////////////////////////////////////////////////////////////////////////
1361/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1362
1364 Double_t *step) const
1365{
1366 for (Int_t i = 0; i < vecsize; i++)
1367 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1368}
1369
1370////////////////////////////////////////////////////////////////////////////////
1371/// Compute safe distance from each of the points in the input array.
1372/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1373/// Output: Safety values
1374
1376{
1377 for (Int_t i = 0; i < vecsize; i++)
1378 safe[i] = Safety(&points[3 * i], inside[i]);
1379}
#define c(i)
Definition RSha256.hxx:101
#define e(i)
Definition RSha256.hxx:103
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float ymin
float xmax
float ymax
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
Scalar Mag2() const
Magnitute squared ( r^2 in spherical coordinate)
Scalar Dot(const DisplacementVector3D< OtherCoords, Tag > &v) const
Return the scalar (dot) product of two displacement vectors.
Generic 3D primitive description class.
Definition TBuffer3D.h:18
UInt_t NbPnts() const
Definition TBuffer3D.h:80
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
@ kBoundingBox
Definition TBuffer3D.h:51
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Bool_t fLocalFrame
Definition TBuffer3D.h:90
void SetAABoundingBox(const Double_t origin[3], const Double_t halfLengths[3])
Set fBBVertex in kBoundingBox section to a axis aligned (local) BB using supplied origin and box half...
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:113
Double_t fBBVertex[8][3]
Definition TBuffer3D.h:108
Box class.
Definition TGeoBBox.h:18
ROOT::Math::XYZVector GetWorldCenter(const TGeoMatrix *m) const
Compute the world-space center of a placed TGeoBBox.
Definition TGeoBBox.cxx:835
void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const override
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
Compute distance from inside point to surface of the box.
Definition TGeoBBox.cxx:353
Double_t fDX
Definition TGeoBBox.h:21
const char * GetAxisName(Int_t iaxis) const override
Returns name of axis IAXIS.
Definition TGeoBBox.cxx:557
void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Compute distance from array of input points having directions specified by dirs. Store output in dist...
void GetBoundingCylinder(Double_t *param) const override
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoBBox.cxx:599
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) const override
Computes normal to closest surface from POINT.
Definition TGeoBBox.cxx:218
virtual Bool_t GetPointsOnFacet(Int_t index, Int_t npoints, Double_t *array) const
Fills array with n random points located on the surface of indexed facet.
Definition TGeoBBox.cxx:636
void SetPoints(Double_t *points) const override
Fill box points.
void SetBoxDimensions(Double_t dx, Double_t dy, Double_t dz, Double_t *origin=nullptr)
Set parameters of the box.
Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
Compute distance from outside point to surface of the box.
Definition TGeoBBox.cxx:431
Double_t Capacity() const override
Computes capacity of the shape in [length^3].
Definition TGeoBBox.cxx:210
virtual Double_t GetFacetArea(Int_t index=0) const
Get area in internal units of the facet with a given index.
Definition TGeoBBox.cxx:613
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute closest distance from point px,py to each corner.
Definition TGeoBBox.cxx:263
Bool_t CouldBeCrossed(const Double_t *point, const Double_t *dir) const override
Decides fast if the bounding box could be crossed by a vector.
Definition TGeoBBox.cxx:233
void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const override
Check the inside status for each of the points in the array.
void SetBoxPoints(Double_t *points) const
Fill box vertices to an array.
Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const override
Fills array with n random points located on the line segments of the shape mesh.
Definition TGeoBBox.cxx:702
Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const override
Fills real parameters of a positioned box inside this one. Returns 0 if successful.
Definition TGeoBBox.cxx:762
Double_t fOrigin[3]
Definition TGeoBBox.h:24
Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const override
Get range of shape for a given axis.
Definition TGeoBBox.cxx:570
TBuffer3D * MakeBuffer3D() const override
Creates a TBuffer3D describing this shape.
Definition TGeoBBox.cxx:909
static Bool_t IsSeparatingAxis(const ROOT::Math::XYZVector &L, const ROOT::Math::XYZVector &D, const ROOT::Math::XYZVector &Ax, const ROOT::Math::XYZVector &Ay, const ROOT::Math::XYZVector &Az, const ROOT::Math::XYZVector &Bx, const ROOT::Math::XYZVector &By, const ROOT::Math::XYZVector &Bz, const ROOT::Math::XYZVector &dA, const ROOT::Math::XYZVector &dB, Double_t tol)
Test whether a given axis is a separating axis for two oriented boxes.
Definition TGeoBBox.cxx:886
void InspectShape() const override
Prints shape parameters.
Definition TGeoBBox.cxx:845
static bool MayIntersect(const TGeoBBox *boxA, const TGeoMatrix *mA, const TGeoBBox *boxB, const TGeoMatrix *mB)
Fast "may-intersect" test for two placed TGeoBBox objects.
Definition TGeoBBox.cxx:952
void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Compute distance from array of input points having directions specified by dirs. Store output in dist...
Bool_t Contains(const Double_t *point) const override
Test if point is inside this shape.
Definition TGeoBBox.cxx:322
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
Computes the closest distance from given point to this shape.
void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const override
Returns numbers of vertices, segments and polygons composing the shape mesh.
Definition TGeoBBox.cxx:819
void Sizeof3D() const override
TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const override
In case shape has some negative parameters, these has to be computed in order to fit the mother.
Definition TGeoBBox.cxx:803
Double_t fDY
Definition TGeoBBox.h:22
void SetSegsAndPols(TBuffer3D &buffer) const override
Fills TBuffer3D structure for segments and polygons.
void SetDimensions(Double_t *param) override
Set dimensions based on the array of parameters param[0] - half-length in x param[1] - half-length in...
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const override
Fills a static 3D buffer and returns a reference.
Double_t fDZ
Definition TGeoBBox.h:23
void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize) override
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
~TGeoBBox() override
Destructor.
Definition TGeoBBox.cxx:205
void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const override
Compute safe distance from each of the points in the input array.
void ComputeBBox() override
Compute bounding box - nothing to do in this case.
Definition TGeoBBox.cxx:317
TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) override
Divide this box shape belonging to volume "voldiv" into ndiv equal volumes called divname,...
Definition TGeoBBox.cxx:276
TGeoBBox()
Default constructor.
Definition TGeoBBox.cxx:162
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Geometrical transformation package.
Definition TGeoMatrix.h:39
Node containing an offset.
Definition TGeoNode.h:185
base finder class for patterns. A pattern is specifying a division type
a Y axis divison pattern
a Z axis divison pattern
Base abstract class for all shapes.
Definition TGeoShape.h:25
static Double_t Big()
Definition TGeoShape.h:95
Int_t GetBasicColor() const
Get the basic color (0-7).
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.
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,...
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
const char * GetName() const override
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoRunTimeShape
Definition TGeoShape.h:40
static Double_t Tolerance()
Definition TGeoShape.h:98
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:177
Volume families.
Definition TGeoVolume.h:267
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:202
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:881
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1088
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:558
Basic string class.
Definition TString.h:138
const char * Data() const
Definition TString.h:384
DisplacementVector3D< Cartesian3D< double >, DefaultCoordinateSystemTag > XYZVector
3D Vector based on the cartesian coordinates x,y,z in double precision
Definition Vector3Dfwd.h:44
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:993
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:249
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:197
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:329
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:122
TMarker m
Definition textangle.C:8