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