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