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 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 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 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
160
161////////////////////////////////////////////////////////////////////////////////
162/// Default constructor
163
165{
167 fDX = fDY = fDZ = 0;
168 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
169}
170
171////////////////////////////////////////////////////////////////////////////////
172/// Constructor where half-lengths are provided.
173
175 :TGeoShape("")
176{
178 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
179 SetBoxDimensions(dx, dy, dz, origin);
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Constructor with shape name.
184
185TGeoBBox::TGeoBBox(const char *name, Double_t dx, Double_t dy, Double_t dz, Double_t *origin)
187{
189 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
190 SetBoxDimensions(dx, dy, dz, origin);
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// Constructor based on the array of parameters.
195/// - param[0] - half-length in x
196/// - param[1] - half-length in y
197/// - param[2] - half-length in z
198
200 :TGeoShape("")
201{
203 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
204 SetDimensions(param);
205}
206
207////////////////////////////////////////////////////////////////////////////////
208/// Destructor
209
211{
212}
213
214////////////////////////////////////////////////////////////////////////////////
215/// Check if 2 positioned boxes overlap.
216
217Bool_t TGeoBBox::AreOverlapping(const TGeoBBox *box1, const TGeoMatrix *mat1, const TGeoBBox *box2, const TGeoMatrix *mat2)
218{
219 Double_t master[3];
220 Double_t local[3];
221 Double_t ldir1[3], ldir2[3];
222 const Double_t *o1 = box1->GetOrigin();
223 const Double_t *o2 = box2->GetOrigin();
224 // Convert center of first box to the local frame of second
225 mat1->LocalToMaster(o1, master);
226 mat2->MasterToLocal(master, local);
227 if (TGeoBBox::Contains(local,box2->GetDX(),box2->GetDY(),box2->GetDZ(),o2)) return kTRUE;
228 Double_t distsq = (local[0]-o2[0])*(local[0]-o2[0]) +
229 (local[1]-o2[1])*(local[1]-o2[1]) +
230 (local[2]-o2[2])*(local[2]-o2[2]);
231 // Compute distance between box centers and compare with max value
232 Double_t rmaxsq = (box1->GetDX()+box2->GetDX())*(box1->GetDX()+box2->GetDX()) +
233 (box1->GetDY()+box2->GetDY())*(box1->GetDY()+box2->GetDY()) +
234 (box1->GetDZ()+box2->GetDZ())*(box1->GetDZ()+box2->GetDZ());
235 if (distsq > rmaxsq + TGeoShape::Tolerance()) return kFALSE;
236 // We are still not sure: shoot a ray from the center of "1" towards the
237 // center of 2.
238 Double_t dir[3];
239 mat1->LocalToMaster(o1, ldir1);
240 mat2->LocalToMaster(o2, ldir2);
241 distsq = 1./TMath::Sqrt(distsq);
242 dir[0] = (ldir2[0]-ldir1[0])*distsq;
243 dir[1] = (ldir2[1]-ldir1[1])*distsq;
244 dir[2] = (ldir2[2]-ldir1[2])*distsq;
245 mat1->MasterToLocalVect(dir, ldir1);
246 mat2->MasterToLocalVect(dir, ldir2);
247 // Distance to exit from o1
248 Double_t dist1 = TGeoBBox::DistFromInside(o1,ldir1,box1->GetDX(),box1->GetDY(),box1->GetDZ(),o1);
249 // Distance to enter from o2
250 Double_t dist2 = TGeoBBox::DistFromOutside(local,ldir2,box2->GetDX(),box2->GetDY(),box2->GetDZ(),o2);
251 if (dist1 > dist2) return kTRUE;
252 return kFALSE;
253}
254
255////////////////////////////////////////////////////////////////////////////////
256/// Computes capacity of the shape in [length^3].
257
259{
260 return (8.*fDX*fDY*fDZ);
261}
262
263////////////////////////////////////////////////////////////////////////////////
264/// Computes normal to closest surface from POINT.
265
266void TGeoBBox::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
267{
268 memset(norm,0,3*sizeof(Double_t));
269 Double_t saf[3];
270 Int_t i;
271 saf[0]=TMath::Abs(TMath::Abs(point[0]-fOrigin[0])-fDX);
272 saf[1]=TMath::Abs(TMath::Abs(point[1]-fOrigin[1])-fDY);
273 saf[2]=TMath::Abs(TMath::Abs(point[2]-fOrigin[2])-fDZ);
274 i = TMath::LocMin(3,saf);
275 norm[i] = (dir[i]>0)?1:(-1);
276}
277
278////////////////////////////////////////////////////////////////////////////////
279/// Decides fast if the bounding box could be crossed by a vector.
280
281Bool_t TGeoBBox::CouldBeCrossed(const Double_t *point, const Double_t *dir) const
282{
283 Double_t mind = fDX;
284 if (fDY<mind) mind=fDY;
285 if (fDZ<mind) mind=fDZ;
286 Double_t dx = fOrigin[0]-point[0];
287 Double_t dy = fOrigin[1]-point[1];
288 Double_t dz = fOrigin[2]-point[2];
289 Double_t do2 = dx*dx+dy*dy+dz*dz;
290 if (do2<=(mind*mind)) return kTRUE;
291 Double_t rmax2 = fDX*fDX+fDY*fDY+fDZ*fDZ;
292 if (do2<=rmax2) return kTRUE;
293 // inside bounding sphere
294 Double_t doct = dx*dir[0]+dy*dir[1]+dz*dir[2];
295 // leaving ray
296 if (doct<=0) return kFALSE;
297 Double_t dirnorm=dir[0]*dir[0]+dir[1]*dir[1]+dir[2]*dir[2];
298 if ((doct*doct)>=(do2-rmax2)*dirnorm) return kTRUE;
299 return kFALSE;
300}
301
302////////////////////////////////////////////////////////////////////////////////
303/// Compute closest distance from point px,py to each corner.
304
306{
307 const Int_t numPoints = 8;
308 return ShapeDistancetoPrimitive(numPoints, px, py);
309}
310
311////////////////////////////////////////////////////////////////////////////////
312/// Divide this box shape belonging to volume "voldiv" into ndiv equal volumes
313/// called divname, from start position with the given step. Returns pointer
314/// to created division cell volume. In case a wrong division axis is supplied,
315/// returns pointer to volume to be divided.
316
317TGeoVolume *TGeoBBox::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
318 Double_t start, Double_t step)
319{
320 TGeoShape *shape; //--- shape to be created
321 TGeoVolume *vol; //--- division volume to be created
322 TGeoVolumeMulti *vmulti; //--- generic divided volume
323 TGeoPatternFinder *finder; //--- finder to be attached
324 TString opt = ""; //--- option to be attached
325 Double_t end = start+ndiv*step;
326 switch (iaxis) {
327 case 1: //--- divide on X
328 shape = new TGeoBBox(step/2., fDY, fDZ);
329 finder = new TGeoPatternX(voldiv, ndiv, start, end);
330 opt = "X";
331 break;
332 case 2: //--- divide on Y
333 shape = new TGeoBBox(fDX, step/2., fDZ);
334 finder = new TGeoPatternY(voldiv, ndiv, start, end);
335 opt = "Y";
336 break;
337 case 3: //--- divide on Z
338 shape = new TGeoBBox(fDX, fDY, step/2.);
339 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
340 opt = "Z";
341 break;
342 default:
343 Error("Divide", "Wrong axis type for division");
344 return 0;
345 }
346 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
347 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
348 vmulti->AddVolume(vol);
349 voldiv->SetFinder(finder);
350 finder->SetDivIndex(voldiv->GetNdaughters());
351 for (Int_t ic=0; ic<ndiv; ic++) {
352 voldiv->AddNodeOffset(vol, ic, start+step/2.+ic*step, opt.Data());
353 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
354 }
355 return vmulti;
356}
357
358////////////////////////////////////////////////////////////////////////////////
359/// Compute bounding box - nothing to do in this case.
360
362{
363}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Test if point is inside this shape.
367
369{
370 if (TMath::Abs(point[2]-fOrigin[2]) > fDZ) return kFALSE;
371 if (TMath::Abs(point[0]-fOrigin[0]) > fDX) return kFALSE;
372 if (TMath::Abs(point[1]-fOrigin[1]) > fDY) return kFALSE;
373 return kTRUE;
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// Static method to check if point[3] is located inside a box of having dx, dy, dz
378/// as half-lengths.
379
380Bool_t TGeoBBox::Contains(const Double_t *point, Double_t dx, Double_t dy, Double_t dz, const Double_t *origin)
381{
382 if (TMath::Abs(point[2]-origin[2]) > dz) return kFALSE;
383 if (TMath::Abs(point[0]-origin[0]) > dx) return kFALSE;
384 if (TMath::Abs(point[1]-origin[1]) > dy) return kFALSE;
385 return kTRUE;
386}
387
388////////////////////////////////////////////////////////////////////////////////
389/// Compute distance from inside point to surface of the box.
390/// Boundary safe algorithm.
391
392Double_t TGeoBBox::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
393{
394 Double_t s,smin,saf[6];
395 Double_t newpt[3];
396 Int_t i;
397 for (i=0; i<3; i++) newpt[i] = point[i] - fOrigin[i];
398 saf[0] = fDX+newpt[0];
399 saf[1] = fDX-newpt[0];
400 saf[2] = fDY+newpt[1];
401 saf[3] = fDY-newpt[1];
402 saf[4] = fDZ+newpt[2];
403 saf[5] = fDZ-newpt[2];
404 if (iact<3 && safe) {
405 smin = saf[0];
406 // compute safe distance
407 for (i=1;i<6;i++) if (saf[i] < smin) smin = saf[i];
408 *safe = smin;
409 if (smin<0) *safe = 0.0;
410 if (iact==0) return TGeoShape::Big();
411 if (iact==1 && step<*safe) return TGeoShape::Big();
412 }
413 // compute distance to surface
414 smin=TGeoShape::Big();
415 for (i=0; i<3; i++) {
416 if (dir[i]!=0) {
417 s = (dir[i]>0)?(saf[(i<<1)+1]/dir[i]):(-saf[i<<1]/dir[i]);
418 if (s < 0) return 0.0;
419 if (s < smin) smin = s;
420 }
421 }
422 return smin;
423}
424
425////////////////////////////////////////////////////////////////////////////////
426/// Compute distance from inside point to surface of the box.
427/// Boundary safe algorithm.
428
430 Double_t dx, Double_t dy, Double_t dz, const Double_t *origin, Double_t /*stepmax*/)
431{
432 Double_t s,smin,saf[6];
433 Double_t newpt[3];
434 Int_t i;
435 for (i=0; i<3; i++) newpt[i] = point[i] - origin[i];
436 saf[0] = dx+newpt[0];
437 saf[1] = dx-newpt[0];
438 saf[2] = dy+newpt[1];
439 saf[3] = dy-newpt[1];
440 saf[4] = dz+newpt[2];
441 saf[5] = dz-newpt[2];
442 // compute distance to surface
443 smin=TGeoShape::Big();
444 for (i=0; i<3; i++) {
445 if (dir[i]!=0) {
446 s = (dir[i]>0)?(saf[(i<<1)+1]/dir[i]):(-saf[i<<1]/dir[i]);
447 if (s < 0) return 0.0;
448 if (s < smin) smin = s;
449 }
450 }
451 return smin;
452}
453
454////////////////////////////////////////////////////////////////////////////////
455/// Compute distance from outside point to surface of the box.
456/// Boundary safe algorithm.
457
458Double_t TGeoBBox::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
459{
460 Bool_t in = kTRUE;
461 Double_t saf[3];
462 Double_t par[3];
463 Double_t newpt[3];
464 Int_t i,j;
465 for (i=0; i<3; i++) newpt[i] = point[i] - fOrigin[i];
466 par[0] = fDX;
467 par[1] = fDY;
468 par[2] = fDZ;
469 for (i=0; i<3; i++) {
470 saf[i] = TMath::Abs(newpt[i])-par[i];
471 if (saf[i]>=step) return TGeoShape::Big();
472 if (in && saf[i]>0) in=kFALSE;
473 }
474 if (iact<3 && safe) {
475 // compute safe distance
476 if (in) {
477 *safe = 0.0;
478 } else {
479 *safe = saf[0];
480 if (saf[1] > *safe) *safe = saf[1];
481 if (saf[2] > *safe) *safe = saf[2];
482 }
483 if (iact==0) return TGeoShape::Big();
484 if (iact==1 && step<*safe) return TGeoShape::Big();
485 }
486 // compute distance from point to box
487 // protection in case point is actually inside box
488 if (in) {
489 j = 0;
490 Double_t ss = saf[0];
491 if (saf[1]>ss) {
492 ss = saf[1];
493 j = 1;
494 }
495 if (saf[2]>ss) j = 2;
496 if (newpt[j]*dir[j]>0) return TGeoShape::Big(); // in fact exiting
497 return 0.0;
498 }
499 for (i=0; i<3; i++) {
500 if (saf[i]<0) continue;
501 if (newpt[i]*dir[i] >= 0) continue;
502 Double_t snxt = saf[i]/TMath::Abs(dir[i]);
503 Int_t ibreak = 0;
504 for (j=0; j<3; j++) {
505 if (j==i) continue;
506 Double_t coord=newpt[j]+snxt*dir[j];
507 if (TMath::Abs(coord)>par[j]) {
508 ibreak=1;
509 break;
510 }
511 }
512 if (!ibreak) return snxt;
513 }
514 return TGeoShape::Big();
515}
516
517////////////////////////////////////////////////////////////////////////////////
518/// Compute distance from outside point to surface of the box.
519/// Boundary safe algorithm.
520
522 Double_t dx, Double_t dy, Double_t dz, const Double_t *origin, Double_t stepmax)
523{
524 Bool_t in = kTRUE;
525 Double_t saf[3];
526 Double_t par[3];
527 Double_t newpt[3];
528 Int_t i,j;
529 for (i=0; i<3; i++) newpt[i] = point[i] - origin[i];
530 par[0] = dx;
531 par[1] = dy;
532 par[2] = dz;
533 for (i=0; i<3; i++) {
534 saf[i] = TMath::Abs(newpt[i])-par[i];
535 if (saf[i]>=stepmax) return TGeoShape::Big();
536 if (in && saf[i]>0) in=kFALSE;
537 }
538 // In case point is inside return ZERO
539 if (in) return 0.0;
540 Double_t coord, snxt=TGeoShape::Big();
541 Int_t ibreak=0;
542 for (i=0; i<3; i++) {
543 if (saf[i]<0) continue;
544 if (newpt[i]*dir[i] >= 0) continue;
545 snxt = saf[i]/TMath::Abs(dir[i]);
546 ibreak = 0;
547 for (j=0; j<3; j++) {
548 if (j==i) continue;
549 coord=newpt[j]+snxt*dir[j];
550 if (TMath::Abs(coord)>par[j]) {
551 ibreak=1;
552 break;
553 }
554 }
555 if (!ibreak) return snxt;
556 }
557 return TGeoShape::Big();
558}
559
560////////////////////////////////////////////////////////////////////////////////
561/// Returns name of axis IAXIS.
562
563const char *TGeoBBox::GetAxisName(Int_t iaxis) const
564{
565 switch (iaxis) {
566 case 1:
567 return "X";
568 case 2:
569 return "Y";
570 case 3:
571 return "Z";
572 default:
573 return "UNDEFINED";
574 }
575}
576
577////////////////////////////////////////////////////////////////////////////////
578/// Get range of shape for a given axis.
579
581{
582 xlo = 0;
583 xhi = 0;
584 Double_t dx = 0;
585 switch (iaxis) {
586 case 1:
587 xlo = fOrigin[0]-fDX;
588 xhi = fOrigin[0]+fDX;
589 dx = 2*fDX;
590 return dx;
591 case 2:
592 xlo = fOrigin[1]-fDY;
593 xhi = fOrigin[1]+fDY;
594 dx = 2*fDY;
595 return dx;
596 case 3:
597 xlo = fOrigin[2]-fDZ;
598 xhi = fOrigin[2]+fDZ;
599 dx = 2*fDZ;
600 return dx;
601 }
602 return dx;
603}
604
605////////////////////////////////////////////////////////////////////////////////
606/// Fill vector param[4] with the bounding cylinder parameters. The order
607/// is the following : Rmin, Rmax, Phi1, Phi2
608
610{
611 param[0] = 0.; // Rmin
612 param[1] = fDX*fDX+fDY*fDY; // Rmax
613 param[2] = 0.; // Phi1
614 param[3] = 360.; // Phi2
615}
616
617////////////////////////////////////////////////////////////////////////////////
618/// Get area in internal units of the facet with a given index.
619/// Possible index values:
620/// - 0 - all facets together
621/// - 1 to 6 - facet index from bottom to top Z
622
624{
625 Double_t area = 0.;
626 switch (index) {
627 case 0:
628 area = 8.*(fDX*fDY + fDX*fDZ + fDY*fDZ);
629 return area;
630 case 1:
631 case 6:
632 area = 4.*fDX*fDY;
633 return area;
634 case 2:
635 case 4:
636 area = 4.*fDX*fDZ;
637 return area;
638 case 3:
639 case 5:
640 area = 4.*fDY*fDZ;
641 return area;
642 }
643 return area;
644}
645
646////////////////////////////////////////////////////////////////////////////////
647/// Fills array with n random points located on the surface of indexed facet.
648/// The output array must be provided with a length of minimum 3*npoints. Returns
649/// true if operation succeeded.
650/// Possible index values:
651/// - 0 - all facets together
652/// - 1 to 6 - facet index from bottom to top Z
653
655{
656 if (index<0 || index>6) return kFALSE;
657 Double_t surf[6];
658 Double_t area = 0.;
659 if (index==0) {
660 for (Int_t isurf=0; isurf<6; isurf++) {
661 surf[isurf] = TGeoBBox::GetFacetArea(isurf+1);
662 if (isurf>0) surf[isurf] += surf[isurf-1];
663 }
664 area = surf[5];
665 }
666
667 for (Int_t i=0; i<npoints; i++) {
668 // Generate randomly a surface index if needed.
669 Double_t *point = &array[3*i];
670 Int_t surfindex = index;
671 if (surfindex==0) {
672 Double_t val = area*gRandom->Rndm();
673 surfindex = 2+TMath::BinarySearch(6, surf, val);
674 if (surfindex>6) surfindex=6;
675 }
676 switch (surfindex) {
677 case 1:
678 point[0] = -fDX + 2*fDX*gRandom->Rndm();
679 point[1] = -fDY + 2*fDY*gRandom->Rndm();
680 point[2] = -fDZ;
681 break;
682 case 2:
683 point[0] = -fDX + 2*fDX*gRandom->Rndm();
684 point[1] = -fDY;
685 point[2] = -fDZ + 2*fDZ*gRandom->Rndm();
686 break;
687 case 3:
688 point[0] = -fDX;
689 point[1] = -fDY + 2*fDY*gRandom->Rndm();
690 point[2] = -fDZ + 2*fDZ*gRandom->Rndm();
691 break;
692 case 4:
693 point[0] = -fDX + 2*fDX*gRandom->Rndm();
694 point[1] = fDY;
695 point[2] = -fDZ + 2*fDZ*gRandom->Rndm();
696 break;
697 case 5:
698 point[0] = fDX;
699 point[1] = -fDY + 2*fDY*gRandom->Rndm();
700 point[2] = -fDZ + 2*fDZ*gRandom->Rndm();
701 break;
702 case 6:
703 point[0] = -fDX + 2*fDX*gRandom->Rndm();
704 point[1] = -fDY + 2*fDY*gRandom->Rndm();
705 point[2] = fDZ;
706 break;
707 }
708 }
709 return kTRUE;
710}
711
712////////////////////////////////////////////////////////////////////////////////
713/// Fills array with n random points located on the line segments of the shape mesh.
714/// The output array must be provided with a length of minimum 3*npoints. Returns
715/// true if operation is implemented.
716
718{
719 if (npoints<GetNmeshVertices()) {
720 Error("GetPointsOnSegments", "You should require at least %d points", GetNmeshVertices());
721 return kFALSE;
722 }
724 Int_t npnts = buff.NbPnts();
725 Int_t nsegs = buff.NbSegs();
726 // Copy buffered points in the array
727 memcpy(array, buff.fPnts, 3*npnts*sizeof(Double_t));
728 Int_t ipoints = npoints - npnts;
729 Int_t icrt = 3*npnts;
730 Int_t nperseg = (Int_t)(Double_t(ipoints)/nsegs);
731 Double_t *p0, *p1;
732 Double_t x,y,z, dx,dy,dz;
733 for (Int_t i=0; i<nsegs; i++) {
734 p0 = &array[3*buff.fSegs[3*i+1]];
735 p1 = &array[3*buff.fSegs[3*i+2]];
736 if (i==(nsegs-1)) nperseg = ipoints;
737 dx = (p1[0]-p0[0])/(nperseg+1);
738 dy = (p1[1]-p0[1])/(nperseg+1);
739 dz = (p1[2]-p0[2])/(nperseg+1);
740 for (Int_t j=0; j<nperseg; j++) {
741 x = p0[0] + (j+1)*dx;
742 y = p0[1] + (j+1)*dy;
743 z = p0[2] + (j+1)*dz;
744 array[icrt++] = x; array[icrt++] = y; array[icrt++] = z;
745 ipoints--;
746 }
747 }
748 return kTRUE;
749}
750
751////////////////////////////////////////////////////////////////////////////////
752/// Fills real parameters of a positioned box inside this one. Returns 0 if successful.
753
754Int_t TGeoBBox::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
755{
756 dx=dy=dz=0;
757 if (mat->IsRotation()) {
758 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
759 return 1; // ### rotation not accepted ###
760 }
761 //--> translate the origin of the parametrized box to the frame of this box.
762 Double_t origin[3];
763 mat->LocalToMaster(parambox->GetOrigin(), origin);
764 if (!Contains(origin)) {
765 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
766 return 1; // ### wrong matrix ###
767 }
768 //--> now we have to get the valid range for all parametrized axis
769 Double_t xlo=0, xhi=0;
770 Double_t dd[3];
771 dd[0] = parambox->GetDX();
772 dd[1] = parambox->GetDY();
773 dd[2] = parambox->GetDZ();
774 for (Int_t iaxis=0; iaxis<3; iaxis++) {
775 if (dd[iaxis]>=0) continue;
776 TGeoBBox::GetAxisRange(iaxis+1, xlo, xhi);
777 //-> compute best fitting parameter
778 dd[iaxis] = TMath::Min(origin[iaxis]-xlo, xhi-origin[iaxis]);
779 if (dd[iaxis]<0) {
780 Error("GetFittingBox", "wrong matrix");
781 return 1;
782 }
783 }
784 dx = dd[0];
785 dy = dd[1];
786 dz = dd[2];
787 return 0;
788}
789
790////////////////////////////////////////////////////////////////////////////////
791/// In case shape has some negative parameters, these has to be computed
792/// in order to fit the mother
793
795{
796 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
797 Double_t dx, dy, dz;
798 Int_t ierr = mother->GetFittingBox(this, mat, dx, dy, dz);
799 if (ierr) {
800 Error("GetMakeRuntimeShape", "cannot fit this to mother");
801 return 0;
802 }
803 return (new TGeoBBox(dx, dy, dz));
804}
805
806////////////////////////////////////////////////////////////////////////////////
807/// Returns numbers of vertices, segments and polygons composing the shape mesh.
808
809void TGeoBBox::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
810{
811 nvert = 8;
812 nsegs = 12;
813 npols = 6;
814}
815
816////////////////////////////////////////////////////////////////////////////////
817/// Prints shape parameters
818
820{
821 printf("*** Shape %s: TGeoBBox ***\n", GetName());
822 printf(" dX = %11.5f\n", fDX);
823 printf(" dY = %11.5f\n", fDY);
824 printf(" dZ = %11.5f\n", fDZ);
825 printf(" origin: x=%11.5f y=%11.5f z=%11.5f\n", fOrigin[0], fOrigin[1], fOrigin[2]);
826}
827
828////////////////////////////////////////////////////////////////////////////////
829/// Creates a TBuffer3D describing *this* shape.
830/// Coordinates are in local reference frame.
831
833{
834 TBuffer3D* buff = new TBuffer3D(TBuffer3DTypes::kGeneric, 8, 24, 12, 36, 6, 36);
835 if (buff)
836 {
837 SetPoints(buff->fPnts);
838 SetSegsAndPols(*buff);
839 }
840
841 return buff;
842}
843
844////////////////////////////////////////////////////////////////////////////////
845/// Fills TBuffer3D structure for segments and polygons.
846
848{
850
851 buff.fSegs[ 0] = c ; buff.fSegs[ 1] = 0 ; buff.fSegs[ 2] = 1 ;
852 buff.fSegs[ 3] = c+1 ; buff.fSegs[ 4] = 1 ; buff.fSegs[ 5] = 2 ;
853 buff.fSegs[ 6] = c+1 ; buff.fSegs[ 7] = 2 ; buff.fSegs[ 8] = 3 ;
854 buff.fSegs[ 9] = c ; buff.fSegs[10] = 3 ; buff.fSegs[11] = 0 ;
855 buff.fSegs[12] = c+2 ; buff.fSegs[13] = 4 ; buff.fSegs[14] = 5 ;
856 buff.fSegs[15] = c+2 ; buff.fSegs[16] = 5 ; buff.fSegs[17] = 6 ;
857 buff.fSegs[18] = c+3 ; buff.fSegs[19] = 6 ; buff.fSegs[20] = 7 ;
858 buff.fSegs[21] = c+3 ; buff.fSegs[22] = 7 ; buff.fSegs[23] = 4 ;
859 buff.fSegs[24] = c ; buff.fSegs[25] = 0 ; buff.fSegs[26] = 4 ;
860 buff.fSegs[27] = c+2 ; buff.fSegs[28] = 1 ; buff.fSegs[29] = 5 ;
861 buff.fSegs[30] = c+1 ; buff.fSegs[31] = 2 ; buff.fSegs[32] = 6 ;
862 buff.fSegs[33] = c+3 ; buff.fSegs[34] = 3 ; buff.fSegs[35] = 7 ;
863
864 buff.fPols[ 0] = c ; buff.fPols[ 1] = 4 ; buff.fPols[ 2] = 0 ;
865 buff.fPols[ 3] = 9 ; buff.fPols[ 4] = 4 ; buff.fPols[ 5] = 8 ;
866 buff.fPols[ 6] = c+1 ; buff.fPols[ 7] = 4 ; buff.fPols[ 8] = 1 ;
867 buff.fPols[ 9] = 10 ; buff.fPols[10] = 5 ; buff.fPols[11] = 9 ;
868 buff.fPols[12] = c ; buff.fPols[13] = 4 ; buff.fPols[14] = 2 ;
869 buff.fPols[15] = 11 ; buff.fPols[16] = 6 ; buff.fPols[17] = 10 ;
870 buff.fPols[18] = c+1 ; buff.fPols[19] = 4 ; buff.fPols[20] = 3 ;
871 buff.fPols[21] = 8 ; buff.fPols[22] = 7 ; buff.fPols[23] = 11 ;
872 buff.fPols[24] = c+2 ; buff.fPols[25] = 4 ; buff.fPols[26] = 0 ;
873 buff.fPols[27] = 3 ; buff.fPols[28] = 2 ; buff.fPols[29] = 1 ;
874 buff.fPols[30] = c+3 ; buff.fPols[31] = 4 ; buff.fPols[32] = 4 ;
875 buff.fPols[33] = 5 ; buff.fPols[34] = 6 ; buff.fPols[35] = 7 ;
876}
877
878////////////////////////////////////////////////////////////////////////////////
879/// Computes the closest distance from given point to this shape.
880
882{
883 Double_t safe, safy, safz;
884 if (in) {
885 safe = fDX - TMath::Abs(point[0]-fOrigin[0]);
886 safy = fDY - TMath::Abs(point[1]-fOrigin[1]);
887 safz = fDZ - TMath::Abs(point[2]-fOrigin[2]);
888 if (safy < safe) safe = safy;
889 if (safz < safe) safe = safz;
890 } else {
891 safe = -fDX + TMath::Abs(point[0]-fOrigin[0]);
892 safy = -fDY + TMath::Abs(point[1]-fOrigin[1]);
893 safz = -fDZ + TMath::Abs(point[2]-fOrigin[2]);
894 if (safy > safe) safe = safy;
895 if (safz > safe) safe = safz;
896 }
897 return safe;
898}
899
900////////////////////////////////////////////////////////////////////////////////
901/// Save a primitive as a C++ statement(s) on output stream "out".
902
903void TGeoBBox::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
904{
906 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
907 out << " dx = " << fDX << ";" << std::endl;
908 out << " dy = " << fDY << ";" << std::endl;
909 out << " dz = " << fDZ << ";" << std::endl;
913 out << " origin[0] = " << fOrigin[0] << ";" << std::endl;
914 out << " origin[1] = " << fOrigin[1] << ";" << std::endl;
915 out << " origin[2] = " << fOrigin[2] << ";" << std::endl;
916 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz,origin);" << std::endl;
917 } else {
918 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz);" << std::endl;
919 }
921}
922
923////////////////////////////////////////////////////////////////////////////////
924/// Set parameters of the box.
925
927{
928 fDX = dx;
929 fDY = dy;
930 fDZ = dz;
931 if (origin) {
932 fOrigin[0] = origin[0];
933 fOrigin[1] = origin[1];
934 fOrigin[2] = origin[2];
935 }
939 if ((fDX<0) || (fDY<0) || (fDZ<0)) SetShapeBit(kGeoRunTimeShape);
940}
941
942////////////////////////////////////////////////////////////////////////////////
943/// Set dimensions based on the array of parameters
944/// param[0] - half-length in x
945/// param[1] - half-length in y
946/// param[2] - half-length in z
947
949{
950 if (!param) {
951 Error("SetDimensions", "null parameters");
952 return;
953 }
954 fDX = param[0];
955 fDY = param[1];
956 fDZ = param[2];
960 if ((fDX<0) || (fDY<0) || (fDZ<0)) SetShapeBit(kGeoRunTimeShape);
961}
962
963////////////////////////////////////////////////////////////////////////////////
964/// Fill box vertices to an array.
965
967{
969}
970
971////////////////////////////////////////////////////////////////////////////////
972/// Fill box points.
973
975{
976 if (!points) return;
977 Double_t xmin,xmax,ymin,ymax,zmin,zmax;
978 xmin = -fDX+fOrigin[0];
979 xmax = fDX+fOrigin[0];
980 ymin = -fDY+fOrigin[1];
981 ymax = fDY+fOrigin[1];
982 zmin = -fDZ+fOrigin[2];
983 zmax = fDZ+fOrigin[2];
984 points[ 0] = xmin; points[ 1] = ymin; points[ 2] = zmin;
985 points[ 3] = xmin; points[ 4] = ymax; points[ 5] = zmin;
986 points[ 6] = xmax; points[ 7] = ymax; points[ 8] = zmin;
987 points[ 9] = xmax; points[10] = ymin; points[11] = zmin;
988 points[12] = xmin; points[13] = ymin; points[14] = zmax;
989 points[15] = xmin; points[16] = ymax; points[17] = zmax;
990 points[18] = xmax; points[19] = ymax; points[20] = zmax;
991 points[21] = xmax; points[22] = ymin; points[23] = zmax;
992}
993
994////////////////////////////////////////////////////////////////////////////////
995/// Fill box points.
996
998{
999 if (!points) return;
1000 Double_t xmin,xmax,ymin,ymax,zmin,zmax;
1001 xmin = -fDX+fOrigin[0];
1002 xmax = fDX+fOrigin[0];
1003 ymin = -fDY+fOrigin[1];
1004 ymax = fDY+fOrigin[1];
1005 zmin = -fDZ+fOrigin[2];
1006 zmax = fDZ+fOrigin[2];
1007 points[ 0] = xmin; points[ 1] = ymin; points[ 2] = zmin;
1008 points[ 3] = xmin; points[ 4] = ymax; points[ 5] = zmin;
1009 points[ 6] = xmax; points[ 7] = ymax; points[ 8] = zmin;
1010 points[ 9] = xmax; points[10] = ymin; points[11] = zmin;
1011 points[12] = xmin; points[13] = ymin; points[14] = zmax;
1012 points[15] = xmin; points[16] = ymax; points[17] = zmax;
1013 points[18] = xmax; points[19] = ymax; points[20] = zmax;
1014 points[21] = xmax; points[22] = ymin; points[23] = zmax;
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018////// fill size of this 3-D object
1019//// TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
1020//// if (painter) painter->AddSize3D(8, 12, 6);
1021
1023{
1024}
1025
1026////////////////////////////////////////////////////////////////////////////////
1027/// Fills a static 3D buffer and returns a reference.
1028
1029const TBuffer3D & TGeoBBox::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1030{
1031 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1032
1033 FillBuffer3D(buffer, reqSections, localFrame);
1034
1035 // TODO: A box itself has has nothing more as already described
1036 // by bounding box. How will viewer interpret?
1037 if (reqSections & TBuffer3D::kRawSizes) {
1038 if (buffer.SetRawSizes(8, 3*8, 12, 3*12, 6, 6*6)) {
1040 }
1041 }
1042 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1043 SetPoints(buffer.fPnts);
1044 if (!buffer.fLocalFrame) {
1045 TransformPoints(buffer.fPnts, buffer.NbPnts());
1046 }
1047
1048 SetSegsAndPols(buffer);
1050 }
1051
1052 return buffer;
1053}
1054
1055////////////////////////////////////////////////////////////////////////////////
1056/// Fills the supplied buffer, with sections in desired frame
1057/// See TBuffer3D.h for explanation of sections, frame etc.
1058
1059void TGeoBBox::FillBuffer3D(TBuffer3D & buffer, Int_t reqSections, Bool_t localFrame) const
1060{
1061 TGeoShape::FillBuffer3D(buffer, reqSections, localFrame);
1062
1063 if (reqSections & TBuffer3D::kBoundingBox) {
1064 Double_t halfLengths[3] = { fDX, fDY, fDZ };
1065 buffer.SetAABoundingBox(fOrigin, halfLengths);
1066
1067 if (!buffer.fLocalFrame) {
1068 TransformPoints(buffer.fBBVertex[0], 8);
1069 }
1071 }
1072}
1073
1074////////////////////////////////////////////////////////////////////////////////
1075/// Check the inside status for each of the points in the array.
1076/// Input: Array of point coordinates + vector size
1077/// Output: Array of Booleans for the inside of each point
1078
1079void TGeoBBox::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1080{
1081 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
1082}
1083
1084////////////////////////////////////////////////////////////////////////////////
1085/// Compute the normal for an array o points so that norm.dot.dir is positive
1086/// Input: Arrays of point coordinates and directions + vector size
1087/// Output: Array of normal directions
1088
1089void TGeoBBox::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1090{
1091 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
1092}
1093
1094////////////////////////////////////////////////////////////////////////////////
1095/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1096
1097void TGeoBBox::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1098{
1099 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
1100}
1101
1102////////////////////////////////////////////////////////////////////////////////
1103/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1104
1105void TGeoBBox::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
1106{
1107 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
1108}
1109
1110////////////////////////////////////////////////////////////////////////////////
1111/// Compute safe distance from each of the points in the input array.
1112/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1113/// Output: Safety values
1114
1115void TGeoBBox::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1116{
1117 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
1118}
#define c(i)
Definition RSha256.hxx:101
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
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
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
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:114
UInt_t NbPnts() const
Definition TBuffer3D.h:80
UInt_t NbSegs() const
Definition TBuffer3D.h:81
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
Int_t * fSegs
Definition TBuffer3D.h:113
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:112
Double_t fBBVertex[8][3]
Definition TBuffer3D.h:107
Box class.
Definition TGeoBBox.h:18
virtual const Double_t * GetOrigin() const
Definition TGeoBBox.h:77
virtual void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
Returns numbers of vertices, segments and polygons composing the shape mesh.
Definition TGeoBBox.cxx:809
Double_t fDX
Definition TGeoBBox.h:21
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:654
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
Computes the closest distance from given point to this shape.
Definition TGeoBBox.cxx:881
void SetBoxDimensions(Double_t dx, Double_t dy, Double_t dz, Double_t *origin=nullptr)
Set parameters of the box.
Definition TGeoBBox.cxx:926
virtual Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
Fills real parameters of a positioned box inside this one. Returns 0 if successful.
Definition TGeoBBox.cxx:754
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this box shape belonging to volume "voldiv" into ndiv equal volumes called divname,...
Definition TGeoBBox.cxx:317
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
virtual Double_t GetDX() const
Definition TGeoBBox.h:74
virtual Double_t GetFacetArea(Int_t index=0) const
Get area in internal units of the facet with a given index.
Definition TGeoBBox.cxx:623
void SetBoxPoints(Double_t *points) const
Fill box vertices to an array.
Definition TGeoBBox.cxx:966
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition TGeoBBox.cxx:580
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
Compute distance from outside point to surface of the box.
Definition TGeoBBox.cxx:458
virtual Double_t GetDZ() const
Definition TGeoBBox.h:76
virtual void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
Check the inside status for each of the points in the array.
virtual void Sizeof3D() const
virtual TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const
In case shape has some negative parameters, these has to be computed in order to fit the mother.
Definition TGeoBBox.cxx:794
virtual Int_t GetNmeshVertices() const
Definition TGeoBBox.h:73
virtual Double_t GetDY() const
Definition TGeoBBox.h:75
virtual void ComputeBBox()
Compute bounding box - nothing to do in this case.
Definition TGeoBBox.cxx:361
virtual void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
Compute safe distance from each of the points in the input array.
static Bool_t AreOverlapping(const TGeoBBox *box1, const TGeoMatrix *mat1, const TGeoBBox *box2, const TGeoMatrix *mat2)
Check if 2 positioned boxes overlap.
Definition TGeoBBox.cxx:217
Double_t fOrigin[3]
Definition TGeoBBox.h:24
virtual void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual ~TGeoBBox()
Destructor.
Definition TGeoBBox.cxx:210
virtual void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
Compute the normal for an array o points so that norm.dot.dir is positive Input: Arrays of point coor...
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const
Compute distance from inside point to surface of the box.
Definition TGeoBBox.cxx:392
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoBBox.cxx:258
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute closest distance from point px,py to each corner.
Definition TGeoBBox.cxx:305
virtual void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const
Compute distance from array of input points having directions specified by dirs. Store output in dist...
virtual Bool_t CouldBeCrossed(const Double_t *point, const Double_t *dir) const
Decides fast if the bounding box could be crossed by a vector.
Definition TGeoBBox.cxx:281
Double_t fDY
Definition TGeoBBox.h:22
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Fills a static 3D buffer and returns a reference.
virtual Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const
Fills array with n random points located on the line segments of the shape mesh.
Definition TGeoBBox.cxx:717
virtual const char * GetAxisName(Int_t iaxis) const
Returns name of axis IAXIS.
Definition TGeoBBox.cxx:563
virtual void SetSegsAndPols(TBuffer3D &buffer) const
Fills TBuffer3D structure for segments and polygons.
Definition TGeoBBox.cxx:847
virtual Bool_t Contains(const Double_t *point) const
Test if point is inside this shape.
Definition TGeoBBox.cxx:368
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoBBox.cxx:903
Double_t fDZ
Definition TGeoBBox.h:23
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
Fills the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections...
virtual TBuffer3D * MakeBuffer3D() const
Creates a TBuffer3D describing this shape.
Definition TGeoBBox.cxx:832
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Computes normal to closest surface from POINT.
Definition TGeoBBox.cxx:266
TGeoBBox()
Default constructor.
Definition TGeoBBox.cxx:164
virtual void SetDimensions(Double_t *param)
Set dimensions based on the array of parameters param[0] - half-length in x param[1] - half-length in...
Definition TGeoBBox.cxx:948
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoBBox.cxx:609
virtual void SetPoints(Double_t *points) const
Fill box points.
Definition TGeoBBox.cxx:974
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Geometrical transformation package.
Definition TGeoMatrix.h:41
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void MasterToLocalVect(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
Bool_t IsRotation() const
Definition TGeoMatrix.h:68
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
Node containing an offset.
Definition TGeoNode.h:184
Base finder class for patterns.
void SetDivIndex(Int_t index)
Base abstract class for all shapes.
Definition TGeoShape.h:26
static Double_t Big()
Definition TGeoShape.h:88
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.
virtual Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const =0
virtual const char * GetName() const
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoRunTimeShape
Definition TGeoShape.h:41
static Double_t Tolerance()
Definition TGeoShape.h:91
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
Volume families.
Definition TGeoVolume.h:255
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:49
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:174
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:232
Int_t GetNdaughters() const
Definition TGeoVolume.h:351
TObjArray * GetNodes()
Definition TGeoVolume.h:168
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:552
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:980
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
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:347
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123