00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "Riostream.h"
00013
00014 #include "TGeoManager.h"
00015 #include "TGeoVolume.h"
00016 #include "TGeoArb8.h"
00017 #include "TGeoMatrix.h"
00018 #include "TMath.h"
00019
00020 ClassImp(TGeoArb8)
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104 TGeoArb8::TGeoArb8()
00105 {
00106
00107 fDz = 0;
00108 fTwist = 0;
00109 for (Int_t i=0; i<8; i++) {
00110 fXY[i][0] = 0.0;
00111 fXY[i][1] = 0.0;
00112 }
00113 SetShapeBit(kGeoArb8);
00114 }
00115
00116
00117 TGeoArb8::TGeoArb8(Double_t dz, Double_t *vertices)
00118 :TGeoBBox(0,0,0)
00119 {
00120
00121
00122 fDz = dz;
00123 fTwist = 0;
00124 SetShapeBit(kGeoArb8);
00125 if (vertices) {
00126 for (Int_t i=0; i<8; i++) {
00127 fXY[i][0] = vertices[2*i];
00128 fXY[i][1] = vertices[2*i+1];
00129 }
00130 ComputeTwist();
00131 ComputeBBox();
00132 } else {
00133 for (Int_t i=0; i<8; i++) {
00134 fXY[i][0] = 0.0;
00135 fXY[i][1] = 0.0;
00136 }
00137 }
00138 }
00139
00140
00141 TGeoArb8::TGeoArb8(const char *name, Double_t dz, Double_t *vertices)
00142 :TGeoBBox(name, 0,0,0)
00143 {
00144
00145
00146 fDz = dz;
00147 fTwist = 0;
00148 SetShapeBit(kGeoArb8);
00149 if (vertices) {
00150 for (Int_t i=0; i<8; i++) {
00151 fXY[i][0] = vertices[2*i];
00152 fXY[i][1] = vertices[2*i+1];
00153 }
00154 ComputeTwist();
00155 ComputeBBox();
00156 } else {
00157 for (Int_t i=0; i<8; i++) {
00158 fXY[i][0] = 0.0;
00159 fXY[i][1] = 0.0;
00160 }
00161 }
00162 }
00163
00164
00165 TGeoArb8::TGeoArb8(const TGeoArb8& ga8) :
00166 TGeoBBox(ga8),
00167 fDz(ga8.fDz),
00168 fTwist(ga8.fTwist)
00169 {
00170
00171 for(Int_t i=0; i<8; i++) {
00172 fXY[i][0]=ga8.fXY[i][0];
00173 fXY[i][1]=ga8.fXY[i][1];
00174 }
00175 }
00176
00177
00178 TGeoArb8& TGeoArb8::operator=(const TGeoArb8& ga8)
00179 {
00180
00181 if(this!=&ga8) {
00182 TGeoBBox::operator=(ga8);
00183 fDz=ga8.fDz;
00184 fTwist=ga8.fTwist;
00185 for(Int_t i=0; i<8; i++) {
00186 fXY[i][0]=ga8.fXY[i][0];
00187 fXY[i][1]=ga8.fXY[i][1];
00188 }
00189 }
00190 return *this;
00191 }
00192
00193
00194 TGeoArb8::~TGeoArb8()
00195 {
00196
00197 if (fTwist) delete [] fTwist;
00198 }
00199
00200
00201 Double_t TGeoArb8::Capacity() const
00202 {
00203
00204 Int_t i,j;
00205 Double_t capacity = 0;
00206 for (i=0; i<4; i++) {
00207 j = (i+1)%4;
00208 capacity += 0.25*fDz*((fXY[i][0]+fXY[i+4][0])*(fXY[j][1]+fXY[j+4][1]) -
00209 (fXY[j][0]+fXY[j+4][0])*(fXY[i][1]+fXY[i+4][1]) +
00210 (1./3)*((fXY[i+4][0]-fXY[i][0])*(fXY[j+4][1]-fXY[j][1]) -
00211 (fXY[j][0]-fXY[j+4][0])*(fXY[i][1]-fXY[i+4][1])));
00212 }
00213 return TMath::Abs(capacity);
00214 }
00215
00216
00217 void TGeoArb8::ComputeBBox()
00218 {
00219
00220 Double_t xmin, xmax, ymin, ymax;
00221 xmin = xmax = fXY[0][0];
00222 ymin = ymax = fXY[0][1];
00223
00224 for (Int_t i=1; i<8; i++) {
00225 if (xmin>fXY[i][0]) xmin=fXY[i][0];
00226 if (xmax<fXY[i][0]) xmax=fXY[i][0];
00227 if (ymin>fXY[i][1]) ymin=fXY[i][1];
00228 if (ymax<fXY[i][1]) ymax=fXY[i][1];
00229 }
00230 fDX = 0.5*(xmax-xmin);
00231 fDY = 0.5*(ymax-ymin);
00232 fDZ = fDz;
00233 fOrigin[0] = 0.5*(xmax+xmin);
00234 fOrigin[1] = 0.5*(ymax+ymin);
00235 fOrigin[2] = 0;
00236 SetShapeBit(kGeoClosedShape);
00237 }
00238
00239
00240 void TGeoArb8::ComputeTwist()
00241 {
00242
00243
00244
00245 Double_t twist[4];
00246 Bool_t twisted = kFALSE;
00247 Double_t dx1, dy1, dx2, dy2;
00248 for (Int_t i=0; i<4; i++) {
00249 dx1 = fXY[(i+1)%4][0]-fXY[i][0];
00250 dy1 = fXY[(i+1)%4][1]-fXY[i][1];
00251 if (TMath::Abs(dx1)<TGeoShape::Tolerance() && TMath::Abs(dy1)<TGeoShape::Tolerance()) {
00252 twist[i] = 0;
00253 continue;
00254 }
00255 dx2 = fXY[4+(i+1)%4][0]-fXY[4+i][0];
00256 dy2 = fXY[4+(i+1)%4][1]-fXY[4+i][1];
00257 if (TMath::Abs(dx2)<TGeoShape::Tolerance() && TMath::Abs(dy2)<TGeoShape::Tolerance()) {
00258 twist[i] = 0;
00259 continue;
00260 }
00261 twist[i] = dy1*dx2 - dx1*dy2;
00262 if (TMath::Abs(twist[i])<TGeoShape::Tolerance()) {
00263 twist[i] = 0;
00264 continue;
00265 }
00266 twist[i] = TMath::Sign(1.,twist[i]);
00267 twisted = kTRUE;
00268 }
00269 if (twisted) {
00270 if (fTwist) delete [] fTwist;
00271 fTwist = new Double_t[4];
00272 memcpy(fTwist, &twist[0], 4*sizeof(Double_t));
00273 }
00274 Double_t sum1 = 0.;
00275 Double_t sum2 = 0.;
00276 Int_t j;
00277 for (Int_t i=0; i<4; i++) {
00278 j = (i+1)%4;
00279 sum1 += fXY[i][0]*fXY[j][1]-fXY[j][0]*fXY[i][1];
00280 sum2 += fXY[i+4][0]*fXY[j+4][1]-fXY[j+4][0]*fXY[i+4][1];
00281 }
00282 if (sum1*sum2 < -TGeoShape::Tolerance()) {
00283 Fatal("ComputeTwist", "Shape %s type Arb8: Lower/upper faces defined with opposite clockwise", GetName());
00284 return;
00285 }
00286 if (sum1>0.) {
00287 Error("ComputeTwist", "Shape %s type Arb8: Vertices must be defined clockwise in XY planes. Re-ordering...", GetName());
00288 Double_t xtemp, ytemp;
00289 xtemp = fXY[1][0];
00290 ytemp = fXY[1][1];
00291 fXY[1][0] = fXY[3][0];
00292 fXY[1][1] = fXY[3][1];
00293 fXY[3][0] = xtemp;
00294 fXY[3][1] = ytemp;
00295 xtemp = fXY[5][0];
00296 ytemp = fXY[5][1];
00297 fXY[5][0] = fXY[7][0];
00298 fXY[5][1] = fXY[7][1];
00299 fXY[7][0] = xtemp;
00300 fXY[7][1] = ytemp;
00301 }
00302
00303 Bool_t illegal_cross = kFALSE;
00304 illegal_cross = TGeoShape::IsSegCrossing(fXY[0][0],fXY[0][1],fXY[1][0],fXY[1][1],
00305 fXY[2][0],fXY[2][1],fXY[3][0],fXY[3][1]);
00306 if (!illegal_cross)
00307 illegal_cross = TGeoShape::IsSegCrossing(fXY[4][0],fXY[4][1],fXY[5][0],fXY[5][1],
00308 fXY[6][0],fXY[6][1],fXY[7][0],fXY[7][1]);
00309 if (illegal_cross) {
00310 Error("ComputeTwist", "Shape %s type Arb8: Malformed polygon with crossing opposite segments", GetName());
00311 InspectShape();
00312 }
00313 }
00314
00315
00316 Double_t TGeoArb8::GetTwist(Int_t iseg) const
00317 {
00318
00319 if (!fTwist) return 0.;
00320 if (iseg<0 || iseg>3) return 0.;
00321 return fTwist[iseg];
00322 }
00323
00324
00325 Double_t TGeoArb8::GetClosestEdge(Double_t *point, Double_t *vert, Int_t &isegment) const
00326 {
00327
00328
00329
00330 isegment = 0;
00331 Int_t isegmin = 0;
00332 Int_t i1, i2;
00333 Double_t p1[2], p2[2];
00334 Double_t lsq, ssq, dx, dy, dpx, dpy, u;
00335 Double_t umin = -1.;
00336 Double_t safe=1E30;
00337 for (i1=0; i1<4; i1++) {
00338 if (TGeoShape::IsSameWithinTolerance(safe,0)) {
00339 isegment = isegmin;
00340 return umin;
00341 }
00342 i2 = (i1+1)%4;
00343 p1[0] = vert[2*i1];
00344 p1[1] = vert[2*i1+1];
00345 p2[0] = vert[2*i2];
00346 p2[1] = vert[2*i2+1];
00347 dx = p2[0] - p1[0];
00348 dy = p2[1] - p1[1];
00349 dpx = point[0] - p1[0];
00350 dpy = point[1] - p1[1];
00351 lsq = dx*dx + dy*dy;
00352
00353 if (TGeoShape::IsSameWithinTolerance(lsq,0)) {
00354 ssq = dpx*dpx + dpy*dpy;
00355 if (ssq < safe) {
00356 safe = ssq;
00357 isegmin = i1;
00358 umin = -1;
00359 }
00360 continue;
00361 }
00362
00363 u = (dpx*dx + dpy*dy)/lsq;
00364
00365 if (u>1) {
00366
00367 dpx = point[0]-p2[0];
00368 dpy = point[1]-p2[1];
00369 u = -1.;
00370 } else {
00371 if (u>=0) {
00372
00373 dpx -= u*dx;
00374 dpy -= u*dy;
00375 } else {
00376
00377 u = -1.;
00378 }
00379 }
00380 ssq = dpx*dpx + dpy*dpy;
00381 if (ssq < safe) {
00382 safe = ssq;
00383 isegmin = i1;
00384 umin = u;
00385 }
00386 }
00387 isegment = isegmin;
00388
00389 return umin;
00390 }
00391
00392
00393 void TGeoArb8::ComputeNormal(Double_t *point, Double_t *dir, Double_t *norm)
00394 {
00395
00396 Double_t safc;
00397 Double_t x0, y0, z0, x1, y1, z1, x2, y2, z2;
00398 Double_t ax, ay, az, bx, by, bz;
00399 Double_t fn;
00400 safc = fDz-TMath::Abs(point[2]);
00401 if (safc<1E-4) {
00402 memset(norm,0,3*sizeof(Double_t));
00403 norm[2] = (dir[2]>0)?1:(-1);
00404 return;
00405 }
00406 Double_t vert[8];
00407 SetPlaneVertices(point[2], vert);
00408
00409 Int_t iseg;
00410 Double_t frac = GetClosestEdge(point, vert, iseg);
00411 if (frac<0) frac = 0.;
00412 Int_t jseg = (iseg+1)%4;
00413 x0 = vert[2*iseg];
00414 y0 = vert[2*iseg+1];
00415 z0 = point[2];
00416 x2 = vert[2*jseg];
00417 y2 = vert[2*jseg+1];
00418 z2 = point[2];
00419 x0 += frac*(x2-x0);
00420 y0 += frac*(y2-y0);
00421 x1 = fXY[iseg+4][0];
00422 y1 = fXY[iseg+4][1];
00423 z1 = fDz;
00424 x1 += frac*(fXY[jseg+4][0]-x1);
00425 y1 += frac*(fXY[jseg+4][1]-y1);
00426 ax = x1-x0;
00427 ay = y1-y0;
00428 az = z1-z0;
00429 bx = x2-x0;
00430 by = y2-y0;
00431 bz = z2-z0;
00432
00433
00434 norm[0] = ay*bz-az*by;
00435 norm[1] = az*bx-ax*bz;
00436 norm[2] = ax*by-ay*bx;
00437 fn = TMath::Sqrt(norm[0]*norm[0]+norm[1]*norm[1]+norm[2]*norm[2]);
00438
00439 if (fn<1E-10) {
00440 norm[0] = 1.;
00441 norm[1] = 0.;
00442 norm[2] = 0.;
00443 } else {
00444 norm[0] /= fn;
00445 norm[1] /= fn;
00446 norm[2] /= fn;
00447 }
00448
00449 if (dir[0]*norm[0]+dir[1]*norm[1]+dir[2]*norm[2] < 0) {
00450 norm[0] = -norm[0];
00451 norm[1] = -norm[1];
00452 norm[2] = -norm[2];
00453 }
00454 }
00455
00456
00457 Bool_t TGeoArb8::Contains(Double_t *point) const
00458 {
00459
00460
00461 if (TMath::Abs(point[2]) > fDz) return kFALSE;
00462
00463 Double_t poly[8];
00464
00465
00466 Double_t cf = 0.5*(fDz-point[2])/fDz;
00467 Int_t i;
00468 for (i=0; i<4; i++) {
00469 poly[2*i] = fXY[i+4][0]+cf*(fXY[i][0]-fXY[i+4][0]);
00470 poly[2*i+1] = fXY[i+4][1]+cf*(fXY[i][1]-fXY[i+4][1]);
00471 }
00472 return InsidePolygon(point[0],point[1],poly);
00473 }
00474
00475
00476 Double_t TGeoArb8::DistToPlane(Double_t *point, Double_t *dir, Int_t ipl, Bool_t in) const
00477 {
00478
00479
00480
00481
00482
00483 Double_t xa,xb,xc,xd;
00484 Double_t ya,yb,yc,yd;
00485 Double_t eps = 100.*TGeoShape::Tolerance();
00486 Int_t j = (ipl+1)%4;
00487 xa=fXY[ipl][0];
00488 ya=fXY[ipl][1];
00489 xb=fXY[ipl+4][0];
00490 yb=fXY[ipl+4][1];
00491 xc=fXY[j][0];
00492 yc=fXY[j][1];
00493 xd=fXY[4+j][0];
00494 yd=fXY[4+j][1];
00495 Double_t dz2 =0.5/fDz;
00496 Double_t tx1 =dz2*(xb-xa);
00497 Double_t ty1 =dz2*(yb-ya);
00498 Double_t tx2 =dz2*(xd-xc);
00499 Double_t ty2 =dz2*(yd-yc);
00500 Double_t dzp =fDz+point[2];
00501 Double_t xs1 =xa+tx1*dzp;
00502 Double_t ys1 =ya+ty1*dzp;
00503 Double_t xs2 =xc+tx2*dzp;
00504 Double_t ys2 =yc+ty2*dzp;
00505 Double_t dxs =xs2-xs1;
00506 Double_t dys =ys2-ys1;
00507 Double_t dtx =tx2-tx1;
00508 Double_t dty =ty2-ty1;
00509 Double_t a=(dtx*dir[1]-dty*dir[0]+(tx1*ty2-tx2*ty1)*dir[2])*dir[2];
00510 Double_t b=dxs*dir[1]-dys*dir[0]+(dtx*point[1]-dty*point[0]+ty2*xs1-ty1*xs2
00511 +tx1*ys2-tx2*ys1)*dir[2];
00512 Double_t c=dxs*point[1]-dys*point[0]+xs1*ys2-xs2*ys1;
00513 Double_t s=TGeoShape::Big();
00514 Double_t x1,x2,y1,y2,xp,yp,zi;
00515 if (TMath::Abs(a)<eps) {
00516 if (TMath::Abs(b)<eps) return TGeoShape::Big();
00517 s=-c/b;
00518 if (s>eps) {
00519 if (in) return s;
00520 zi=point[2]+s*dir[2];
00521 if (TMath::Abs(zi)<fDz) {
00522 x1=xs1+tx1*dir[2]*s;
00523 x2=xs2+tx2*dir[2]*s;
00524 xp=point[0]+s*dir[0];
00525 y1=ys1+ty1*dir[2]*s;
00526 y2=ys2+ty2*dir[2]*s;
00527 yp=point[1]+s*dir[1];
00528 zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
00529 if (zi<=0) return s;
00530 }
00531 }
00532 return TGeoShape::Big();
00533 }
00534 Double_t d=b*b-4*a*c;
00535 if (d>=0) {
00536 if (a>0) s=0.5*(-b-TMath::Sqrt(d))/a;
00537 else s=0.5*(-b+TMath::Sqrt(d))/a;
00538 if (s>eps) {
00539 if (in) return s;
00540 zi=point[2]+s*dir[2];
00541 if (TMath::Abs(zi)<fDz) {
00542 x1=xs1+tx1*dir[2]*s;
00543 x2=xs2+tx2*dir[2]*s;
00544 xp=point[0]+s*dir[0];
00545 y1=ys1+ty1*dir[2]*s;
00546 y2=ys2+ty2*dir[2]*s;
00547 yp=point[1]+s*dir[1];
00548 zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
00549 if (zi<=0) return s;
00550 }
00551 }
00552 if (a>0) s=0.5*(-b+TMath::Sqrt(d))/a;
00553 else s=0.5*(-b-TMath::Sqrt(d))/a;
00554 if (s>eps) {
00555 if (in) return s;
00556 zi=point[2]+s*dir[2];
00557 if (TMath::Abs(zi)<fDz) {
00558 x1=xs1+tx1*dir[2]*s;
00559 x2=xs2+tx2*dir[2]*s;
00560 xp=point[0]+s*dir[0];
00561 y1=ys1+ty1*dir[2]*s;
00562 y2=ys2+ty2*dir[2]*s;
00563 yp=point[1]+s*dir[1];
00564 zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
00565 if (zi<=0) return s;
00566 }
00567 }
00568 }
00569 return TGeoShape::Big();
00570 }
00571
00572
00573 Double_t TGeoArb8::DistFromOutside(Double_t *point, Double_t *dir, Int_t , Double_t step, Double_t * ) const
00574 {
00575
00576 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
00577 if (sdist>=step) return TGeoShape::Big();
00578 Double_t dist[5];
00579
00580 Int_t i;
00581 for (i=0; i<4; i++) {
00582 dist[i]=DistToPlane(point, dir, i, kFALSE);
00583 }
00584
00585 dist[4]=TGeoShape::Big();
00586 if (TMath::Abs(point[2])>fDz) {
00587 if (!TGeoShape::IsSameWithinTolerance(dir[2],0)) {
00588 Double_t pt[3];
00589 if (point[2]>0) {
00590 dist[4] = (fDz-point[2])/dir[2];
00591 pt[2]=fDz;
00592 } else {
00593 dist[4] = (-fDz-point[2])/dir[2];
00594 pt[2]=-fDz;
00595 }
00596 if (dist[4]<0) {
00597 dist[4]=TGeoShape::Big();
00598 } else {
00599 for (Int_t j=0; j<2; j++) pt[j]=point[j]+dist[4]*dir[j];
00600 if (!Contains(&pt[0])) dist[4]=TGeoShape::Big();
00601 }
00602 }
00603 }
00604 Double_t distmin = dist[0];
00605 for (i=1;i<5;i++) if (dist[i] < distmin) distmin = dist[i];
00606 return distmin;
00607 }
00608
00609
00610 Double_t TGeoArb8::DistFromInside(Double_t *point, Double_t *dir, Int_t , Double_t , Double_t * ) const
00611 {
00612
00613 #ifdef OLDALGORITHM
00614 Int_t i;
00615 Double_t dist[6];
00616 dist[0]=dist[1]=TGeoShape::Big();
00617 if (dir[2]<0) {
00618 dist[0]=(-fDz-point[2])/dir[2];
00619 } else {
00620 if (dir[2]>0) dist[1]=(fDz-point[2])/dir[2];
00621 }
00622 for (i=0; i<4; i++) {
00623 dist[i+2]=DistToPlane(point, dir, i, kTRUE);
00624 }
00625
00626 Double_t distmin = dist[0];
00627 for (i=1;i<6;i++) if (dist[i] < distmin) distmin = dist[i];
00628 return distmin;
00629 #else
00630
00631
00632
00633
00634
00635 Double_t distmin;
00636 Bool_t lateral_cross = kFALSE;
00637 if (dir[2]<0) {
00638 distmin=(-fDz-point[2])/dir[2];
00639 } else {
00640 if (dir[2]>0) distmin =(fDz-point[2])/dir[2];
00641 else distmin = TGeoShape::Big();
00642 }
00643 Double_t dz2 =0.5/fDz;
00644 Double_t xa,xb,xc,xd;
00645 Double_t ya,yb,yc,yd;
00646 Double_t eps = 100.*TGeoShape::Tolerance();
00647 for (Int_t ipl=0;ipl<4;ipl++) {
00648 Int_t j = (ipl+1)%4;
00649 xa=fXY[ipl][0];
00650 ya=fXY[ipl][1];
00651 xb=fXY[ipl+4][0];
00652 yb=fXY[ipl+4][1];
00653 xc=fXY[j][0];
00654 yc=fXY[j][1];
00655 xd=fXY[4+j][0];
00656 yd=fXY[4+j][1];
00657
00658 Double_t tx1 =dz2*(xb-xa);
00659 Double_t ty1 =dz2*(yb-ya);
00660 Double_t tx2 =dz2*(xd-xc);
00661 Double_t ty2 =dz2*(yd-yc);
00662 Double_t dzp =fDz+point[2];
00663 Double_t xs1 =xa+tx1*dzp;
00664 Double_t ys1 =ya+ty1*dzp;
00665 Double_t xs2 =xc+tx2*dzp;
00666 Double_t ys2 =yc+ty2*dzp;
00667 Double_t dxs =xs2-xs1;
00668 Double_t dys =ys2-ys1;
00669 Double_t dtx =tx2-tx1;
00670 Double_t dty =ty2-ty1;
00671 Double_t a=(dtx*dir[1]-dty*dir[0]+(tx1*ty2-tx2*ty1)*dir[2])*dir[2];
00672 Double_t b=dxs*dir[1]-dys*dir[0]+(dtx*point[1]-dty*point[0]+ty2*xs1-ty1*xs2
00673 +tx1*ys2-tx2*ys1)*dir[2];
00674 Double_t c=dxs*point[1]-dys*point[0]+xs1*ys2-xs2*ys1;
00675 Double_t s=TGeoShape::Big();
00676 if (TMath::Abs(a)<eps) {
00677 if (TMath::Abs(b)<eps) continue;
00678 s=-c/b;
00679 if (s>eps && s < distmin) {
00680 distmin =s;
00681 lateral_cross=kTRUE;
00682 }
00683 continue;
00684 }
00685 Double_t d=b*b-4*a*c;
00686 if (d>=0.) {
00687 if (a>0) s=0.5*(-b-TMath::Sqrt(d))/a;
00688 else s=0.5*(-b+TMath::Sqrt(d))/a;
00689 if (s>eps) {
00690 if (s < distmin) {
00691 distmin = s;
00692 lateral_cross = kTRUE;
00693 }
00694 } else {
00695 if (a>0) s=0.5*(-b+TMath::Sqrt(d))/a;
00696 else s=0.5*(-b-TMath::Sqrt(d))/a;
00697 if (s>eps && s < distmin) {
00698 distmin =s;
00699 lateral_cross = kTRUE;
00700 }
00701 }
00702 }
00703 }
00704
00705 if (!lateral_cross) {
00706
00707 if (distmin > 1.E10) return TGeoShape::Tolerance();
00708 Double_t pt[2];
00709 pt[0] = point[0]+distmin*dir[0];
00710 pt[1] = point[1]+distmin*dir[1];
00711
00712 Double_t poly[8];
00713 Int_t i = 0;
00714 if (dir[2]>0.) i=4;
00715 for (Int_t j=0; j<4; j++) {
00716 poly[2*j] = fXY[j+i][0];
00717 poly[2*j+1] = fXY[j+i][1];
00718 }
00719 if (!InsidePolygon(pt[0],pt[1],poly)) return TGeoShape::Tolerance();
00720 }
00721 return distmin;
00722 #endif
00723 }
00724
00725
00726 TGeoVolume *TGeoArb8::Divide(TGeoVolume *voldiv, const char * , Int_t , Int_t ,
00727 Double_t , Double_t )
00728 {
00729
00730 Error("Divide", "Division of an arbitrary trapezoid not implemented");
00731 return voldiv;
00732 }
00733
00734
00735 Double_t TGeoArb8::GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
00736 {
00737
00738 xlo = 0;
00739 xhi = 0;
00740 Double_t dx = 0;
00741 if (iaxis==3) {
00742 xlo = -fDz;
00743 xhi = fDz;
00744 dx = xhi-xlo;
00745 return dx;
00746 }
00747 return dx;
00748 }
00749
00750
00751 void TGeoArb8::GetBoundingCylinder(Double_t *param) const
00752 {
00753
00754
00755
00756 Double_t rmaxsq = 0;
00757 Double_t rsq;
00758 Int_t i;
00759 for (i=0; i<8; i++) {
00760 rsq = fXY[i][0]*fXY[i][0] + fXY[i][1]*fXY[i][1];
00761 rmaxsq = TMath::Max(rsq, rmaxsq);
00762 }
00763 param[0] = 0.;
00764 param[1] = rmaxsq;
00765 param[2] = 0.;
00766 param[3] = 360.;
00767 }
00768
00769
00770 Int_t TGeoArb8::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
00771 {
00772
00773 dx=dy=dz=0;
00774 if (mat->IsRotation()) {
00775 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
00776 return 1;
00777 }
00778
00779 Double_t origin[3];
00780 mat->LocalToMaster(parambox->GetOrigin(), origin);
00781 if (!Contains(origin)) {
00782 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
00783 return 1;
00784 }
00785
00786 Double_t dd[3];
00787 dd[0] = parambox->GetDX();
00788 dd[1] = parambox->GetDY();
00789 dd[2] = parambox->GetDZ();
00790
00791 if (dd[2]<0) {
00792 dd[2] = TMath::Min(origin[2]+fDz, fDz-origin[2]);
00793 if (dd[2]<0) {
00794 Error("GetFittingBox", "wrong matrix");
00795 return 1;
00796 }
00797 }
00798 if (dd[0]>=0 && dd[1]>=0) {
00799 dx = dd[0];
00800 dy = dd[1];
00801 dz = dd[2];
00802 return 0;
00803 }
00804
00805 Double_t upper[8];
00806 Double_t lower[8];
00807 SetPlaneVertices(origin[2]-dd[2], lower);
00808 SetPlaneVertices(origin[2]+dd[2], upper);
00809 Double_t ddmin=TGeoShape::Big();
00810 for (Int_t iaxis=0; iaxis<2; iaxis++) {
00811 if (dd[iaxis]>=0) continue;
00812 ddmin=TGeoShape::Big();
00813 for (Int_t ivert=0; ivert<4; ivert++) {
00814 ddmin = TMath::Min(ddmin, TMath::Abs(origin[iaxis]-lower[2*ivert+iaxis]));
00815 ddmin = TMath::Min(ddmin, TMath::Abs(origin[iaxis]-upper[2*ivert+iaxis]));
00816 }
00817 dd[iaxis] = ddmin;
00818 }
00819 dx = dd[0];
00820 dy = dd[1];
00821 dz = dd[2];
00822 return 0;
00823 }
00824
00825
00826 void TGeoArb8::GetPlaneNormal(Double_t *p1, Double_t *p2, Double_t *p3, Double_t *norm)
00827 {
00828
00829 Double_t cross = 0.;
00830 Double_t v1[3], v2[3];
00831 Int_t i;
00832 for (i=0; i<3; i++) {
00833 v1[i] = p2[i] - p1[i];
00834 v2[i] = p3[i] - p1[i];
00835 }
00836 norm[0] = v1[1]*v2[2]-v1[2]*v2[1];
00837 cross += norm[0]*norm[0];
00838 norm[1] = v1[2]*v2[0]-v1[0]*v2[2];
00839 cross += norm[1]*norm[1];
00840 norm[2] = v1[0]*v2[1]-v1[1]*v2[0];
00841 cross += norm[2]*norm[2];
00842 if (TMath::Abs(cross) < TGeoShape::Tolerance()) return;
00843 cross = 1./TMath::Sqrt(cross);
00844 for (i=0; i<3; i++) norm[i] *= cross;
00845 }
00846
00847
00848 Bool_t TGeoArb8::GetPointsOnFacet(Int_t , Int_t , Double_t * ) const
00849 {
00850
00851
00852
00853
00854
00855
00856 return kFALSE;
00857
00858
00859
00860
00861
00862
00863
00864
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890
00891
00892
00893
00894
00895
00896
00897
00898
00899
00900
00901
00902
00903
00904
00905 }
00906
00907
00908 Bool_t TGeoArb8::InsidePolygon(Double_t x, Double_t y, Double_t *pts)
00909 {
00910
00911 Int_t i,j;
00912 Double_t x1,y1,x2,y2;
00913 Double_t cross;
00914 for (i=0; i<4; i++) {
00915 j = (i+1)%4;
00916 x1 = pts[i<<1];
00917 y1 = pts[(i<<1)+1];
00918 x2 = pts[j<<1];
00919 y2 = pts[(j<<1)+1];
00920 cross = (x-x1)*(y2-y1)-(y-y1)*(x2-x1);
00921 if (cross<0) return kFALSE;
00922 }
00923 return kTRUE;
00924 }
00925
00926
00927 void TGeoArb8::InspectShape() const
00928 {
00929
00930 printf("*** Shape %s: TGeoArb8 ***\n", GetName());
00931 if (IsTwisted()) printf(" = TWISTED\n");
00932 for (Int_t ip=0; ip<8; ip++) {
00933 printf(" point #%i : x=%11.5f y=%11.5f z=%11.5f\n",
00934 ip, fXY[ip][0], fXY[ip][1], fDz*((ip<4)?-1:1));
00935 }
00936 printf(" Bounding box:\n");
00937 TGeoBBox::InspectShape();
00938 }
00939
00940
00941 Double_t TGeoArb8::Safety(Double_t *point, Bool_t in) const
00942 {
00943
00944 Double_t safz = fDz-TMath::Abs(point[2]);
00945 if (!in) safz = -safz;
00946 Int_t iseg;
00947 Double_t safe = TGeoShape::Big();
00948 Double_t lsq, ssq, dx, dy, dpx, dpy, u;
00949 if (IsTwisted()) {
00950 if (!in) {
00951 if (!TGeoBBox::Contains(point)) return TGeoBBox::Safety(point,kFALSE);
00952 }
00953
00954
00955 Double_t vert[8];
00956 Double_t *p1, *p2;
00957 Int_t isegmin=0;
00958 Double_t umin = 0.;
00959 SetPlaneVertices (point[2], vert);
00960 for (iseg=0; iseg<4; iseg++) {
00961 if (safe<TGeoShape::Tolerance()) return 0.;
00962 p1 = &vert[2*iseg];
00963 p2 = &vert[2*((iseg+1)%4)];
00964 dx = p2[0] - p1[0];
00965 dy = p2[1] - p1[1];
00966 dpx = point[0] - p1[0];
00967 dpy = point[1] - p1[1];
00968
00969 lsq = dx*dx + dy*dy;
00970 u = (dpx*dx + dpy*dy)/lsq;
00971 if (u>1) {
00972 dpx = point[0]-p2[0];
00973 dpy = point[1]-p2[1];
00974 } else {
00975 if (u>=0) {
00976 dpx -= u*dx;
00977 dpy -= u*dy;
00978 }
00979 }
00980 ssq = dpx*dpx + dpy*dpy;
00981 if (ssq < safe) {
00982 isegmin = iseg;
00983 umin = u;
00984 safe = ssq;
00985 }
00986 }
00987 if (umin<0) umin = 0.;
00988 if (umin>1) {
00989 isegmin = (isegmin+1)%4;
00990 umin = 0.;
00991 }
00992 Int_t i1 = isegmin;
00993 Int_t i2 = (isegmin+1)%4;
00994 Double_t dx1 = fXY[i2][0]-fXY[i1][0];
00995 Double_t dx2 = fXY[i2+4][0]-fXY[i1+4][0];
00996 Double_t dy1 = fXY[i2][1]-fXY[i1][1];
00997 Double_t dy2 = fXY[i2+4][1]-fXY[i1+4][1];
00998 dx = dx1 + umin*(dx2-dx1);
00999 dy = dy1 + umin*(dy2-dy1);
01000 safe *= 1.- 4.*fDz*fDz/(dx*dx+dy*dy+4.*fDz*fDz);
01001 safe = TMath::Sqrt(safe);
01002 return safe;
01003 }
01004
01005 Double_t saf[5];
01006 saf[0] = safz;
01007
01008 for (iseg=0; iseg<4; iseg++) saf[iseg+1] = SafetyToFace(point,iseg,in);
01009 if (in) safe = saf[TMath::LocMin(5, saf)];
01010 else safe = saf[TMath::LocMax(5, saf)];
01011 if (safe<0) return 0.;
01012 return safe;
01013 }
01014
01015
01016 Double_t TGeoArb8::SafetyToFace(Double_t *point, Int_t iseg, Bool_t in) const
01017 {
01018
01019
01020 Double_t vertices[12];
01021 Int_t ipln = (iseg+1)%4;
01022
01023 vertices[0] = fXY[iseg][0];
01024 vertices[1] = fXY[iseg][1];
01025 vertices[2] = -fDz;
01026
01027 vertices[3] = fXY[ipln][0];
01028 vertices[4] = fXY[ipln][1];
01029 vertices[5] = -fDz;
01030
01031 vertices[6] = fXY[ipln+4][0];
01032 vertices[7] = fXY[ipln+4][1];
01033 vertices[8] = fDz;
01034
01035 vertices[9] = fXY[iseg+4][0];
01036 vertices[10] = fXY[iseg+4][1];
01037 vertices[11] = fDz;
01038 Double_t safe;
01039 Double_t norm[3];
01040 Double_t *p1, *p2, *p3;
01041 p1 = &vertices[0];
01042 p2 = &vertices[9];
01043 p3 = &vertices[6];
01044 if (IsSamePoint(p2,p3)) {
01045 p3 = &vertices[3];
01046 if (IsSamePoint(p1,p3)) return -TGeoShape::Big();
01047 }
01048 GetPlaneNormal(p1,p2,p3,norm);
01049 safe = (point[0]-p1[0])*norm[0]+(point[1]-p1[1])*norm[1]+(point[2]-p1[2])*norm[2];
01050 if (in) return (-safe);
01051 return safe;
01052 }
01053
01054
01055 void TGeoArb8::SavePrimitive(ostream &out, Option_t * )
01056 {
01057
01058 if (TObject::TestBit(kGeoSavePrimitive)) return;
01059 out << " // Shape: " << GetName() << " type: " << ClassName() << endl;
01060 out << " dz = " << fDz << ";" << endl;
01061 out << " vert[0] = " << fXY[0][0] << ";" << endl;
01062 out << " vert[1] = " << fXY[0][1] << ";" << endl;
01063 out << " vert[2] = " << fXY[1][0] << ";" << endl;
01064 out << " vert[3] = " << fXY[1][1] << ";" << endl;
01065 out << " vert[4] = " << fXY[2][0] << ";" << endl;
01066 out << " vert[5] = " << fXY[2][1] << ";" << endl;
01067 out << " vert[6] = " << fXY[3][0] << ";" << endl;
01068 out << " vert[7] = " << fXY[3][1] << ";" << endl;
01069 out << " vert[8] = " << fXY[4][0] << ";" << endl;
01070 out << " vert[9] = " << fXY[4][1] << ";" << endl;
01071 out << " vert[10] = " << fXY[5][0] << ";" << endl;
01072 out << " vert[11] = " << fXY[5][1] << ";" << endl;
01073 out << " vert[12] = " << fXY[6][0] << ";" << endl;
01074 out << " vert[13] = " << fXY[6][1] << ";" << endl;
01075 out << " vert[14] = " << fXY[7][0] << ";" << endl;
01076 out << " vert[15] = " << fXY[7][1] << ";" << endl;
01077 out << " TGeoShape *" << GetPointerName() << " = new TGeoArb8(\"" << GetName() << "\", dz,vert);" << endl;
01078 TObject::SetBit(TGeoShape::kGeoSavePrimitive);
01079 }
01080
01081
01082 void TGeoArb8::SetPlaneVertices(Double_t zpl, Double_t *vertices) const
01083 {
01084
01085 Double_t cf = 0.5*(fDz-zpl)/fDz;
01086 for (Int_t i=0; i<4; i++) {
01087 vertices[2*i] = fXY[i+4][0]+cf*(fXY[i][0]-fXY[i+4][0]);
01088 vertices[2*i+1] = fXY[i+4][1]+cf*(fXY[i][1]-fXY[i+4][1]);
01089 }
01090 }
01091
01092
01093 void TGeoArb8::SetDimensions(Double_t *param)
01094 {
01095
01096
01097
01098
01099
01100 fDz = param[0];
01101 for (Int_t i=0; i<8; i++) {
01102 fXY[i][0] = param[2*i+1];
01103 fXY[i][1] = param[2*i+2];
01104 }
01105 ComputeTwist();
01106 ComputeBBox();
01107 }
01108
01109
01110 void TGeoArb8::SetPoints(Double_t *points) const
01111 {
01112
01113 for (Int_t i=0; i<8; i++) {
01114 points[3*i] = fXY[i][0];
01115 points[3*i+1] = fXY[i][1];
01116 points[3*i+2] = (i<4)?-fDz:fDz;
01117 }
01118 }
01119
01120
01121 void TGeoArb8::SetPoints(Float_t *points) const
01122 {
01123
01124 for (Int_t i=0; i<8; i++) {
01125 points[3*i] = fXY[i][0];
01126 points[3*i+1] = fXY[i][1];
01127 points[3*i+2] = (i<4)?-fDz:fDz;
01128 }
01129 }
01130
01131
01132 void TGeoArb8::SetVertex(Int_t vnum, Double_t x, Double_t y)
01133 {
01134
01135 if (vnum<0 || vnum >7) {
01136 Error("SetVertex", "Invalid vertex number");
01137 return;
01138 }
01139 fXY[vnum][0] = x;
01140 fXY[vnum][1] = y;
01141 if (vnum == 7) {
01142 ComputeTwist();
01143 ComputeBBox();
01144 }
01145 }
01146
01147
01148 void TGeoArb8::Sizeof3D() const
01149 {
01150
01151 TGeoBBox::Sizeof3D();
01152 }
01153
01154
01155 void TGeoArb8::Streamer(TBuffer &R__b)
01156 {
01157
01158 if (R__b.IsReading()) {
01159 R__b.ReadClassBuffer(TGeoArb8::Class(), this);
01160 ComputeTwist();
01161 } else {
01162 R__b.WriteClassBuffer(TGeoArb8::Class(), this);
01163 }
01164 }
01165
01166 ClassImp(TGeoTrap)
01167
01168
01169 TGeoTrap::TGeoTrap()
01170 {
01171
01172 fDz = 0;
01173 fTheta = 0;
01174 fPhi = 0;
01175 fH1 = fH2 = fBl1 = fBl2 = fTl1 = fTl2 = fAlpha1 = fAlpha2 = 0;
01176 }
01177
01178
01179 TGeoTrap::TGeoTrap(Double_t dz, Double_t theta, Double_t phi)
01180 :TGeoArb8("", 0, 0)
01181 {
01182
01183 fDz = dz;
01184 fTheta = theta;
01185 fPhi = phi;
01186 fH1 = fH2 = fBl1 = fBl2 = fTl1 = fTl2 = fAlpha1 = fAlpha2 = 0;
01187 }
01188
01189
01190 TGeoTrap::TGeoTrap(Double_t dz, Double_t theta, Double_t phi, Double_t h1,
01191 Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
01192 Double_t tl2, Double_t alpha2)
01193 :TGeoArb8("", 0, 0)
01194 {
01195
01196 fDz = dz;
01197 fTheta = theta;
01198 fPhi = phi;
01199 fH1 = h1;
01200 fH2 = h2;
01201 fBl1 = bl1;
01202 fBl2 = bl2;
01203 fTl1 = tl1;
01204 fTl2 = tl2;
01205 fAlpha1 = alpha1;
01206 fAlpha2 = alpha2;
01207 Double_t tx = TMath::Tan(theta*TMath::DegToRad())*TMath::Cos(phi*TMath::DegToRad());
01208 Double_t ty = TMath::Tan(theta*TMath::DegToRad())*TMath::Sin(phi*TMath::DegToRad());
01209 Double_t ta1 = TMath::Tan(alpha1*TMath::DegToRad());
01210 Double_t ta2 = TMath::Tan(alpha2*TMath::DegToRad());
01211 fXY[0][0] = -dz*tx-h1*ta1-bl1; fXY[0][1] = -dz*ty-h1;
01212 fXY[1][0] = -dz*tx+h1*ta1-tl1; fXY[1][1] = -dz*ty+h1;
01213 fXY[2][0] = -dz*tx+h1*ta1+tl1; fXY[2][1] = -dz*ty+h1;
01214 fXY[3][0] = -dz*tx-h1*ta1+bl1; fXY[3][1] = -dz*ty-h1;
01215 fXY[4][0] = dz*tx-h2*ta2-bl2; fXY[4][1] = dz*ty-h2;
01216 fXY[5][0] = dz*tx+h2*ta2-tl2; fXY[5][1] = dz*ty+h2;
01217 fXY[6][0] = dz*tx+h2*ta2+tl2; fXY[6][1] = dz*ty+h2;
01218 fXY[7][0] = dz*tx-h2*ta2+bl2; fXY[7][1] = dz*ty-h2;
01219 ComputeTwist();
01220 if ((dz<0) || (h1<0) || (bl1<0) || (tl1<0) ||
01221 (h2<0) || (bl2<0) || (tl2<0)) {
01222 SetShapeBit(kGeoRunTimeShape);
01223 }
01224 else TGeoArb8::ComputeBBox();
01225 }
01226
01227
01228 TGeoTrap::TGeoTrap(const char *name, Double_t dz, Double_t theta, Double_t phi, Double_t h1,
01229 Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
01230 Double_t tl2, Double_t alpha2)
01231 :TGeoArb8(name, 0, 0)
01232 {
01233
01234 SetName(name);
01235 fDz = dz;
01236 fTheta = theta;
01237 fPhi = phi;
01238 fH1 = h1;
01239 fH2 = h2;
01240 fBl1 = bl1;
01241 fBl2 = bl2;
01242 fTl1 = tl1;
01243 fTl2 = tl2;
01244 fAlpha1 = alpha1;
01245 fAlpha2 = alpha2;
01246 for (Int_t i=0; i<8; i++) {
01247 fXY[i][0] = 0.0;
01248 fXY[i][1] = 0.0;
01249 }
01250 Double_t tx = TMath::Tan(theta*TMath::DegToRad())*TMath::Cos(phi*TMath::DegToRad());
01251 Double_t ty = TMath::Tan(theta*TMath::DegToRad())*TMath::Sin(phi*TMath::DegToRad());
01252 Double_t ta1 = TMath::Tan(alpha1*TMath::DegToRad());
01253 Double_t ta2 = TMath::Tan(alpha2*TMath::DegToRad());
01254 fXY[0][0] = -dz*tx-h1*ta1-bl1; fXY[0][1] = -dz*ty-h1;
01255 fXY[1][0] = -dz*tx+h1*ta1-tl1; fXY[1][1] = -dz*ty+h1;
01256 fXY[2][0] = -dz*tx+h1*ta1+tl1; fXY[2][1] = -dz*ty+h1;
01257 fXY[3][0] = -dz*tx-h1*ta1+bl1; fXY[3][1] = -dz*ty-h1;
01258 fXY[4][0] = dz*tx-h2*ta2-bl2; fXY[4][1] = dz*ty-h2;
01259 fXY[5][0] = dz*tx+h2*ta2-tl2; fXY[5][1] = dz*ty+h2;
01260 fXY[6][0] = dz*tx+h2*ta2+tl2; fXY[6][1] = dz*ty+h2;
01261 fXY[7][0] = dz*tx-h2*ta2+bl2; fXY[7][1] = dz*ty-h2;
01262 ComputeTwist();
01263 if ((dz<0) || (h1<0) || (bl1<0) || (tl1<0) ||
01264 (h2<0) || (bl2<0) || (tl2<0)) {
01265 SetShapeBit(kGeoRunTimeShape);
01266 }
01267 else TGeoArb8::ComputeBBox();
01268 }
01269
01270
01271 TGeoTrap::~TGeoTrap()
01272 {
01273
01274 }
01275
01276
01277 Double_t TGeoTrap::DistFromInside(Double_t *point, Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
01278 {
01279
01280 if (iact<3 && safe) {
01281
01282 *safe = Safety(point, kTRUE);
01283 if (iact==0) return TGeoShape::Big();
01284 if (iact==1 && step<*safe) return TGeoShape::Big();
01285 }
01286
01287
01288
01289
01290
01291
01292
01293
01294 Double_t distmin;
01295 if (dir[2]<0) {
01296 distmin=(-fDz-point[2])/dir[2];
01297 } else {
01298 if (dir[2]>0) distmin =(fDz-point[2])/dir[2];
01299 else distmin = TGeoShape::Big();
01300 }
01301 Double_t xa,xb,xc;
01302 Double_t ya,yb,yc;
01303 for (Int_t ipl=0;ipl<4;ipl++) {
01304 Int_t j = (ipl+1)%4;
01305 xa=fXY[ipl][0];
01306 ya=fXY[ipl][1];
01307 xb=fXY[ipl+4][0];
01308 yb=fXY[ipl+4][1];
01309 xc=fXY[j][0];
01310 yc=fXY[j][1];
01311 Double_t ax,ay,az;
01312 ax = xb-xa;
01313 ay = yb-ya;
01314 az = 2.*fDz;
01315 Double_t bx,by;
01316 bx = xc-xa;
01317 by = yc-ya;
01318 Double_t ddotn = -dir[0]*az*by + dir[1]*az*bx+dir[2]*(ax*by-ay*bx);
01319 if (ddotn<=0) continue;
01320 Double_t saf = -(point[0]-xa)*az*by + (point[1]-ya)*az*bx + (point[2]+fDz)*(ax*by-ay*bx);
01321 if (saf>=0.0) return 0.0;
01322 Double_t s = -saf/ddotn;
01323 if (s<distmin) distmin=s;
01324 }
01325 return distmin;
01326 }
01327
01328
01329 Double_t TGeoTrap::DistFromOutside(Double_t *point, Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
01330 {
01331
01332 if (iact<3 && safe) {
01333
01334 *safe = Safety(point, kFALSE);
01335 if (iact==0) return TGeoShape::Big();
01336 if (iact==1 && step<*safe) return TGeoShape::Big();
01337 }
01338
01339 Double_t sdist = TGeoBBox::DistFromOutside(point,dir, fDX, fDY, fDZ, fOrigin, step);
01340 if (sdist>=step) return TGeoShape::Big();
01341
01342 Bool_t in = kTRUE;
01343 Double_t pts[8];
01344 Double_t snxt;
01345 Double_t xnew, ynew, znew;
01346 Int_t i,j;
01347 if (point[2]<-fDz+TGeoShape::Tolerance()) {
01348 if (dir[2]<=0) return TGeoShape::Big();
01349 in = kFALSE;
01350 snxt = -(fDz+point[2])/dir[2];
01351 xnew = point[0] + snxt*dir[0];
01352 ynew = point[1] + snxt*dir[1];
01353 for (i=0;i<4;i++) {
01354 j = i<<1;
01355 pts[j] = fXY[i][0];
01356 pts[j+1] = fXY[i][1];
01357 }
01358 if (InsidePolygon(xnew,ynew,pts)) return snxt;
01359 } else if (point[2]>fDz-TGeoShape::Tolerance()) {
01360 if (dir[2]>=0) return TGeoShape::Big();
01361 in = kFALSE;
01362 snxt = (fDz-point[2])/dir[2];
01363 xnew = point[0] + snxt*dir[0];
01364 ynew = point[1] + snxt*dir[1];
01365 for (i=0;i<4;i++) {
01366 j = i<<1;
01367 pts[j] = fXY[i+4][0];
01368 pts[j+1] = fXY[i+4][1];
01369 }
01370 if (InsidePolygon(xnew,ynew,pts)) return snxt;
01371 }
01372 snxt = TGeoShape::Big();
01373
01374
01375
01376 Double_t dz2 =0.5/fDz;
01377 Double_t xa,xb,xc,xd;
01378 Double_t ya,yb,yc,yd;
01379 Double_t ax,ay,az;
01380 Double_t bx,by;
01381 Double_t ddotn, saf;
01382 Double_t safmin = TGeoShape::Big();
01383 Bool_t exiting = kFALSE;
01384 for (i=0; i<4; i++) {
01385 j = (i+1)%4;
01386 xa=fXY[i][0];
01387 ya=fXY[i][1];
01388 xb=fXY[i+4][0];
01389 yb=fXY[i+4][1];
01390 xc=fXY[j][0];
01391 yc=fXY[j][1];
01392 xd=fXY[4+j][0];
01393 yd=fXY[4+j][1];
01394 ax = xb-xa;
01395 ay = yb-ya;
01396 az = 2.*fDz;
01397 bx = xc-xa;
01398 by = yc-ya;
01399 ddotn = -dir[0]*az*by + dir[1]*az*bx+dir[2]*(ax*by-ay*bx);
01400 saf = (point[0]-xa)*az*by - (point[1]-ya)*az*bx - (point[2]+fDz)*(ax*by-ay*bx);
01401
01402 if (saf<=0) {
01403
01404 in = kFALSE;
01405 if (ddotn>=0) return TGeoShape::Big();
01406 snxt = saf/ddotn;
01407 znew = point[2]+snxt*dir[2];
01408 if (TMath::Abs(znew)<=fDz) {
01409 xnew = point[0]+snxt*dir[0];
01410 ynew = point[1]+snxt*dir[1];
01411 Double_t tx1 =dz2*(xb-xa);
01412 Double_t ty1 =dz2*(yb-ya);
01413 Double_t tx2 =dz2*(xd-xc);
01414 Double_t ty2 =dz2*(yd-yc);
01415 Double_t dzp =fDz+znew;
01416 Double_t xs1 =xa+tx1*dzp;
01417 Double_t ys1 =ya+ty1*dzp;
01418 Double_t xs2 =xc+tx2*dzp;
01419 Double_t ys2 =yc+ty2*dzp;
01420 if (TMath::Abs(xs1-xs2)>TMath::Abs(ys1-ys2)) {
01421 if ((xnew-xs1)*(xs2-xnew)>=0) return snxt;
01422 } else {
01423 if ((ynew-ys1)*(ys2-ynew)>=0) return snxt;
01424 }
01425 }
01426 } else {
01427 if (saf<safmin) {
01428 safmin = saf;
01429 if (ddotn>=0) exiting = kTRUE;
01430 else exiting = kFALSE;
01431 }
01432 }
01433 }
01434 if (!in) return TGeoShape::Big();
01435 if (exiting) return TGeoShape::Big();
01436 return 0.0;
01437 }
01438
01439
01440 TGeoVolume *TGeoTrap::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
01441 Double_t start, Double_t step)
01442 {
01443
01444
01445
01446
01447
01448 TGeoShape *shape;
01449 TGeoVolume *vol;
01450 TGeoVolumeMulti *vmulti;
01451 TGeoPatternFinder *finder;
01452 TString opt = "";
01453 if (iaxis!=3) {
01454 Error("Divide", "cannot divide trapezoids on other axis than Z");
01455 return 0;
01456 }
01457 Double_t end = start+ndiv*step;
01458 Double_t points_lo[8];
01459 Double_t points_hi[8];
01460 finder = new TGeoPatternTrapZ(voldiv, ndiv, start, end);
01461 voldiv->SetFinder(finder);
01462 finder->SetDivIndex(voldiv->GetNdaughters());
01463 opt = "Z";
01464 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
01465 Double_t txz = ((TGeoPatternTrapZ*)finder)->GetTxz();
01466 Double_t tyz = ((TGeoPatternTrapZ*)finder)->GetTyz();
01467 Double_t zmin, zmax, ox,oy,oz;
01468 for (Int_t idiv=0; idiv<ndiv; idiv++) {
01469 zmin = start+idiv*step;
01470 zmax = start+(idiv+1)*step;
01471 oz = start+idiv*step+step/2;
01472 ox = oz*txz;
01473 oy = oz*tyz;
01474 SetPlaneVertices(zmin, &points_lo[0]);
01475 SetPlaneVertices(zmax, &points_hi[0]);
01476 shape = new TGeoTrap(step/2, fTheta, fPhi);
01477 for (Int_t vert1=0; vert1<4; vert1++)
01478 ((TGeoArb8*)shape)->SetVertex(vert1, points_lo[2*vert1]-ox, points_lo[2*vert1+1]-oy);
01479 for (Int_t vert2=0; vert2<4; vert2++)
01480 ((TGeoArb8*)shape)->SetVertex(vert2+4, points_hi[2*vert2]-ox, points_hi[2*vert2+1]-oy);
01481 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
01482 vmulti->AddVolume(vol);
01483 voldiv->AddNodeOffset(vol, idiv, oz, opt.Data());
01484 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
01485 }
01486 return vmulti;
01487 }
01488
01489
01490 TGeoShape *TGeoTrap::GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix * ) const
01491 {
01492
01493
01494 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
01495 if (mother->IsRunTimeShape()) {
01496 Error("GetMakeRuntimeShape", "invalid mother");
01497 return 0;
01498 }
01499 Double_t dz, h1, bl1, tl1, h2, bl2, tl2;
01500 if (fDz<0) dz=((TGeoTrap*)mother)->GetDz();
01501 else dz=fDz;
01502
01503 if (fH1<0) h1 = ((TGeoTrap*)mother)->GetH1();
01504 else h1 = fH1;
01505
01506 if (fH2<0) h2 = ((TGeoTrap*)mother)->GetH2();
01507 else h2 = fH2;
01508
01509 if (fBl1<0) bl1 = ((TGeoTrap*)mother)->GetBl1();
01510 else bl1 = fBl1;
01511
01512 if (fBl2<0) bl2 = ((TGeoTrap*)mother)->GetBl2();
01513 else bl2 = fBl2;
01514
01515 if (fTl1<0) tl1 = ((TGeoTrap*)mother)->GetTl1();
01516 else tl1 = fTl1;
01517
01518 if (fTl2<0) tl2 = ((TGeoTrap*)mother)->GetTl2();
01519 else tl2 = fTl2;
01520
01521 return (new TGeoTrap(dz, fTheta, fPhi, h1, bl1, tl1, fAlpha1, h2, bl2, tl2, fAlpha2));
01522 }
01523
01524
01525 Double_t TGeoTrap::Safety(Double_t *point, Bool_t in) const
01526 {
01527
01528 Double_t safe = TGeoShape::Big();
01529 Double_t saf[5];
01530 Double_t norm[3];
01531 Int_t i,j;
01532 Double_t x0, y0, z0=-fDz, x1, y1, z1=fDz, x2, y2;
01533 Double_t ax, ay, az=z1-z0, bx, by;
01534 Double_t fn;
01535
01536 for (i=0; i<4; i++) {
01537 if (in) saf[i] = TGeoShape::Big();
01538 else saf[i] = 0.;
01539 x0 = fXY[i][0];
01540 y0 = fXY[i][1];
01541 x1 = fXY[i+4][0];
01542 y1 = fXY[i+4][1];
01543 ax = x1-x0;
01544 ay = y1-y0;
01545 az = z1-z0;
01546 j = (i+1)%4;
01547 x2 = fXY[j][0];
01548 y2 = fXY[j][1];
01549 bx = x2-x0;
01550 by = y2-y0;
01551 if (TMath::Abs(bx)<TGeoShape::Tolerance() && TMath::Abs(by)<TGeoShape::Tolerance()) {
01552 x2 = fXY[4+j][0];
01553 y2 = fXY[4+j][1];
01554 bx = x2-x1;
01555 by = y2-y1;
01556 if (TMath::Abs(bx)<TGeoShape::Tolerance() && TMath::Abs(by)<TGeoShape::Tolerance()) continue;
01557 }
01558 norm[0] = -az*by;
01559 norm[1] = az*bx;
01560 norm[2] = ax*by-ay*bx;
01561 fn = TMath::Sqrt(norm[0]*norm[0]+norm[1]*norm[1]+norm[2]*norm[2]);
01562 if (fn<1E-10) continue;
01563 saf[i] = (x0-point[0])*norm[0]+(y0-point[1])*norm[1]+(-fDz-point[2])*norm[2];
01564 if (in) {
01565 saf[i]=TMath::Abs(saf[i])/fn;
01566 } else {
01567 saf[i] = -saf[i]/fn;
01568 }
01569 }
01570 saf[4] = fDz-TMath::Abs(point[2]);
01571 if (in) {
01572 safe = saf[0];
01573 for (j=1;j<5;j++) if (saf[j] <safe) safe = saf[j];
01574 } else {
01575 saf[4]=-saf[4];
01576 safe = saf[0];
01577 for (j=1;j<5;j++) if (saf[j] >safe) safe = saf[j];
01578 }
01579 return safe;
01580 }
01581
01582
01583 void TGeoTrap::SavePrimitive(ostream &out, Option_t * )
01584 {
01585
01586 if (TObject::TestBit(kGeoSavePrimitive)) return;
01587 out << " // Shape: " << GetName() << " type: " << ClassName() << endl;
01588 out << " dz = " << fDz << ";" << endl;
01589 out << " theta = " << fTheta << ";" << endl;
01590 out << " phi = " << fPhi << ";" << endl;
01591 out << " h1 = " << fH1<< ";" << endl;
01592 out << " bl1 = " << fBl1<< ";" << endl;
01593 out << " tl1 = " << fTl1<< ";" << endl;
01594 out << " alpha1 = " << fAlpha1 << ";" << endl;
01595 out << " h2 = " << fH2 << ";" << endl;
01596 out << " bl2 = " << fBl2<< ";" << endl;
01597 out << " tl2 = " << fTl2<< ";" << endl;
01598 out << " alpha2 = " << fAlpha2 << ";" << endl;
01599 out << " TGeoShape *" << GetPointerName() << " = new TGeoTrap(\"" << GetName() << "\", dz,theta,phi,h1,bl1,tl1,alpha1,h2,bl2,tl2,alpha2);" << endl;
01600 TObject::SetBit(TGeoShape::kGeoSavePrimitive);
01601 }
01602
01603
01604 void TGeoTrap::SetDimensions(Double_t *param)
01605 {
01606
01607
01608
01609
01610
01611
01612
01613
01614
01615
01616
01617
01618 fDz = param[0];
01619 fTheta = param[1];
01620 fPhi = param[2];
01621 fH1 = param[3];
01622 fH2 = param[7];
01623 fBl1 = param[4];
01624 fBl2 = param[8];
01625 fTl1 = param[5];
01626 fTl2 = param[9];
01627 fAlpha1 = param[6];
01628 fAlpha2 = param[10];
01629 Double_t tx = TMath::Tan(fTheta*TMath::DegToRad())*TMath::Cos(fPhi*TMath::DegToRad());
01630 Double_t ty = TMath::Tan(fTheta*TMath::DegToRad())*TMath::Sin(fPhi*TMath::DegToRad());
01631 Double_t ta1 = TMath::Tan(fAlpha1*TMath::DegToRad());
01632 Double_t ta2 = TMath::Tan(fAlpha2*TMath::DegToRad());
01633 fXY[0][0] = -fDz*tx-fH1*ta1-fBl1; fXY[0][1] = -fDz*ty-fH1;
01634 fXY[1][0] = -fDz*tx+fH1*ta1-fTl1; fXY[1][1] = -fDz*ty+fH1;
01635 fXY[2][0] = -fDz*tx+fH1*ta1+fTl1; fXY[2][1] = -fDz*ty+fH1;
01636 fXY[3][0] = -fDz*tx-fH1*ta1+fBl1; fXY[3][1] = -fDz*ty-fH1;
01637 fXY[4][0] = fDz*tx-fH2*ta2-fBl2; fXY[4][1] = fDz*ty-fH2;
01638 fXY[5][0] = fDz*tx+fH2*ta2-fTl2; fXY[5][1] = fDz*ty+fH2;
01639 fXY[6][0] = fDz*tx+fH2*ta2+fTl2; fXY[6][1] = fDz*ty+fH2;
01640 fXY[7][0] = fDz*tx-fH2*ta2+fBl2; fXY[7][1] = fDz*ty-fH2;
01641 ComputeTwist();
01642 if ((fDz<0) || (fH1<0) || (fBl1<0) || (fTl1<0) ||
01643 (fH2<0) || (fBl2<0) || (fTl2<0)) {
01644 SetShapeBit(kGeoRunTimeShape);
01645 }
01646 else TGeoArb8::ComputeBBox();
01647 }
01648
01649 ClassImp(TGeoGtra)
01650
01651
01652 TGeoGtra::TGeoGtra()
01653 {
01654
01655 fTwistAngle = 0;
01656 }
01657
01658
01659 TGeoGtra::TGeoGtra(Double_t dz, Double_t theta, Double_t phi, Double_t twist, Double_t h1,
01660 Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
01661 Double_t tl2, Double_t alpha2)
01662 :TGeoTrap(dz, theta, phi, h1, bl1, tl1, alpha1, h2, bl2, tl2, alpha2)
01663 {
01664
01665 fTheta = theta;
01666 fPhi = phi;
01667 fH1 = h1;
01668 fH2 = h2;
01669 fBl1 = bl1;
01670 fBl2 = bl2;
01671 fTl1 = tl1;
01672 fTl2 = tl2;
01673 fAlpha1 = alpha1;
01674 fAlpha2 = alpha2;
01675 Double_t x, y, dx, dy, dx1, dx2, th, ph, al1, al2;
01676 al1 = alpha1*TMath::DegToRad();
01677 al2 = alpha2*TMath::DegToRad();
01678 th = theta*TMath::DegToRad();
01679 ph = phi*TMath::DegToRad();
01680 dx = 2*dz*TMath::Sin(th)*TMath::Cos(ph);
01681 dy = 2*dz*TMath::Sin(th)*TMath::Sin(ph);
01682 fDz = dz;
01683 dx1 = 2*h1*TMath::Tan(al1);
01684 dx2 = 2*h2*TMath::Tan(al2);
01685
01686 fTwistAngle = twist;
01687
01688 Int_t i;
01689 for (i=0; i<8; i++) {
01690 fXY[i][0] = 0.0;
01691 fXY[i][1] = 0.0;
01692 }
01693
01694 fXY[0][0] = -bl1; fXY[0][1] = -h1;
01695 fXY[1][0] = -tl1+dx1; fXY[1][1] = h1;
01696 fXY[2][0] = tl1+dx1; fXY[2][1] = h1;
01697 fXY[3][0] = bl1; fXY[3][1] = -h1;
01698 fXY[4][0] = -bl2+dx; fXY[4][1] = -h2+dy;
01699 fXY[5][0] = -tl2+dx+dx2; fXY[5][1] = h2+dy;
01700 fXY[6][0] = tl2+dx+dx2; fXY[6][1] = h2+dy;
01701 fXY[7][0] = bl2+dx; fXY[7][1] = -h2+dy;
01702 for (i=4; i<8; i++) {
01703 x = fXY[i][0];
01704 y = fXY[i][1];
01705 fXY[i][0] = x*TMath::Cos(twist*TMath::DegToRad()) + y*TMath::Sin(twist*TMath::DegToRad());
01706 fXY[i][1] = -x*TMath::Sin(twist*TMath::DegToRad()) + y*TMath::Cos(twist*TMath::DegToRad());
01707 }
01708 ComputeTwist();
01709 if ((dz<0) || (h1<0) || (bl1<0) || (tl1<0) ||
01710 (h2<0) || (bl2<0) || (tl2<0)) SetShapeBit(kGeoRunTimeShape);
01711 else TGeoArb8::ComputeBBox();
01712 }
01713
01714
01715 TGeoGtra::TGeoGtra(const char *name, Double_t dz, Double_t theta, Double_t phi, Double_t twist, Double_t h1,
01716 Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
01717 Double_t tl2, Double_t alpha2)
01718 :TGeoTrap(name, dz, theta, phi, h1, bl1, tl1, alpha1, h2, bl2, tl2, alpha2)
01719 {
01720
01721 SetName(name);
01722 fTheta = theta;
01723 fPhi = phi;
01724 fH1 = h1;
01725 fH2 = h2;
01726 fBl1 = bl1;
01727 fBl2 = bl2;
01728 fTl1 = tl1;
01729 fTl2 = tl2;
01730 fAlpha1 = alpha1;
01731 fAlpha2 = alpha2;
01732 Double_t x, y, dx, dy, dx1, dx2, th, ph, al1, al2;
01733 al1 = alpha1*TMath::DegToRad();
01734 al2 = alpha2*TMath::DegToRad();
01735 th = theta*TMath::DegToRad();
01736 ph = phi*TMath::DegToRad();
01737 dx = 2*dz*TMath::Sin(th)*TMath::Cos(ph);
01738 dy = 2*dz*TMath::Sin(th)*TMath::Sin(ph);
01739 fDz = dz;
01740 dx1 = 2*h1*TMath::Tan(al1);
01741 dx2 = 2*h2*TMath::Tan(al2);
01742
01743 fTwistAngle = twist;
01744
01745 Int_t i;
01746 for (i=0; i<8; i++) {
01747 fXY[i][0] = 0.0;
01748 fXY[i][1] = 0.0;
01749 }
01750
01751 fXY[0][0] = -bl1; fXY[0][1] = -h1;
01752 fXY[1][0] = -tl1+dx1; fXY[1][1] = h1;
01753 fXY[2][0] = tl1+dx1; fXY[2][1] = h1;
01754 fXY[3][0] = bl1; fXY[3][1] = -h1;
01755 fXY[4][0] = -bl2+dx; fXY[4][1] = -h2+dy;
01756 fXY[5][0] = -tl2+dx+dx2; fXY[5][1] = h2+dy;
01757 fXY[6][0] = tl2+dx+dx2; fXY[6][1] = h2+dy;
01758 fXY[7][0] = bl2+dx; fXY[7][1] = -h2+dy;
01759 for (i=4; i<8; i++) {
01760 x = fXY[i][0];
01761 y = fXY[i][1];
01762 fXY[i][0] = x*TMath::Cos(twist*TMath::DegToRad()) + y*TMath::Sin(twist*TMath::DegToRad());
01763 fXY[i][1] = -x*TMath::Sin(twist*TMath::DegToRad()) + y*TMath::Cos(twist*TMath::DegToRad());
01764 }
01765 ComputeTwist();
01766 if ((dz<0) || (h1<0) || (bl1<0) || (tl1<0) ||
01767 (h2<0) || (bl2<0) || (tl2<0)) SetShapeBit(kGeoRunTimeShape);
01768 else TGeoArb8::ComputeBBox();
01769 }
01770
01771
01772 TGeoGtra::~TGeoGtra()
01773 {
01774
01775 }
01776
01777
01778 Double_t TGeoGtra::DistFromInside(Double_t *point, Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
01779 {
01780
01781 if (iact<3 && safe) {
01782
01783 *safe = Safety(point, kTRUE);
01784 if (iact==0) return TGeoShape::Big();
01785 if (iact==1 && step<*safe) return TGeoShape::Big();
01786 }
01787
01788 return TGeoArb8::DistFromInside(point, dir, iact, step, safe);
01789 }
01790
01791
01792 Double_t TGeoGtra::DistFromOutside(Double_t *point, Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
01793 {
01794
01795 if (iact<3 && safe) {
01796
01797 *safe = Safety(point, kTRUE);
01798 if (iact==0) return TGeoShape::Big();
01799 if (iact==1 && step<*safe) return TGeoShape::Big();
01800 }
01801
01802 return TGeoArb8::DistFromOutside(point, dir, iact, step, safe);
01803 }
01804
01805
01806 TGeoShape *TGeoGtra::GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix * ) const
01807 {
01808
01809
01810 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
01811 if (mother->IsRunTimeShape()) {
01812 Error("GetMakeRuntimeShape", "invalid mother");
01813 return 0;
01814 }
01815 Double_t dz, h1, bl1, tl1, h2, bl2, tl2;
01816 if (fDz<0) dz=((TGeoTrap*)mother)->GetDz();
01817 else dz=fDz;
01818 if (fH1<0)
01819 h1 = ((TGeoTrap*)mother)->GetH1();
01820 else
01821 h1 = fH1;
01822 if (fH2<0)
01823 h2 = ((TGeoTrap*)mother)->GetH2();
01824 else
01825 h2 = fH2;
01826 if (fBl1<0)
01827 bl1 = ((TGeoTrap*)mother)->GetBl1();
01828 else
01829 bl1 = fBl1;
01830 if (fBl2<0)
01831 bl2 = ((TGeoTrap*)mother)->GetBl2();
01832 else
01833 bl2 = fBl2;
01834 if (fTl1<0)
01835 tl1 = ((TGeoTrap*)mother)->GetTl1();
01836 else
01837 tl1 = fTl1;
01838 if (fTl2<0)
01839 tl2 = ((TGeoTrap*)mother)->GetTl2();
01840 else
01841 tl2 = fTl2;
01842 return (new TGeoGtra(dz, fTheta, fPhi, fTwistAngle ,h1, bl1, tl1, fAlpha1, h2, bl2, tl2, fAlpha2));
01843 }
01844
01845
01846 void TGeoGtra::SavePrimitive(ostream &out, Option_t * )
01847 {
01848
01849 if (TObject::TestBit(kGeoSavePrimitive)) return;
01850 out << " // Shape: " << GetName() << " type: " << ClassName() << endl;
01851 out << " dz = " << fDz << ";" << endl;
01852 out << " theta = " << fTheta << ";" << endl;
01853 out << " phi = " << fPhi << ";" << endl;
01854 out << " twist = " << fTwistAngle << ";" << endl;
01855 out << " h1 = " << fH1<< ";" << endl;
01856 out << " bl1 = " << fBl1<< ";" << endl;
01857 out << " tl1 = " << fTl1<< ";" << endl;
01858 out << " alpha1 = " << fAlpha1 << ";" << endl;
01859 out << " h2 = " << fH2 << ";" << endl;
01860 out << " bl2 = " << fBl2<< ";" << endl;
01861 out << " tl2 = " << fTl2<< ";" << endl;
01862 out << " alpha2 = " << fAlpha2 << ";" << endl;
01863 out << " TGeoShape *" << GetPointerName() << " = new TGeoGtra(\"" << GetName() << "\", dz,theta,phi,twist,h1,bl1,tl1,alpha1,h2,bl2,tl2,alpha2);" << endl;
01864 TObject::SetBit(TGeoShape::kGeoSavePrimitive);
01865 }
01866
01867
01868 void TGeoGtra::SetDimensions(Double_t *param)
01869 {
01870
01871
01872
01873
01874
01875
01876
01877
01878
01879
01880
01881
01882
01883 fDz = param[0];
01884 fTheta = param[1];
01885 fPhi = param[2];
01886 fH1 = param[3];
01887 fH2 = param[7];
01888 fBl1 = param[4];
01889 fBl2 = param[8];
01890 fTl1 = param[5];
01891 fTl2 = param[9];
01892 fAlpha1 = param[6];
01893 fAlpha2 = param[10];
01894 fTwistAngle = param[11];
01895 Double_t x, y, dx, dy, dx1, dx2, th, ph, al1, al2;
01896 al1 = fAlpha1*TMath::DegToRad();
01897 al2 = fAlpha2*TMath::DegToRad();
01898 th = fTheta*TMath::DegToRad();
01899 ph = fPhi*TMath::DegToRad();
01900 dx = 2*fDz*TMath::Sin(th)*TMath::Cos(ph);
01901 dy = 2*fDz*TMath::Sin(th)*TMath::Sin(ph);
01902 dx1 = 2*fH1*TMath::Tan(al1);
01903 dx2 = 2*fH2*TMath::Tan(al2);
01904 Int_t i;
01905 for (i=0; i<8; i++) {
01906 fXY[i][0] = 0.0;
01907 fXY[i][1] = 0.0;
01908 }
01909
01910 fXY[0][0] = -fBl1; fXY[0][1] = -fH1;
01911 fXY[1][0] = -fTl1+dx1; fXY[1][1] = fH1;
01912 fXY[2][0] = fTl1+dx1; fXY[2][1] = fH1;
01913 fXY[3][0] = fBl1; fXY[3][1] = -fH1;
01914 fXY[4][0] = -fBl2+dx; fXY[4][1] = -fH2+dy;
01915 fXY[5][0] = -fTl2+dx+dx2; fXY[5][1] = fH2+dy;
01916 fXY[6][0] = fTl2+dx+dx2; fXY[6][1] = fH2+dy;
01917 fXY[7][0] = fBl2+dx; fXY[7][1] = -fH2+dy;
01918 for (i=4; i<8; i++) {
01919 x = fXY[i][0];
01920 y = fXY[i][1];
01921 fXY[i][0] = x*TMath::Cos(fTwistAngle*TMath::DegToRad()) + y*TMath::Sin(fTwistAngle*TMath::DegToRad());
01922 fXY[i][1] = -x*TMath::Sin(fTwistAngle*TMath::DegToRad()) + y*TMath::Cos(fTwistAngle*TMath::DegToRad());
01923 }
01924 ComputeTwist();
01925 if ((fDz<0) || (fH1<0) || (fBl1<0) || (fTl1<0) ||
01926 (fH2<0) || (fBl2<0) || (fTl2<0)) SetShapeBit(kGeoRunTimeShape);
01927 else TGeoArb8::ComputeBBox();
01928 }