Logo ROOT   6.18/05
Reference Guide
TGeoTrd2.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 31/01/02
3// TGeoTrd2::Contains() and DistFromInside() 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 TGeoTrd2
14\ingroup Geometry_classes
15A trapezoid with both x and y lengths varying with z. It
16has 5 parameters, the half lengths in x at -dz and +dz, the half
17lengths in y at -dz and +dz, and the half length in z (dz).
18
19Begin_Macro(source)
20{
21 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
22 new TGeoManager("trd2", "poza9");
23 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
24 TGeoMedium *med = new TGeoMedium("MED",1,mat);
25 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
26 gGeoManager->SetTopVolume(top);
27 TGeoVolume *vol = gGeoManager->MakeTrd2("Trd2",med, 10,20,30,10,40);
28 vol->SetLineWidth(2);
29 gGeoManager->CloseGeometry();
30 gGeoManager->SetNsegments(80);
31 top->Draw();
32 TView *view = gPad->GetView();
33 view->ShowAxis();
34}
35End_Macro
36*/
37
38#include "Riostream.h"
39
40#include "TGeoManager.h"
41#include "TGeoMatrix.h"
42#include "TGeoVolume.h"
43#include "TGeoTrd2.h"
44#include "TMath.h"
45
47
48////////////////////////////////////////////////////////////////////////////////
49/// dummy ctor
50
52{
54 fDz = fDx1 = fDx2 = fDy1 = fDy2 = 0;
55}
56
57////////////////////////////////////////////////////////////////////////////////
58/// constructor.
59
61 :TGeoBBox(0,0,0)
62{
64 fDx1 = dx1;
65 fDx2 = dx2;
66 fDy1 = dy1;
67 fDy2 = dy2;
68 fDz = dz;
69 if ((fDx1<0) || (fDx2<0) || (fDy1<0) || (fDy2<0) || (fDz<0)) {
71 printf("trd2 : dx1=%f, dx2=%f, dy1=%f, dy2=%f, dz=%f\n",
72 dx1,dx2,dy1,dy2,dz);
73 }
74 else ComputeBBox();
75}
76
77////////////////////////////////////////////////////////////////////////////////
78/// constructor.
79
80TGeoTrd2::TGeoTrd2(const char * name, Double_t dx1, Double_t dx2, Double_t dy1, Double_t dy2, Double_t dz)
81 :TGeoBBox(name, 0,0,0)
82{
84 fDx1 = dx1;
85 fDx2 = dx2;
86 fDy1 = dy1;
87 fDy2 = dy2;
88 fDz = dz;
89 if ((fDx1<0) || (fDx2<0) || (fDy1<0) || (fDy2<0) || (fDz<0)) {
91 printf("trd2 : dx1=%f, dx2=%f, dy1=%f, dy2=%f, dz=%f\n",
92 dx1,dx2,dy1,dy2,dz);
93 }
94 else ComputeBBox();
95}
96
97////////////////////////////////////////////////////////////////////////////////
98/// ctor with an array of parameters
99/// - param[0] = dx1
100/// - param[1] = dx2
101/// - param[2] = dy1
102/// - param[3] = dy2
103/// - param[4] = dz
104
106 :TGeoBBox(0,0,0)
107{
109 SetDimensions(param);
110 if ((fDx1<0) || (fDx2<0) || (fDy1<0) || (fDy2<0) || (fDz<0)) SetShapeBit(kGeoRunTimeShape);
111 else ComputeBBox();
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// destructor
116
118{
119}
120
121////////////////////////////////////////////////////////////////////////////////
122/// Computes capacity of the shape in [length^3]
123
125{
126 Double_t capacity = 2*(fDx1+fDx2)*(fDy1+fDy2)*fDz +
127 (2./3.)*(fDx1-fDx2)*(fDy1-fDy2)*fDz;
128 return capacity;
129}
130
131////////////////////////////////////////////////////////////////////////////////
132/// compute bounding box for a trd2
133
135{
138 fDZ = fDz;
139 memset(fOrigin, 0, 3*sizeof(Double_t));
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// Compute normal to closest surface from POINT.
144
145void TGeoTrd2::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
146{
147 Double_t safe, safemin;
148 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
149 Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
150 // check Z facettes
151 safe = safemin = TMath::Abs(fDz-TMath::Abs(point[2]));
152 norm[0] = norm[1] = 0;
153 norm[2] = (dir[2]>=0)?1:-1;
154 if (safe<TGeoShape::Tolerance()) return;
155 // check X facettes
156 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
157 if (distx>=0) {
158 safe=TMath::Abs(distx-TMath::Abs(point[0]))*calf;
159 if (safe<safemin) {
160 safemin = safe;
161 norm[0] = (point[0]>0)?calf:(-calf);
162 norm[1] = 0;
163 norm[2] = calf*fx;
164 Double_t dot = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
165 if (dot<0) {
166 norm[0]=-norm[0];
167 norm[2]=-norm[2];
168 }
169 if (safe<TGeoShape::Tolerance()) return;
170 }
171 }
172
173 Double_t fy = 0.5*(fDy1-fDy2)/fDz;
174 calf = 1./TMath::Sqrt(1.0+fy*fy);
175
176 // check Y facettes
177 distx = 0.5*(fDy1+fDy2)-fy*point[2];
178 if (distx>=0) {
179 safe=TMath::Abs(distx-TMath::Abs(point[1]))*calf;
180 if (safe<safemin) {
181 norm[0] = 0;
182 norm[1] = (point[1]>0)?calf:(-calf);
183 norm[2] = calf*fy;
184 Double_t dot = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
185 if (dot<0) {
186 norm[1]=-norm[1];
187 norm[2]=-norm[2];
188 }
189 }
190 }
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// test if point is inside this shape
195/// check Z range
196
198{
199 if (TMath::Abs(point[2]) > fDz) return kFALSE;
200 // then y
201 Double_t dy = 0.5*(fDy2*(point[2]+fDz)+fDy1*(fDz-point[2]))/fDz;
202 if (TMath::Abs(point[1]) > dy) return kFALSE;
203 // then x
204 Double_t dx = 0.5*(fDx2*(point[2]+fDz)+fDx1*(fDz-point[2]))/fDz;
205 if (TMath::Abs(point[0]) > dx) return kFALSE;
206 return kTRUE;
207}
208
209////////////////////////////////////////////////////////////////////////////////
210/// Compute distance from inside point to surface of the trd2
211/// Boundary safe algorithm
212
213Double_t TGeoTrd2::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
214{
215 Double_t snxt = TGeoShape::Big();
216 if (iact<3 && safe) {
217 // compute safe distance
218 *safe = Safety(point, kTRUE);
219 if (iact==0) return TGeoShape::Big();
220 if (iact==1 && step<*safe) return TGeoShape::Big();
221 }
222
223 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
224 Double_t fy = 0.5*(fDy1-fDy2)/fDz;
225 Double_t cn;
226
227 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
228 Double_t disty = 0.5*(fDy1+fDy2)-fy*point[2];
229
230 //--- Compute distance to this shape
231 // first check if Z facettes are crossed
232 Double_t dist[3];
233 for (Int_t i=0; i<3; i++) dist[i]=TGeoShape::Big();
234 if (dir[2]<0) {
235 dist[0]=-(point[2]+fDz)/dir[2];
236 } else if (dir[2]>0) {
237 dist[0]=(fDz-point[2])/dir[2];
238 }
239 if (dist[0]<=0) return 0.0;
240 // now check X facettes
241 cn = -dir[0]+fx*dir[2];
242 if (cn>0) {
243 dist[1] = point[0]+distx;
244 if (dist[1]<=0) return 0.0;
245 dist[1] /= cn;
246 }
247 cn = dir[0]+fx*dir[2];
248 if (cn>0) {
249 Double_t s = distx-point[0];
250 if (s<=0) return 0.0;
251 s /= cn;
252 if (s<dist[1]) dist[1] = s;
253 }
254 // now check Y facettes
255 cn = -dir[1]+fy*dir[2];
256 if (cn>0) {
257 dist[2] = point[1]+disty;
258 if (dist[2]<=0) return 0.0;
259 dist[2] /= cn;
260 }
261 cn = dir[1]+fy*dir[2];
262 if (cn>0) {
263 Double_t s = disty-point[1];
264 if (s<=0) return 0.0;
265 s /= cn;
266 if (s<dist[2]) dist[2] = s;
267 }
268 snxt = dist[TMath::LocMin(3,dist)];
269 return snxt;
270}
271
272////////////////////////////////////////////////////////////////////////////////
273/// Compute distance from outside point to surface of the trd2
274/// Boundary safe algorithm
275
276Double_t TGeoTrd2::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
277{
278 Double_t snxt = TGeoShape::Big();
279 if (iact<3 && safe) {
280 // compute safe distance
281 *safe = Safety(point, kFALSE);
282 if (iact==0) return TGeoShape::Big();
283 if (iact==1 && step<*safe) return TGeoShape::Big();
284 }
285 // find a visible face
286 Double_t xnew,ynew,znew;
287 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
288 Double_t fy = 0.5*(fDy1-fDy2)/fDz;
289 Double_t cn;
290 // check visibility of X faces
291 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
292 Double_t disty = 0.5*(fDy1+fDy2)-fy*point[2];
293 Bool_t in = kTRUE;
294 Double_t safx = distx-TMath::Abs(point[0]);
295 Double_t safy = disty-TMath::Abs(point[1]);
296 Double_t safz = fDz-TMath::Abs(point[2]);
297 //--- Compute distance to this shape
298 // first check if Z facettes are crossed
299 if (point[2]<=-fDz) {
300 cn = -dir[2];
301 if (cn>=0) return TGeoShape::Big();
302 in = kFALSE;
303 snxt = (fDz+point[2])/cn;
304 // find extrapolated X and Y
305 xnew = point[0]+snxt*dir[0];
306 if (TMath::Abs(xnew) < fDx1) {
307 ynew = point[1]+snxt*dir[1];
308 if (TMath::Abs(ynew) < fDy1) return snxt;
309 }
310 } else if (point[2]>=fDz) {
311 cn = dir[2];
312 if (cn>=0) return TGeoShape::Big();
313 in = kFALSE;
314 snxt = (fDz-point[2])/cn;
315 // find extrapolated X and Y
316 xnew = point[0]+snxt*dir[0];
317 if (TMath::Abs(xnew) < fDx2) {
318 ynew = point[1]+snxt*dir[1];
319 if (TMath::Abs(ynew) < fDy2) return snxt;
320 }
321 }
322 // check if X facettes are crossed
323 if (point[0]<=-distx) {
324 cn = -dir[0]+fx*dir[2];
325 if (cn>=0) return TGeoShape::Big();
326 in = kFALSE;
327 snxt = (point[0]+distx)/cn;
328 // find extrapolated Y and Z
329 znew = point[2]+snxt*dir[2];
330 if (TMath::Abs(znew) < fDz) {
331 Double_t dy = 0.5*(fDy1+fDy2)-fy*znew;
332 ynew = point[1]+snxt*dir[1];
333 if (TMath::Abs(ynew) < dy) return snxt;
334 }
335 }
336 if (point[0]>=distx) {
337 cn = dir[0]+fx*dir[2];
338 if (cn>=0) return TGeoShape::Big();
339 in = kFALSE;
340 snxt = (distx-point[0])/cn;
341 // find extrapolated Y and Z
342 znew = point[2]+snxt*dir[2];
343 if (TMath::Abs(znew) < fDz) {
344 Double_t dy = 0.5*(fDy1+fDy2)-fy*znew;
345 ynew = point[1]+snxt*dir[1];
346 if (TMath::Abs(ynew) < dy) return snxt;
347 }
348 }
349 // finally check Y facettes
350 if (point[1]<=-disty) {
351 cn = -dir[1]+fy*dir[2];
352 in = kFALSE;
353 if (cn>=0) return TGeoShape::Big();
354 snxt = (point[1]+disty)/cn;
355 // find extrapolated X and Z
356 znew = point[2]+snxt*dir[2];
357 if (TMath::Abs(znew) < fDz) {
358 Double_t dx = 0.5*(fDx1+fDx2)-fx*znew;
359 xnew = point[0]+snxt*dir[0];
360 if (TMath::Abs(xnew) < dx) return snxt;
361 }
362 }
363 if (point[1]>=disty) {
364 cn = dir[1]+fy*dir[2];
365 if (cn>=0) return TGeoShape::Big();
366 in = kFALSE;
367 snxt = (disty-point[1])/cn;
368 // find extrapolated X and Z
369 znew = point[2]+snxt*dir[2];
370 if (TMath::Abs(znew) < fDz) {
371 Double_t dx = 0.5*(fDx1+fDx2)-fx*znew;
372 xnew = point[0]+snxt*dir[0];
373 if (TMath::Abs(xnew) < dx) return snxt;
374 }
375 }
376 if (!in) return TGeoShape::Big();
377 // Point actually inside
378 if (safz<safx && safz<safy) {
379 if (point[2]*dir[2]>=0) return TGeoShape::Big();
380 return 0.0;
381 }
382 if (safy<safx) {
383 cn = TMath::Sign(1.0,point[1])*dir[1]+fy*dir[2];
384 if (cn>=0) return TGeoShape::Big();
385 return 0.0;
386 }
387 cn = TMath::Sign(1.0,point[0])*dir[0]+fx*dir[2];
388 if (cn>=0) return TGeoShape::Big();
389 return 0.0;
390}
391
392////////////////////////////////////////////////////////////////////////////////
393/// Get range of shape for a given axis.
394
396{
397 xlo = 0;
398 xhi = 0;
399 Double_t dx = 0;
400 switch (iaxis) {
401 case 3:
402 xlo = -fDz;
403 xhi = fDz;
404 dx = xhi-xlo;
405 return dx;
406 }
407 return dx;
408}
409
410////////////////////////////////////////////////////////////////////////////////
411/// get the most visible corner from outside point and the normals
412
413void TGeoTrd2::GetVisibleCorner(const Double_t *point, Double_t *vertex, Double_t *normals) const
414{
415 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
416 Double_t fy = 0.5*(fDy1-fDy2)/fDz;
417 Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
418 Double_t salf = calf*fx;
419 Double_t cbet = 1./TMath::Sqrt(1.0+fy*fy);
420 Double_t sbet = cbet*fy;
421 // check visibility of X,Y faces
422 Double_t distx = fDx1-fx*(fDz+point[2]);
423 Double_t disty = fDy1-fy*(fDz+point[2]);
424 memset(normals, 0, 9*sizeof(Double_t));
425 TGeoTrd2 *trd2 = (TGeoTrd2*)this;
426 if (point[0]>distx) {
427 // hi x face visible
428 trd2->SetShapeBit(kGeoVisX);
429 normals[0]=calf;
430 normals[2]=salf;
431 } else {
433 normals[0]=-calf;
434 normals[2]=salf;
435 }
436 if (point[1]>disty) {
437 // hi y face visible
438 trd2->SetShapeBit(kGeoVisY);
439 normals[4]=cbet;
440 normals[5]=sbet;
441 } else {
443 normals[4]=-cbet;
444 normals[5]=sbet;
445 }
446 if (point[2]>fDz) {
447 // hi z face visible
448 trd2->SetShapeBit(kGeoVisZ);
449 normals[8]=1;
450 } else {
452 normals[8]=-1;
453 }
455}
456
457////////////////////////////////////////////////////////////////////////////////
458/// get the opposite corner of the intersected face
459
460void TGeoTrd2::GetOppositeCorner(const Double_t * /*point*/, Int_t inorm, Double_t *vertex, Double_t *normals) const
461{
462 TGeoTrd2 *trd2 = (TGeoTrd2*)this;
463 if (inorm != 0) {
464 // change x face
466 normals[0]=-normals[0];
467 }
468 if (inorm != 1) {
469 // change y face
471 normals[4]=-normals[4];
472 }
473 if (inorm != 2) {
474 // hi z face visible
476 normals[8]=-normals[8];
477 }
479}
480
481////////////////////////////////////////////////////////////////////////////////
482/// Divide this trd2 shape belonging to volume "voldiv" into ndiv volumes
483/// called divname, from start position with the given step. Only Z divisions
484/// are supported. For Z divisions just return the pointer to the volume to be
485/// divided. In case a wrong division axis is supplied, returns pointer to
486/// volume that was divided.
487
488TGeoVolume *TGeoTrd2::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
489 Double_t start, Double_t step)
490{
491 TGeoShape *shape; //--- shape to be created
492 TGeoVolume *vol; //--- division volume to be created
493 TGeoVolumeMulti *vmulti; //--- generic divided volume
494 TGeoPatternFinder *finder; //--- finder to be attached
495 TString opt = ""; //--- option to be attached
496 Double_t zmin, zmax, dx1n, dx2n, dy1n, dy2n;
497 Int_t id;
498 Double_t end = start+ndiv*step;
499 switch (iaxis) {
500 case 1:
501 Warning("Divide", "dividing a Trd2 on X not implemented");
502 return 0;
503 case 2:
504 Warning("Divide", "dividing a Trd2 on Y not implemented");
505 return 0;
506 case 3:
507 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
508 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
509 voldiv->SetFinder(finder);
510 finder->SetDivIndex(voldiv->GetNdaughters());
511 for (id=0; id<ndiv; id++) {
512 zmin = start+id*step;
513 zmax = start+(id+1)*step;
514 dx1n = 0.5*(fDx1*(fDz-zmin)+fDx2*(fDz+zmin))/fDz;
515 dx2n = 0.5*(fDx1*(fDz-zmax)+fDx2*(fDz+zmax))/fDz;
516 dy1n = 0.5*(fDy1*(fDz-zmin)+fDy2*(fDz+zmin))/fDz;
517 dy2n = 0.5*(fDy1*(fDz-zmax)+fDy2*(fDz+zmax))/fDz;
518 shape = new TGeoTrd2(dx1n, dx2n, dy1n, dy2n, step/2.);
519 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
520 vmulti->AddVolume(vol);
521 opt = "Z";
522 voldiv->AddNodeOffset(vol, id, start+step/2+id*step, opt.Data());
523 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
524 }
525 return vmulti;
526 default:
527 Error("Divide", "Wrong axis type for division");
528 return 0;
529 }
530}
531
532////////////////////////////////////////////////////////////////////////////////
533/// Fill vector param[4] with the bounding cylinder parameters. The order
534/// is the following : Rmin, Rmax, Phi1, Phi2
535
537{
539}
540
541////////////////////////////////////////////////////////////////////////////////
542/// Fills real parameters of a positioned box inside this. Returns 0 if successful.
543
544Int_t TGeoTrd2::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
545{
546 dx=dy=dz=0;
547 if (mat->IsRotation()) {
548 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
549 return 1; // ### rotation not accepted ###
550 }
551 //--> translate the origin of the parametrized box to the frame of this box.
552 Double_t origin[3];
553 mat->LocalToMaster(parambox->GetOrigin(), origin);
554 if (!Contains(origin)) {
555 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
556 return 1; // ### wrong matrix ###
557 }
558 //--> now we have to get the valid range for all parametrized axis
559 Double_t dd[3];
560 dd[0] = parambox->GetDX();
561 dd[1] = parambox->GetDY();
562 dd[2] = parambox->GetDZ();
563 //-> check if Z range is fixed
564 if (dd[2]<0) {
565 dd[2] = TMath::Min(origin[2]+fDz, fDz-origin[2]);
566 if (dd[2]<0) {
567 Error("GetFittingBox", "wrong matrix");
568 return 1;
569 }
570 }
571 if (dd[0]>=0 && dd[1]>=0) {
572 dx = dd[0];
573 dy = dd[1];
574 dz = dd[2];
575 return 0;
576 }
577 //-> check now range at Z = origin[2] +/- dd[2]
578 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
579 Double_t fy = 0.5*(fDy1-fDy2)/fDz;
580 Double_t dx0 = 0.5*(fDx1+fDx2);
581 Double_t dy0 = 0.5*(fDy1+fDy2);
582 Double_t z=origin[2]-dd[2];
583 dd[0] = dx0-fx*z-origin[0];
584 dd[1] = dy0-fy*z-origin[1];
585 z=origin[2]+dd[2];
586 dd[0] = TMath::Min(dd[0], dx0-fx*z-origin[0]);
587 dd[1] = TMath::Min(dd[1], dy0-fy*z-origin[1]);
588 if (dd[0]<0 || dd[1]<0) {
589 Error("GetFittingBox", "wrong matrix");
590 return 1;
591 }
592 dx = dd[0];
593 dy = dd[1];
594 dz = dd[2];
595 return 0;
596}
597
598////////////////////////////////////////////////////////////////////////////////
599/// in case shape has some negative parameters, these has to be computed
600/// in order to fit the mother
601
603{
604 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
605 if (!mother->TestShapeBit(kGeoTrd2)) {
606 Error("GetMakeRuntimeShape", "invalid mother");
607 return 0;
608 }
609 Double_t dx1, dx2, dy1, dy2, dz;
610 if (fDx1<0) dx1=((TGeoTrd2*)mother)->GetDx1();
611 else dx1=fDx1;
612 if (fDx2<0) dx2=((TGeoTrd2*)mother)->GetDx2();
613 else dx2=fDx2;
614 if (fDy1<0) dy1=((TGeoTrd2*)mother)->GetDy1();
615 else dy1=fDy1;
616 if (fDy2<0) dy2=((TGeoTrd2*)mother)->GetDy2();
617 else dy2=fDy2;
618 if (fDz<0) dz=((TGeoTrd2*)mother)->GetDz();
619 else dz=fDz;
620
621 return (new TGeoTrd2(dx1, dx2, dy1, dy2, dz));
622}
623
624////////////////////////////////////////////////////////////////////////////////
625/// print shape parameters
626
628{
629 printf("*** Shape %s: TGeoTrd2 ***\n", GetName());
630 printf(" dx1 = %11.5f\n", fDx1);
631 printf(" dx2 = %11.5f\n", fDx2);
632 printf(" dy1 = %11.5f\n", fDy1);
633 printf(" dy2 = %11.5f\n", fDy2);
634 printf(" dz = %11.5f\n", fDz);
635 printf(" Bounding box:\n");
637}
638
639////////////////////////////////////////////////////////////////////////////////
640/// computes the closest distance from given point to this shape, according
641/// to option. The matching point on the shape is stored in spoint.
642
644{
645 Double_t saf[3];
646 //--- Compute safety first
647 // check Z facettes
648 saf[0] = fDz-TMath::Abs(point[2]);
649 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
650 Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
651 // check X facettes
652 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
653 if (distx<0) saf[1]=TGeoShape::Big();
654 else saf[1]=(distx-TMath::Abs(point[0]))*calf;
655
656 Double_t fy = 0.5*(fDy1-fDy2)/fDz;
657 calf = 1./TMath::Sqrt(1.0+fy*fy);
658 // check Y facettes
659 distx = 0.5*(fDy1+fDy2)-fy*point[2];
660 if (distx<0) saf[2]=TGeoShape::Big();
661 else saf[2]=(distx-TMath::Abs(point[1]))*calf;
662
663 if (in) return saf[TMath::LocMin(3,saf)];
664 for (Int_t i=0; i<3; i++) saf[i]=-saf[i];
665 return saf[TMath::LocMax(3,saf)];
666}
667
668////////////////////////////////////////////////////////////////////////////////
669/// Save a primitive as a C++ statement(s) on output stream "out".
670
671void TGeoTrd2::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
672{
674 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
675 out << " dx1 = " << fDx1 << ";" << std::endl;
676 out << " dx2 = " << fDx2 << ";" << std::endl;
677 out << " dy1 = " << fDy1 << ";" << std::endl;
678 out << " dy2 = " << fDy2 << ";" << std::endl;
679 out << " dz = " << fDZ << ";" << std::endl;
680 out << " TGeoShape *" << GetPointerName() << " = new TGeoTrd2(\"" << GetName() << "\", dx1,dx2,dy1,dy2,dz);" << std::endl;
682}
683
684////////////////////////////////////////////////////////////////////////////////
685/// set arb8 params in one step :
686
688{
689 fDx1 = param[0];
690 fDx2 = param[1];
691 fDy1 = param[2];
692 fDy2 = param[3];
693 fDz = param[4];
694 ComputeBBox();
695}
696
697////////////////////////////////////////////////////////////////////////////////
698/// create trd2 mesh points
699
701{
702 if (!points) return;
703 points[ 0] = -fDx1; points[ 1] = -fDy1; points[ 2] = -fDz;
704 points[ 3] = -fDx1; points[ 4] = fDy1; points[ 5] = -fDz;
705 points[ 6] = fDx1; points[ 7] = fDy1; points[ 8] = -fDz;
706 points[ 9] = fDx1; points[10] = -fDy1; points[11] = -fDz;
707 points[12] = -fDx2; points[13] = -fDy2; points[14] = fDz;
708 points[15] = -fDx2; points[16] = fDy2; points[17] = fDz;
709 points[18] = fDx2; points[19] = fDy2; points[20] = fDz;
710 points[21] = fDx2; points[22] = -fDy2; points[23] = fDz;
711}
712
713////////////////////////////////////////////////////////////////////////////////
714/// create trd2 mesh points
715
717{
718 if (!points) return;
719 points[ 0] = -fDx1; points[ 1] = -fDy1; points[ 2] = -fDz;
720 points[ 3] = -fDx1; points[ 4] = fDy1; points[ 5] = -fDz;
721 points[ 6] = fDx1; points[ 7] = fDy1; points[ 8] = -fDz;
722 points[ 9] = fDx1; points[10] = -fDy1; points[11] = -fDz;
723 points[12] = -fDx2; points[13] = -fDy2; points[14] = fDz;
724 points[15] = -fDx2; points[16] = fDy2; points[17] = fDz;
725 points[18] = fDx2; points[19] = fDy2; points[20] = fDz;
726 points[21] = fDx2; points[22] = -fDy2; points[23] = fDz;
727}
728
729////////////////////////////////////////////////////////////////////////////////
730/// set vertex of a corner according to visibility flags
731
733{
734 if (TestShapeBit(kGeoVisX)) {
735 if (TestShapeBit(kGeoVisZ)) {
736 vertex[0] = fDx2;
737 vertex[2] = fDz;
739 } else {
740 vertex[0] = fDx1;
741 vertex[2] = -fDz;
743 }
744 } else {
745 if (TestShapeBit(kGeoVisZ)) {
746 vertex[0] = -fDx2;
747 vertex[2] = fDz;
749 } else {
750 vertex[0] = -fDx1;
751 vertex[2] = -fDz;
753 }
754 }
755}
756
757////////////////////////////////////////////////////////////////////////////////
758/// fill size of this 3-D object
759
761{
763}
764
765////////////////////////////////////////////////////////////////////////////////
766/// Check the inside status for each of the points in the array.
767/// Input: Array of point coordinates + vector size
768/// Output: Array of Booleans for the inside of each point
769
770void TGeoTrd2::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
771{
772 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
773}
774
775////////////////////////////////////////////////////////////////////////////////
776/// Compute the normal for an array o points so that norm.dot.dir is positive
777/// Input: Arrays of point coordinates and directions + vector size
778/// Output: Array of normal directions
779
780void TGeoTrd2::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
781{
782 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
783}
784
785////////////////////////////////////////////////////////////////////////////////
786/// Compute distance from array of input points having directions specified by dirs. Store output in dists
787
788void TGeoTrd2::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
789{
790 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
791}
792
793////////////////////////////////////////////////////////////////////////////////
794/// Compute distance from array of input points having directions specified by dirs. Store output in dists
795
796void TGeoTrd2::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
797{
798 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
799}
800
801////////////////////////////////////////////////////////////////////////////////
802/// Compute safe distance from each of the points in the input array.
803/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
804/// Output: Safety values
805
806void TGeoTrd2::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
807{
808 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
809}
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
XFontStruct * id
Definition: TGX11.cxx:108
char name[80]
Definition: TGX11.cxx:109
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:601
point * points
Definition: X3DBuffer.c:22
Box class.
Definition: TGeoBBox.h:18
virtual const Double_t * GetOrigin() const
Definition: TGeoBBox.h:73
Double_t fDX
Definition: TGeoBBox.h:21
virtual void InspectShape() const
Prints shape parameters.
Definition: TGeoBBox.cxx:793
virtual Double_t GetDX() const
Definition: TGeoBBox.h:70
virtual Double_t GetDZ() const
Definition: TGeoBBox.h:72
virtual void Sizeof3D() const
Definition: TGeoBBox.cxx:996
virtual Double_t GetDY() const
Definition: TGeoBBox.h:71
Double_t fOrigin[3]
Definition: TGeoBBox.h:24
Double_t fDY
Definition: TGeoBBox.h:22
Double_t fDZ
Definition: TGeoBBox.h:23
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition: TGeoBBox.cxx:583
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Geometrical transformation package.
Definition: TGeoMatrix.h:41
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
Definition: TGeoMatrix.cxx:339
Node containing an offset.
Definition: TGeoNode.h:182
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
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
Definition: TGeoShape.cxx:524
const char * GetPointerName() const
Provide a pointer name containing uid.
Definition: TGeoShape.cxx:699
virtual const char * GetName() const
Get the shape name.
Definition: TGeoShape.cxx:248
@ kGeoSavePrimitive
Definition: TGeoShape.h:65
@ kGeoVisX
Definition: TGeoShape.h:38
@ kGeoVisZ
Definition: TGeoShape.h:40
@ kGeoRunTimeShape
Definition: TGeoShape.h:41
@ kGeoVisY
Definition: TGeoShape.h:39
@ kGeoTrd2
Definition: TGeoShape.h:58
static Double_t Tolerance()
Definition: TGeoShape.h:91
Bool_t TestShapeBit(UInt_t f) const
Definition: TGeoShape.h:163
A trapezoid with both x and y lengths varying with z.
Definition: TGeoTrd2.h:18
Double_t fDx2
Definition: TGeoTrd2.h:22
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. Returns 0 if successful.
Definition: TGeoTrd2.cxx:544
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...
Definition: TGeoTrd2.cxx:780
virtual void SetPoints(Double_t *points) const
create trd2 mesh points
Definition: TGeoTrd2.cxx:700
Double_t fDy1
Definition: TGeoTrd2.h:23
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this trd2 shape belonging to volume "voldiv" into ndiv volumes called divname,...
Definition: TGeoTrd2.cxx:488
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: TGeoTrd2.cxx:602
virtual void InspectShape() const
print shape parameters
Definition: TGeoTrd2.cxx:627
virtual Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
Compute distance from inside point to surface of the trd2 Boundary safe algorithm.
Definition: TGeoTrd2.cxx:213
virtual ~TGeoTrd2()
destructor
Definition: TGeoTrd2.cxx:117
Double_t fDx1
Definition: TGeoTrd2.h:21
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...
Definition: TGeoTrd2.cxx:796
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition: TGeoTrd2.cxx:395
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...
Definition: TGeoTrd2.cxx:788
virtual Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0) const
Compute distance from outside point to surface of the trd2 Boundary safe algorithm.
Definition: TGeoTrd2.cxx:276
Double_t fDy2
Definition: TGeoTrd2.h:24
virtual void ComputeBBox()
compute bounding box for a trd2
Definition: TGeoTrd2.cxx:134
virtual void Sizeof3D() const
fill size of this 3-D object
Definition: TGeoTrd2.cxx:760
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this shape check Z range
Definition: TGeoTrd2.cxx:197
void GetOppositeCorner(const Double_t *point, Int_t inorm, Double_t *vertex, Double_t *normals) const
get the opposite corner of the intersected face
Definition: TGeoTrd2.cxx:460
void SetVertex(Double_t *vertex) const
set vertex of a corner according to visibility flags
Definition: TGeoTrd2.cxx:732
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
Definition: TGeoTrd2.cxx:643
virtual void SetDimensions(Double_t *param)
set arb8 params in one step :
Definition: TGeoTrd2.cxx:687
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.
Definition: TGeoTrd2.cxx:770
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition: TGeoTrd2.cxx:145
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TGeoTrd2.cxx:671
void GetVisibleCorner(const Double_t *point, Double_t *vertex, Double_t *normals) const
get the most visible corner from outside point and the normals
Definition: TGeoTrd2.cxx:413
Double_t fDz
Definition: TGeoTrd2.h:25
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition: TGeoTrd2.cxx:124
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition: TGeoTrd2.cxx:536
TGeoTrd2()
dummy ctor
Definition: TGeoTrd2.cxx:51
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.
Definition: TGeoTrd2.cxx:806
Volume families.
Definition: TGeoVolume.h:257
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:53
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:176
void SetFinder(TGeoPatternFinder *finder)
Definition: TGeoVolume.h:234
Int_t GetNdaughters() const
Definition: TGeoVolume.h:350
TObjArray * GetNodes()
Definition: TGeoVolume.h:170
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
static constexpr double s
Long64_t LocMin(Long64_t n, const T *a)
Return index of array with the minimum element.
Definition: TMath.h:960
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
T1 Sign(T1 a, T2 b)
Definition: TMathBase.h:165
Long64_t LocMax(Long64_t n, const T *a)
Return index of array with the maximum element.
Definition: TMath.h:988
Double_t Sqrt(Double_t x)
Definition: TMath.h:679
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
REAL * vertex
Definition: triangle.c:512