// @(#)root/geom:$Id$
// Author: Andrei Gheata   31/01/02
// TGeoTrd2::Contains() and DistFromInside() implemented by Mihaela Gheata

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

//_____________________________________________________________________________
// TGeoTrd2 - a trapezoid with both x and y lengths varying with z. It
//   has 5 parameters, the half lengths in x at -dz and +dz, the half
//  lengths in y at -dz and +dz, and the half length in z (dz).
//
//_____________________________________________________________________________
//Begin_Html
/*
<img src="gif/t_trd2.gif">
*/
//End_Html

//Begin_Html
/*
<img src="gif/t_trd2divZ.gif">
*/
//End_Html

//Begin_Html
/*
<img src="gif/t_trd2divstepZ.gif">
*/
//End_Html

#include "Riostream.h"

#include "TGeoManager.h"
#include "TGeoMatrix.h"
#include "TGeoVolume.h"
#include "TGeoTrd2.h"
#include "TMath.h"

ClassImp(TGeoTrd2)

//_____________________________________________________________________________
TGeoTrd2::TGeoTrd2()
{
   // dummy ctor
   SetShapeBit(kGeoTrd2);
   fDz = fDx1 = fDx2 = fDy1 = fDy2 = 0;
}

//_____________________________________________________________________________
TGeoTrd2::TGeoTrd2(Double_t dx1, Double_t dx2, Double_t dy1, Double_t dy2, Double_t dz)
         :TGeoBBox(0,0,0)
{
// constructor.
   SetShapeBit(kGeoTrd2);
   fDx1 = dx1;
   fDx2 = dx2;
   fDy1 = dy1;
   fDy2 = dy2;
   fDz = dz;
   if ((fDx1<0) || (fDx2<0) || (fDy1<0) || (fDy2<0) || (fDz<0)) {
      SetShapeBit(kGeoRunTimeShape);
      printf("trd2 : dx1=%f, dx2=%f, dy1=%f, dy2=%f, dz=%f\n",
              dx1,dx2,dy1,dy2,dz);
   }
   else ComputeBBox();
}

//_____________________________________________________________________________
TGeoTrd2::TGeoTrd2(const char * name, Double_t dx1, Double_t dx2, Double_t dy1, Double_t dy2, Double_t dz)
         :TGeoBBox(name, 0,0,0)
{
// constructor.
   SetShapeBit(kGeoTrd2);
   fDx1 = dx1;
   fDx2 = dx2;
   fDy1 = dy1;
   fDy2 = dy2;
   fDz = dz;
   if ((fDx1<0) || (fDx2<0) || (fDy1<0) || (fDy2<0) || (fDz<0)) {
      SetShapeBit(kGeoRunTimeShape);
      printf("trd2 : dx1=%f, dx2=%f, dy1=%f, dy2=%f, dz=%f\n",
              dx1,dx2,dy1,dy2,dz);
   }
   else ComputeBBox();
}

//_____________________________________________________________________________
TGeoTrd2::TGeoTrd2(Double_t *param)
         :TGeoBBox(0,0,0)
{
   // ctor with an array of parameters
   // param[0] = dx1
   // param[1] = dx2
   // param[2] = dy1
   // param[3] = dy2
   // param[4] = dz
   SetShapeBit(kGeoTrd2);
   SetDimensions(param);
   if ((fDx1<0) || (fDx2<0) || (fDy1<0) || (fDy2<0) || (fDz<0)) SetShapeBit(kGeoRunTimeShape);
   else ComputeBBox();
}

//_____________________________________________________________________________
TGeoTrd2::~TGeoTrd2()
{
// destructor
}

//_____________________________________________________________________________
Double_t TGeoTrd2::Capacity() const
{
// Computes capacity of the shape in [length^3]
   Double_t capacity = 2*(fDx1+fDx2)*(fDy1+fDy2)*fDz +
                      (2./3.)*(fDx1-fDx2)*(fDy1-fDy2)*fDz;
   return capacity;
}

//_____________________________________________________________________________
void TGeoTrd2::ComputeBBox()
{
// compute bounding box for a trd2
   fDX = TMath::Max(fDx1, fDx2);
   fDY = TMath::Max(fDy1, fDy2);
   fDZ = fDz;
   memset(fOrigin, 0, 3*sizeof(Double_t));
}

//_____________________________________________________________________________
void TGeoTrd2::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
{
// Compute normal to closest surface from POINT.
   Double_t safe, safemin;
   Double_t fx = 0.5*(fDx1-fDx2)/fDz;
   Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
   // check Z facettes
   safe = safemin = TMath::Abs(fDz-TMath::Abs(point[2]));
   norm[0] = norm[1] = 0;
   norm[2] = (dir[2]>=0)?1:-1;
   if (safe<TGeoShape::Tolerance()) return;
   // check X facettes
   Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
   if (distx>=0) {
      safe=TMath::Abs(distx-TMath::Abs(point[0]))*calf;
      if (safe<safemin) {
         safemin = safe;
         norm[0] = (point[0]>0)?calf:(-calf);
         norm[1] = 0;
         norm[2] = calf*fx;
         Double_t dot = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
         if (dot<0) {
            norm[0]=-norm[0];
            norm[2]=-norm[2];
         }
         if (safe<TGeoShape::Tolerance()) return;
      }
   }

   Double_t fy = 0.5*(fDy1-fDy2)/fDz;
   calf = 1./TMath::Sqrt(1.0+fy*fy);

   // check Y facettes
   distx = 0.5*(fDy1+fDy2)-fy*point[2];
   if (distx>=0) {
      safe=TMath::Abs(distx-TMath::Abs(point[1]))*calf;
      if (safe<safemin) {
         norm[0] = 0;
         norm[1] = (point[1]>0)?calf:(-calf);
         norm[2] = calf*fy;
         Double_t dot = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
         if (dot<0) {
            norm[1]=-norm[1];
            norm[2]=-norm[2];
         }
      }
   }
}

//_____________________________________________________________________________
Bool_t TGeoTrd2::Contains(const Double_t *point) const
{
// test if point is inside this shape
   // check Z range
   if (TMath::Abs(point[2]) > fDz) return kFALSE;
   // then y
   Double_t dy = 0.5*(fDy2*(point[2]+fDz)+fDy1*(fDz-point[2]))/fDz;
   if (TMath::Abs(point[1]) > dy) return kFALSE;
   // then x
   Double_t dx = 0.5*(fDx2*(point[2]+fDz)+fDx1*(fDz-point[2]))/fDz;
   if (TMath::Abs(point[0]) > dx) return kFALSE;
   return kTRUE;
}

//_____________________________________________________________________________
Double_t TGeoTrd2::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
{
// Compute distance from inside point to surface of the trd2
// Boundary safe algorithm
   Double_t snxt = TGeoShape::Big();
   if (iact<3 && safe) {
   // compute safe distance
      *safe = Safety(point, kTRUE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }

   Double_t fx = 0.5*(fDx1-fDx2)/fDz;
   Double_t fy = 0.5*(fDy1-fDy2)/fDz;
   Double_t cn;

   Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
   Double_t disty = 0.5*(fDy1+fDy2)-fy*point[2];

   //--- Compute distance to this shape
   // first check if Z facettes are crossed
   Double_t dist[3];
   for (Int_t i=0; i<3; i++) dist[i]=TGeoShape::Big();
   if (dir[2]<0) {
      dist[0]=-(point[2]+fDz)/dir[2];
   } else if (dir[2]>0) {
      dist[0]=(fDz-point[2])/dir[2];
   }
   if (dist[0]<=0) return 0.0;
   // now check X facettes
   cn = -dir[0]+fx*dir[2];
   if (cn>0) {
      dist[1] = point[0]+distx;
      if (dist[1]<=0) return 0.0;
      dist[1] /= cn;
   }
   cn = dir[0]+fx*dir[2];
   if (cn>0) {
      Double_t s = distx-point[0];
      if (s<=0) return 0.0;
      s /= cn;
      if (s<dist[1]) dist[1] = s;
   }
   // now check Y facettes
   cn = -dir[1]+fy*dir[2];
   if (cn>0) {
      dist[2] = point[1]+disty;
      if (dist[2]<=0) return 0.0;
      dist[2] /= cn;
   }
   cn = dir[1]+fy*dir[2];
   if (cn>0) {
      Double_t s = disty-point[1];
      if (s<=0) return 0.0;
      s /= cn;
      if (s<dist[2]) dist[2] = s;
   }
   snxt = dist[TMath::LocMin(3,dist)];
   return snxt;
}

//_____________________________________________________________________________
Double_t TGeoTrd2::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
{
// Compute distance from outside point to surface of the trd2
// Boundary safe algorithm
   Double_t snxt = TGeoShape::Big();
   if (iact<3 && safe) {
   // compute safe distance
      *safe = Safety(point, kFALSE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }
   // find a visible face
   Double_t xnew,ynew,znew;
   Double_t fx = 0.5*(fDx1-fDx2)/fDz;
   Double_t fy = 0.5*(fDy1-fDy2)/fDz;
   Double_t cn;
   // check visibility of X faces
   Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
   Double_t disty = 0.5*(fDy1+fDy2)-fy*point[2];
   Bool_t in = kTRUE;
   Double_t safx = distx-TMath::Abs(point[0]);
   Double_t safy = disty-TMath::Abs(point[1]);
   Double_t safz = fDz-TMath::Abs(point[2]);
   //--- Compute distance to this shape
   // first check if Z facettes are crossed
   if (point[2]<=-fDz) {
      cn = -dir[2];
      if (cn>=0) return TGeoShape::Big();
      in = kFALSE;
      snxt = (fDz+point[2])/cn;
      // find extrapolated X and Y
      xnew = point[0]+snxt*dir[0];
      if (TMath::Abs(xnew) < fDx1) {
         ynew = point[1]+snxt*dir[1];
         if (TMath::Abs(ynew) < fDy1) return snxt;
      }
   } else if (point[2]>=fDz) {
      cn = dir[2];
      if (cn>=0) return TGeoShape::Big();
      in = kFALSE;
      snxt = (fDz-point[2])/cn;
      // find extrapolated X and Y
      xnew = point[0]+snxt*dir[0];
      if (TMath::Abs(xnew) < fDx2) {
         ynew = point[1]+snxt*dir[1];
         if (TMath::Abs(ynew) < fDy2) return snxt;
      }
   }
   // check if X facettes are crossed
   if (point[0]<=-distx) {
      cn = -dir[0]+fx*dir[2];
      if (cn>=0) return TGeoShape::Big();
      in = kFALSE;
      snxt = (point[0]+distx)/cn;
      // find extrapolated Y and Z
      znew = point[2]+snxt*dir[2];
      if (TMath::Abs(znew) < fDz) {
         Double_t dy = 0.5*(fDy1+fDy2)-fy*znew;
         ynew = point[1]+snxt*dir[1];
         if (TMath::Abs(ynew) < dy) return snxt;
      }
   }
   if (point[0]>=distx) {
      cn = dir[0]+fx*dir[2];
      if (cn>=0) return TGeoShape::Big();
      in = kFALSE;
      snxt = (distx-point[0])/cn;
      // find extrapolated Y and Z
      znew = point[2]+snxt*dir[2];
      if (TMath::Abs(znew) < fDz) {
         Double_t dy = 0.5*(fDy1+fDy2)-fy*znew;
         ynew = point[1]+snxt*dir[1];
         if (TMath::Abs(ynew) < dy) return snxt;
      }
   }
   // finally check Y facettes
   if (point[1]<=-disty) {
      cn = -dir[1]+fy*dir[2];
      in = kFALSE;
      if (cn>=0) return TGeoShape::Big();
      snxt = (point[1]+disty)/cn;
      // find extrapolated X and Z
      znew = point[2]+snxt*dir[2];
      if (TMath::Abs(znew) < fDz) {
         Double_t dx = 0.5*(fDx1+fDx2)-fx*znew;
         xnew = point[0]+snxt*dir[0];
         if (TMath::Abs(xnew) < dx) return snxt;
      }
   }
   if (point[1]>=disty) {
      cn = dir[1]+fy*dir[2];
      if (cn>=0) return TGeoShape::Big();
      in = kFALSE;
      snxt = (disty-point[1])/cn;
      // find extrapolated X and Z
      znew = point[2]+snxt*dir[2];
      if (TMath::Abs(znew) < fDz) {
         Double_t dx = 0.5*(fDx1+fDx2)-fx*znew;
         xnew = point[0]+snxt*dir[0];
         if (TMath::Abs(xnew) < dx) return snxt;
      }
   }
   if (!in) return TGeoShape::Big();
   // Point actually inside
   if (safz<safx && safz<safy) {
      if (point[2]*dir[2]>=0) return TGeoShape::Big();
      return 0.0;
   }
   if (safy<safx) {
      cn = TMath::Sign(1.0,point[1])*dir[1]+fy*dir[2];
      if (cn>=0) return TGeoShape::Big();
      return 0.0;
   }
   cn = TMath::Sign(1.0,point[0])*dir[0]+fx*dir[2];
   if (cn>=0) return TGeoShape::Big();
   return 0.0;
}

//_____________________________________________________________________________
Double_t TGeoTrd2::GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
{
// Get range of shape for a given axis.
   xlo = 0;
   xhi = 0;
   Double_t dx = 0;
   switch (iaxis) {
      case 3:
         xlo = -fDz;
         xhi = fDz;
         dx = xhi-xlo;
         return dx;
   }
   return dx;
}

//_____________________________________________________________________________
void TGeoTrd2::GetVisibleCorner(const Double_t *point, Double_t *vertex, Double_t *normals) const
{
// get the most visible corner from outside point and the normals
   Double_t fx = 0.5*(fDx1-fDx2)/fDz;
   Double_t fy = 0.5*(fDy1-fDy2)/fDz;
   Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
   Double_t salf = calf*fx;
   Double_t cbet = 1./TMath::Sqrt(1.0+fy*fy);
   Double_t sbet = cbet*fy;
   // check visibility of X,Y faces
   Double_t distx = fDx1-fx*(fDz+point[2]);
   Double_t disty = fDy1-fy*(fDz+point[2]);
   memset(normals, 0, 9*sizeof(Double_t));
   TGeoTrd2 *trd2 = (TGeoTrd2*)this;
   if (point[0]>distx) {
   // hi x face visible
      trd2->SetShapeBit(kGeoVisX);
      normals[0]=calf;
      normals[2]=salf;
   } else {
      trd2->SetShapeBit(kGeoVisX, kFALSE);
      normals[0]=-calf;
      normals[2]=salf;
   }
   if (point[1]>disty) {
   // hi y face visible
      trd2->SetShapeBit(kGeoVisY);
      normals[4]=cbet;
      normals[5]=sbet;
   } else {
      trd2->SetShapeBit(kGeoVisY, kFALSE);
      normals[4]=-cbet;
      normals[5]=sbet;
   }
   if (point[2]>fDz) {
   // hi z face visible
      trd2->SetShapeBit(kGeoVisZ);
      normals[8]=1;
   } else {
      trd2->SetShapeBit(kGeoVisZ, kFALSE);
      normals[8]=-1;
   }
   SetVertex(vertex);
}

//_____________________________________________________________________________
void TGeoTrd2::GetOppositeCorner(const Double_t * /*point*/, Int_t inorm, Double_t *vertex, Double_t *normals) const
{
// get the opposite corner of the intersected face
   TGeoTrd2 *trd2 = (TGeoTrd2*)this;
   if (inorm != 0) {
   // change x face
      trd2->SetShapeBit(kGeoVisX, !TestShapeBit(kGeoVisX));
      normals[0]=-normals[0];
   }
   if (inorm != 1) {
   // change y face
      trd2->SetShapeBit(kGeoVisY, !TestShapeBit(kGeoVisY));
      normals[4]=-normals[4];
   }
   if (inorm != 2) {
   // hi z face visible
      trd2->SetShapeBit(kGeoVisZ, !TestShapeBit(kGeoVisZ));
      normals[8]=-normals[8];
   }
   SetVertex(vertex);
}

//_____________________________________________________________________________
TGeoVolume *TGeoTrd2::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. Only Z divisions
// are supported. For Z divisions just return the pointer to the volume to be
// divided. In case a wrong division axis is supplied, returns pointer to
// volume that was divided.
   TGeoShape *shape;           //--- shape to be created
   TGeoVolume *vol;            //--- division volume to be created
   TGeoVolumeMulti *vmulti;    //--- generic divided volume
   TGeoPatternFinder *finder;  //--- finder to be attached
   TString opt = "";           //--- option to be attached
   Double_t zmin, zmax, dx1n, dx2n, dy1n, dy2n;
   Int_t id;
   Double_t end = start+ndiv*step;
   switch (iaxis) {
      case 1:
         Warning("Divide", "dividing a Trd2 on X not implemented");
         return 0;
      case 2:
         Warning("Divide", "dividing a Trd2 on Y not implemented");
         return 0;
      case 3:
         finder = new TGeoPatternZ(voldiv, ndiv, start, end);
         vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
         voldiv->SetFinder(finder);
         finder->SetDivIndex(voldiv->GetNdaughters());
         for (id=0; id<ndiv; id++) {
            zmin = start+id*step;
            zmax = start+(id+1)*step;
            dx1n = 0.5*(fDx1*(fDz-zmin)+fDx2*(fDz+zmin))/fDz;
            dx2n = 0.5*(fDx1*(fDz-zmax)+fDx2*(fDz+zmax))/fDz;
            dy1n = 0.5*(fDy1*(fDz-zmin)+fDy2*(fDz+zmin))/fDz;
            dy2n = 0.5*(fDy1*(fDz-zmax)+fDy2*(fDz+zmax))/fDz;
            shape = new TGeoTrd2(dx1n, dx2n, dy1n, dy2n, step/2.);
            vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
            vmulti->AddVolume(vol);
            opt = "Z";
            voldiv->AddNodeOffset(vol, id, start+step/2+id*step, opt.Data());
            ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
         }
         return vmulti;
      default:
         Error("Divide", "Wrong axis type for division");
         return 0;
   }
}

//_____________________________________________________________________________
void TGeoTrd2::GetBoundingCylinder(Double_t *param) const
{
//--- Fill vector param[4] with the bounding cylinder parameters. The order
// is the following : Rmin, Rmax, Phi1, Phi2
   TGeoBBox::GetBoundingCylinder(param);
}

//_____________________________________________________________________________
Int_t TGeoTrd2::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.
   dx=dy=dz=0;
   if (mat->IsRotation()) {
      Error("GetFittingBox", "cannot handle parametrized rotated volumes");
      return 1; // ### rotation not accepted ###
   }
   //--> translate the origin of the parametrized box to the frame of this box.
   Double_t origin[3];
   mat->LocalToMaster(parambox->GetOrigin(), origin);
   if (!Contains(origin)) {
      Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
      return 1; // ### wrong matrix ###
   }
   //--> now we have to get the valid range for all parametrized axis
   Double_t dd[3];
   dd[0] = parambox->GetDX();
   dd[1] = parambox->GetDY();
   dd[2] = parambox->GetDZ();
   //-> check if Z range is fixed
   if (dd[2]<0) {
      dd[2] = TMath::Min(origin[2]+fDz, fDz-origin[2]);
      if (dd[2]<0) {
         Error("GetFittingBox", "wrong matrix");
         return 1;
      }
   }
   if (dd[0]>=0 && dd[1]>=0) {
      dx = dd[0];
      dy = dd[1];
      dz = dd[2];
      return 0;
   }
   //-> check now range at Z = origin[2] +/- dd[2]
   Double_t fx = 0.5*(fDx1-fDx2)/fDz;
   Double_t fy = 0.5*(fDy1-fDy2)/fDz;
   Double_t dx0 = 0.5*(fDx1+fDx2);
   Double_t dy0 = 0.5*(fDy1+fDy2);
   Double_t z=origin[2]-dd[2];
   dd[0] = dx0-fx*z-origin[0];
   dd[1] = dy0-fy*z-origin[1];
   z=origin[2]+dd[2];
   dd[0] = TMath::Min(dd[0], dx0-fx*z-origin[0]);
   dd[1] = TMath::Min(dd[1], dy0-fy*z-origin[1]);
   if (dd[0]<0 || dd[1]<0) {
      Error("GetFittingBox", "wrong matrix");
      return 1;
   }
   dx = dd[0];
   dy = dd[1];
   dz = dd[2];
   return 0;
}

//_____________________________________________________________________________
TGeoShape *TGeoTrd2::GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix * /*mat*/) const
{
// in case shape has some negative parameters, these has to be computed
// in order to fit the mother
   if (!TestShapeBit(kGeoRunTimeShape)) return 0;
   if (!mother->TestShapeBit(kGeoTrd2)) {
      Error("GetMakeRuntimeShape", "invalid mother");
      return 0;
   }
   Double_t dx1, dx2, dy1, dy2, dz;
   if (fDx1<0) dx1=((TGeoTrd2*)mother)->GetDx1();
   else dx1=fDx1;
   if (fDx2<0) dx2=((TGeoTrd2*)mother)->GetDx2();
   else dx2=fDx2;
   if (fDy1<0) dy1=((TGeoTrd2*)mother)->GetDy1();
   else dy1=fDy1;
   if (fDy2<0) dy2=((TGeoTrd2*)mother)->GetDy2();
   else dy2=fDy2;
   if (fDz<0) dz=((TGeoTrd2*)mother)->GetDz();
   else dz=fDz;

   return (new TGeoTrd2(dx1, dx2, dy1, dy2, dz));
}

//_____________________________________________________________________________
void TGeoTrd2::InspectShape() const
{
// print shape parameters
   printf("*** Shape %s: TGeoTrd2 ***\n", GetName());
   printf("    dx1 = %11.5f\n", fDx1);
   printf("    dx2 = %11.5f\n", fDx2);
   printf("    dy1 = %11.5f\n", fDy1);
   printf("    dy2 = %11.5f\n", fDy2);
   printf("    dz  = %11.5f\n", fDz);
   printf(" Bounding box:\n");
   TGeoBBox::InspectShape();
}

//_____________________________________________________________________________
Double_t TGeoTrd2::Safety(const Double_t *point, Bool_t in) const
{
// computes the closest distance from given point to this shape, according
// to option. The matching point on the shape is stored in spoint.
   Double_t saf[3];
   //--- Compute safety first
   // check Z facettes
   saf[0] = fDz-TMath::Abs(point[2]);
   Double_t fx = 0.5*(fDx1-fDx2)/fDz;
   Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
   // check X facettes
   Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
   if (distx<0) saf[1]=TGeoShape::Big();
   else         saf[1]=(distx-TMath::Abs(point[0]))*calf;

   Double_t fy = 0.5*(fDy1-fDy2)/fDz;
   calf = 1./TMath::Sqrt(1.0+fy*fy);
   // check Y facettes
   distx = 0.5*(fDy1+fDy2)-fy*point[2];
   if (distx<0) saf[2]=TGeoShape::Big();
   else         saf[2]=(distx-TMath::Abs(point[1]))*calf;

   if (in) return saf[TMath::LocMin(3,saf)];
   for (Int_t i=0; i<3; i++) saf[i]=-saf[i];
   return saf[TMath::LocMax(3,saf)];
}

//_____________________________________________________________________________
void TGeoTrd2::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
{
// Save a primitive as a C++ statement(s) on output stream "out".
   if (TObject::TestBit(kGeoSavePrimitive)) return;
   out << "   // Shape: " << GetName() << " type: " << ClassName() << std::endl;
   out << "   dx1 = " << fDx1 << ";" << std::endl;
   out << "   dx2 = " << fDx2 << ";" << std::endl;
   out << "   dy1 = " << fDy1 << ";" << std::endl;
   out << "   dy2 = " << fDy2 << ";" << std::endl;
   out << "   dz  = " << fDZ  << ";" << std::endl;
   out << "   TGeoShape *" << GetPointerName() << " = new TGeoTrd2(\"" << GetName() << "\", dx1,dx2,dy1,dy2,dz);" << std::endl;
   TObject::SetBit(TGeoShape::kGeoSavePrimitive);
}

//_____________________________________________________________________________
void TGeoTrd2::SetDimensions(Double_t *param)
{
// set arb8 params in one step :
   fDx1 = param[0];
   fDx2 = param[1];
   fDy1 = param[2];
   fDy2 = param[3];
   fDz  = param[4];
   ComputeBBox();
}

//_____________________________________________________________________________
void TGeoTrd2::SetPoints(Double_t *points) const
{
// create trd2 mesh points
   if (!points) return;
   points[ 0] = -fDx1; points[ 1] = -fDy1; points[ 2] = -fDz;
   points[ 3] = -fDx1; points[ 4] =  fDy1; points[ 5] = -fDz;
   points[ 6] =  fDx1; points[ 7] =  fDy1; points[ 8] = -fDz;
   points[ 9] =  fDx1; points[10] = -fDy1; points[11] = -fDz;
   points[12] = -fDx2; points[13] = -fDy2; points[14] =  fDz;
   points[15] = -fDx2; points[16] =  fDy2; points[17] =  fDz;
   points[18] =  fDx2; points[19] =  fDy2; points[20] =  fDz;
   points[21] =  fDx2; points[22] = -fDy2; points[23] =  fDz;
}

//_____________________________________________________________________________
void TGeoTrd2::SetPoints(Float_t *points) const
{
// create trd2 mesh points
   if (!points) return;
   points[ 0] = -fDx1; points[ 1] = -fDy1; points[ 2] = -fDz;
   points[ 3] = -fDx1; points[ 4] =  fDy1; points[ 5] = -fDz;
   points[ 6] =  fDx1; points[ 7] =  fDy1; points[ 8] = -fDz;
   points[ 9] =  fDx1; points[10] = -fDy1; points[11] = -fDz;
   points[12] = -fDx2; points[13] = -fDy2; points[14] =  fDz;
   points[15] = -fDx2; points[16] =  fDy2; points[17] =  fDz;
   points[18] =  fDx2; points[19] =  fDy2; points[20] =  fDz;
   points[21] =  fDx2; points[22] = -fDy2; points[23] =  fDz;
}

//_____________________________________________________________________________
void TGeoTrd2::SetVertex(Double_t *vertex) const
{
// set vertex of a corner according to visibility flags
   if (TestShapeBit(kGeoVisX)) {
      if (TestShapeBit(kGeoVisZ)) {
         vertex[0] = fDx2;
         vertex[2] = fDz;
         vertex[1] = (TestShapeBit(kGeoVisY))?fDy2:-fDy2;
      } else {
         vertex[0] = fDx1;
         vertex[2] = -fDz;
         vertex[1] = (TestShapeBit(kGeoVisY))?fDy1:-fDy1;
      }
   } else {
      if (TestShapeBit(kGeoVisZ)) {
         vertex[0] = -fDx2;
         vertex[2] = fDz;
         vertex[1] = (TestShapeBit(kGeoVisY))?fDy2:-fDy2;
      } else {
         vertex[0] = -fDx1;
         vertex[2] = -fDz;
         vertex[1] = (TestShapeBit(kGeoVisY))?fDy1:-fDy1;
      }
   }
}

//_____________________________________________________________________________
void TGeoTrd2::Sizeof3D() const
{
// fill size of this 3-D object
   TGeoBBox::Sizeof3D();
}

//_____________________________________________________________________________
void TGeoTrd2::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.
// Input: Array of point coordinates + vector size
// Output: Array of Booleans for the inside of each point
   for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
}

//_____________________________________________________________________________
void TGeoTrd2::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 coordinates and directions + vector size
// Output: Array of normal directions
   for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
}

//_____________________________________________________________________________
void TGeoTrd2::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 dists
   for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
}

//_____________________________________________________________________________
void TGeoTrd2::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 dists
   for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
}

//_____________________________________________________________________________
void TGeoTrd2::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.
// Input: Array of point coordinates, array of statuses for these points, size of the arrays
// Output: Safety values
   for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
}
 TGeoTrd2.cxx:1
 TGeoTrd2.cxx:2
 TGeoTrd2.cxx:3
 TGeoTrd2.cxx:4
 TGeoTrd2.cxx:5
 TGeoTrd2.cxx:6
 TGeoTrd2.cxx:7
 TGeoTrd2.cxx:8
 TGeoTrd2.cxx:9
 TGeoTrd2.cxx:10
 TGeoTrd2.cxx:11
 TGeoTrd2.cxx:12
 TGeoTrd2.cxx:13
 TGeoTrd2.cxx:14
 TGeoTrd2.cxx:15
 TGeoTrd2.cxx:16
 TGeoTrd2.cxx:17
 TGeoTrd2.cxx:18
 TGeoTrd2.cxx:19
 TGeoTrd2.cxx:20
 TGeoTrd2.cxx:21
 TGeoTrd2.cxx:22
 TGeoTrd2.cxx:23
 TGeoTrd2.cxx:24
 TGeoTrd2.cxx:25
 TGeoTrd2.cxx:26
 TGeoTrd2.cxx:27
 TGeoTrd2.cxx:28
 TGeoTrd2.cxx:29
 TGeoTrd2.cxx:30
 TGeoTrd2.cxx:31
 TGeoTrd2.cxx:32
 TGeoTrd2.cxx:33
 TGeoTrd2.cxx:34
 TGeoTrd2.cxx:35
 TGeoTrd2.cxx:36
 TGeoTrd2.cxx:37
 TGeoTrd2.cxx:38
 TGeoTrd2.cxx:39
 TGeoTrd2.cxx:40
 TGeoTrd2.cxx:41
 TGeoTrd2.cxx:42
 TGeoTrd2.cxx:43
 TGeoTrd2.cxx:44
 TGeoTrd2.cxx:45
 TGeoTrd2.cxx:46
 TGeoTrd2.cxx:47
 TGeoTrd2.cxx:48
 TGeoTrd2.cxx:49
 TGeoTrd2.cxx:50
 TGeoTrd2.cxx:51
 TGeoTrd2.cxx:52
 TGeoTrd2.cxx:53
 TGeoTrd2.cxx:54
 TGeoTrd2.cxx:55
 TGeoTrd2.cxx:56
 TGeoTrd2.cxx:57
 TGeoTrd2.cxx:58
 TGeoTrd2.cxx:59
 TGeoTrd2.cxx:60
 TGeoTrd2.cxx:61
 TGeoTrd2.cxx:62
 TGeoTrd2.cxx:63
 TGeoTrd2.cxx:64
 TGeoTrd2.cxx:65
 TGeoTrd2.cxx:66
 TGeoTrd2.cxx:67
 TGeoTrd2.cxx:68
 TGeoTrd2.cxx:69
 TGeoTrd2.cxx:70
 TGeoTrd2.cxx:71
 TGeoTrd2.cxx:72
 TGeoTrd2.cxx:73
 TGeoTrd2.cxx:74
 TGeoTrd2.cxx:75
 TGeoTrd2.cxx:76
 TGeoTrd2.cxx:77
 TGeoTrd2.cxx:78
 TGeoTrd2.cxx:79
 TGeoTrd2.cxx:80
 TGeoTrd2.cxx:81
 TGeoTrd2.cxx:82
 TGeoTrd2.cxx:83
 TGeoTrd2.cxx:84
 TGeoTrd2.cxx:85
 TGeoTrd2.cxx:86
 TGeoTrd2.cxx:87
 TGeoTrd2.cxx:88
 TGeoTrd2.cxx:89
 TGeoTrd2.cxx:90
 TGeoTrd2.cxx:91
 TGeoTrd2.cxx:92
 TGeoTrd2.cxx:93
 TGeoTrd2.cxx:94
 TGeoTrd2.cxx:95
 TGeoTrd2.cxx:96
 TGeoTrd2.cxx:97
 TGeoTrd2.cxx:98
 TGeoTrd2.cxx:99
 TGeoTrd2.cxx:100
 TGeoTrd2.cxx:101
 TGeoTrd2.cxx:102
 TGeoTrd2.cxx:103
 TGeoTrd2.cxx:104
 TGeoTrd2.cxx:105
 TGeoTrd2.cxx:106
 TGeoTrd2.cxx:107
 TGeoTrd2.cxx:108
 TGeoTrd2.cxx:109
 TGeoTrd2.cxx:110
 TGeoTrd2.cxx:111
 TGeoTrd2.cxx:112
 TGeoTrd2.cxx:113
 TGeoTrd2.cxx:114
 TGeoTrd2.cxx:115
 TGeoTrd2.cxx:116
 TGeoTrd2.cxx:117
 TGeoTrd2.cxx:118
 TGeoTrd2.cxx:119
 TGeoTrd2.cxx:120
 TGeoTrd2.cxx:121
 TGeoTrd2.cxx:122
 TGeoTrd2.cxx:123
 TGeoTrd2.cxx:124
 TGeoTrd2.cxx:125
 TGeoTrd2.cxx:126
 TGeoTrd2.cxx:127
 TGeoTrd2.cxx:128
 TGeoTrd2.cxx:129
 TGeoTrd2.cxx:130
 TGeoTrd2.cxx:131
 TGeoTrd2.cxx:132
 TGeoTrd2.cxx:133
 TGeoTrd2.cxx:134
 TGeoTrd2.cxx:135
 TGeoTrd2.cxx:136
 TGeoTrd2.cxx:137
 TGeoTrd2.cxx:138
 TGeoTrd2.cxx:139
 TGeoTrd2.cxx:140
 TGeoTrd2.cxx:141
 TGeoTrd2.cxx:142
 TGeoTrd2.cxx:143
 TGeoTrd2.cxx:144
 TGeoTrd2.cxx:145
 TGeoTrd2.cxx:146
 TGeoTrd2.cxx:147
 TGeoTrd2.cxx:148
 TGeoTrd2.cxx:149
 TGeoTrd2.cxx:150
 TGeoTrd2.cxx:151
 TGeoTrd2.cxx:152
 TGeoTrd2.cxx:153
 TGeoTrd2.cxx:154
 TGeoTrd2.cxx:155
 TGeoTrd2.cxx:156
 TGeoTrd2.cxx:157
 TGeoTrd2.cxx:158
 TGeoTrd2.cxx:159
 TGeoTrd2.cxx:160
 TGeoTrd2.cxx:161
 TGeoTrd2.cxx:162
 TGeoTrd2.cxx:163
 TGeoTrd2.cxx:164
 TGeoTrd2.cxx:165
 TGeoTrd2.cxx:166
 TGeoTrd2.cxx:167
 TGeoTrd2.cxx:168
 TGeoTrd2.cxx:169
 TGeoTrd2.cxx:170
 TGeoTrd2.cxx:171
 TGeoTrd2.cxx:172
 TGeoTrd2.cxx:173
 TGeoTrd2.cxx:174
 TGeoTrd2.cxx:175
 TGeoTrd2.cxx:176
 TGeoTrd2.cxx:177
 TGeoTrd2.cxx:178
 TGeoTrd2.cxx:179
 TGeoTrd2.cxx:180
 TGeoTrd2.cxx:181
 TGeoTrd2.cxx:182
 TGeoTrd2.cxx:183
 TGeoTrd2.cxx:184
 TGeoTrd2.cxx:185
 TGeoTrd2.cxx:186
 TGeoTrd2.cxx:187
 TGeoTrd2.cxx:188
 TGeoTrd2.cxx:189
 TGeoTrd2.cxx:190
 TGeoTrd2.cxx:191
 TGeoTrd2.cxx:192
 TGeoTrd2.cxx:193
 TGeoTrd2.cxx:194
 TGeoTrd2.cxx:195
 TGeoTrd2.cxx:196
 TGeoTrd2.cxx:197
 TGeoTrd2.cxx:198
 TGeoTrd2.cxx:199
 TGeoTrd2.cxx:200
 TGeoTrd2.cxx:201
 TGeoTrd2.cxx:202
 TGeoTrd2.cxx:203
 TGeoTrd2.cxx:204
 TGeoTrd2.cxx:205
 TGeoTrd2.cxx:206
 TGeoTrd2.cxx:207
 TGeoTrd2.cxx:208
 TGeoTrd2.cxx:209
 TGeoTrd2.cxx:210
 TGeoTrd2.cxx:211
 TGeoTrd2.cxx:212
 TGeoTrd2.cxx:213
 TGeoTrd2.cxx:214
 TGeoTrd2.cxx:215
 TGeoTrd2.cxx:216
 TGeoTrd2.cxx:217
 TGeoTrd2.cxx:218
 TGeoTrd2.cxx:219
 TGeoTrd2.cxx:220
 TGeoTrd2.cxx:221
 TGeoTrd2.cxx:222
 TGeoTrd2.cxx:223
 TGeoTrd2.cxx:224
 TGeoTrd2.cxx:225
 TGeoTrd2.cxx:226
 TGeoTrd2.cxx:227
 TGeoTrd2.cxx:228
 TGeoTrd2.cxx:229
 TGeoTrd2.cxx:230
 TGeoTrd2.cxx:231
 TGeoTrd2.cxx:232
 TGeoTrd2.cxx:233
 TGeoTrd2.cxx:234
 TGeoTrd2.cxx:235
 TGeoTrd2.cxx:236
 TGeoTrd2.cxx:237
 TGeoTrd2.cxx:238
 TGeoTrd2.cxx:239
 TGeoTrd2.cxx:240
 TGeoTrd2.cxx:241
 TGeoTrd2.cxx:242
 TGeoTrd2.cxx:243
 TGeoTrd2.cxx:244
 TGeoTrd2.cxx:245
 TGeoTrd2.cxx:246
 TGeoTrd2.cxx:247
 TGeoTrd2.cxx:248
 TGeoTrd2.cxx:249
 TGeoTrd2.cxx:250
 TGeoTrd2.cxx:251
 TGeoTrd2.cxx:252
 TGeoTrd2.cxx:253
 TGeoTrd2.cxx:254
 TGeoTrd2.cxx:255
 TGeoTrd2.cxx:256
 TGeoTrd2.cxx:257
 TGeoTrd2.cxx:258
 TGeoTrd2.cxx:259
 TGeoTrd2.cxx:260
 TGeoTrd2.cxx:261
 TGeoTrd2.cxx:262
 TGeoTrd2.cxx:263
 TGeoTrd2.cxx:264
 TGeoTrd2.cxx:265
 TGeoTrd2.cxx:266
 TGeoTrd2.cxx:267
 TGeoTrd2.cxx:268
 TGeoTrd2.cxx:269
 TGeoTrd2.cxx:270
 TGeoTrd2.cxx:271
 TGeoTrd2.cxx:272
 TGeoTrd2.cxx:273
 TGeoTrd2.cxx:274
 TGeoTrd2.cxx:275
 TGeoTrd2.cxx:276
 TGeoTrd2.cxx:277
 TGeoTrd2.cxx:278
 TGeoTrd2.cxx:279
 TGeoTrd2.cxx:280
 TGeoTrd2.cxx:281
 TGeoTrd2.cxx:282
 TGeoTrd2.cxx:283
 TGeoTrd2.cxx:284
 TGeoTrd2.cxx:285
 TGeoTrd2.cxx:286
 TGeoTrd2.cxx:287
 TGeoTrd2.cxx:288
 TGeoTrd2.cxx:289
 TGeoTrd2.cxx:290
 TGeoTrd2.cxx:291
 TGeoTrd2.cxx:292
 TGeoTrd2.cxx:293
 TGeoTrd2.cxx:294
 TGeoTrd2.cxx:295
 TGeoTrd2.cxx:296
 TGeoTrd2.cxx:297
 TGeoTrd2.cxx:298
 TGeoTrd2.cxx:299
 TGeoTrd2.cxx:300
 TGeoTrd2.cxx:301
 TGeoTrd2.cxx:302
 TGeoTrd2.cxx:303
 TGeoTrd2.cxx:304
 TGeoTrd2.cxx:305
 TGeoTrd2.cxx:306
 TGeoTrd2.cxx:307
 TGeoTrd2.cxx:308
 TGeoTrd2.cxx:309
 TGeoTrd2.cxx:310
 TGeoTrd2.cxx:311
 TGeoTrd2.cxx:312
 TGeoTrd2.cxx:313
 TGeoTrd2.cxx:314
 TGeoTrd2.cxx:315
 TGeoTrd2.cxx:316
 TGeoTrd2.cxx:317
 TGeoTrd2.cxx:318
 TGeoTrd2.cxx:319
 TGeoTrd2.cxx:320
 TGeoTrd2.cxx:321
 TGeoTrd2.cxx:322
 TGeoTrd2.cxx:323
 TGeoTrd2.cxx:324
 TGeoTrd2.cxx:325
 TGeoTrd2.cxx:326
 TGeoTrd2.cxx:327
 TGeoTrd2.cxx:328
 TGeoTrd2.cxx:329
 TGeoTrd2.cxx:330
 TGeoTrd2.cxx:331
 TGeoTrd2.cxx:332
 TGeoTrd2.cxx:333
 TGeoTrd2.cxx:334
 TGeoTrd2.cxx:335
 TGeoTrd2.cxx:336
 TGeoTrd2.cxx:337
 TGeoTrd2.cxx:338
 TGeoTrd2.cxx:339
 TGeoTrd2.cxx:340
 TGeoTrd2.cxx:341
 TGeoTrd2.cxx:342
 TGeoTrd2.cxx:343
 TGeoTrd2.cxx:344
 TGeoTrd2.cxx:345
 TGeoTrd2.cxx:346
 TGeoTrd2.cxx:347
 TGeoTrd2.cxx:348
 TGeoTrd2.cxx:349
 TGeoTrd2.cxx:350
 TGeoTrd2.cxx:351
 TGeoTrd2.cxx:352
 TGeoTrd2.cxx:353
 TGeoTrd2.cxx:354
 TGeoTrd2.cxx:355
 TGeoTrd2.cxx:356
 TGeoTrd2.cxx:357
 TGeoTrd2.cxx:358
 TGeoTrd2.cxx:359
 TGeoTrd2.cxx:360
 TGeoTrd2.cxx:361
 TGeoTrd2.cxx:362
 TGeoTrd2.cxx:363
 TGeoTrd2.cxx:364
 TGeoTrd2.cxx:365
 TGeoTrd2.cxx:366
 TGeoTrd2.cxx:367
 TGeoTrd2.cxx:368
 TGeoTrd2.cxx:369
 TGeoTrd2.cxx:370
 TGeoTrd2.cxx:371
 TGeoTrd2.cxx:372
 TGeoTrd2.cxx:373
 TGeoTrd2.cxx:374
 TGeoTrd2.cxx:375
 TGeoTrd2.cxx:376
 TGeoTrd2.cxx:377
 TGeoTrd2.cxx:378
 TGeoTrd2.cxx:379
 TGeoTrd2.cxx:380
 TGeoTrd2.cxx:381
 TGeoTrd2.cxx:382
 TGeoTrd2.cxx:383
 TGeoTrd2.cxx:384
 TGeoTrd2.cxx:385
 TGeoTrd2.cxx:386
 TGeoTrd2.cxx:387
 TGeoTrd2.cxx:388
 TGeoTrd2.cxx:389
 TGeoTrd2.cxx:390
 TGeoTrd2.cxx:391
 TGeoTrd2.cxx:392
 TGeoTrd2.cxx:393
 TGeoTrd2.cxx:394
 TGeoTrd2.cxx:395
 TGeoTrd2.cxx:396
 TGeoTrd2.cxx:397
 TGeoTrd2.cxx:398
 TGeoTrd2.cxx:399
 TGeoTrd2.cxx:400
 TGeoTrd2.cxx:401
 TGeoTrd2.cxx:402
 TGeoTrd2.cxx:403
 TGeoTrd2.cxx:404
 TGeoTrd2.cxx:405
 TGeoTrd2.cxx:406
 TGeoTrd2.cxx:407
 TGeoTrd2.cxx:408
 TGeoTrd2.cxx:409
 TGeoTrd2.cxx:410
 TGeoTrd2.cxx:411
 TGeoTrd2.cxx:412
 TGeoTrd2.cxx:413
 TGeoTrd2.cxx:414
 TGeoTrd2.cxx:415
 TGeoTrd2.cxx:416
 TGeoTrd2.cxx:417
 TGeoTrd2.cxx:418
 TGeoTrd2.cxx:419
 TGeoTrd2.cxx:420
 TGeoTrd2.cxx:421
 TGeoTrd2.cxx:422
 TGeoTrd2.cxx:423
 TGeoTrd2.cxx:424
 TGeoTrd2.cxx:425
 TGeoTrd2.cxx:426
 TGeoTrd2.cxx:427
 TGeoTrd2.cxx:428
 TGeoTrd2.cxx:429
 TGeoTrd2.cxx:430
 TGeoTrd2.cxx:431
 TGeoTrd2.cxx:432
 TGeoTrd2.cxx:433
 TGeoTrd2.cxx:434
 TGeoTrd2.cxx:435
 TGeoTrd2.cxx:436
 TGeoTrd2.cxx:437
 TGeoTrd2.cxx:438
 TGeoTrd2.cxx:439
 TGeoTrd2.cxx:440
 TGeoTrd2.cxx:441
 TGeoTrd2.cxx:442
 TGeoTrd2.cxx:443
 TGeoTrd2.cxx:444
 TGeoTrd2.cxx:445
 TGeoTrd2.cxx:446
 TGeoTrd2.cxx:447
 TGeoTrd2.cxx:448
 TGeoTrd2.cxx:449
 TGeoTrd2.cxx:450
 TGeoTrd2.cxx:451
 TGeoTrd2.cxx:452
 TGeoTrd2.cxx:453
 TGeoTrd2.cxx:454
 TGeoTrd2.cxx:455
 TGeoTrd2.cxx:456
 TGeoTrd2.cxx:457
 TGeoTrd2.cxx:458
 TGeoTrd2.cxx:459
 TGeoTrd2.cxx:460
 TGeoTrd2.cxx:461
 TGeoTrd2.cxx:462
 TGeoTrd2.cxx:463
 TGeoTrd2.cxx:464
 TGeoTrd2.cxx:465
 TGeoTrd2.cxx:466
 TGeoTrd2.cxx:467
 TGeoTrd2.cxx:468
 TGeoTrd2.cxx:469
 TGeoTrd2.cxx:470
 TGeoTrd2.cxx:471
 TGeoTrd2.cxx:472
 TGeoTrd2.cxx:473
 TGeoTrd2.cxx:474
 TGeoTrd2.cxx:475
 TGeoTrd2.cxx:476
 TGeoTrd2.cxx:477
 TGeoTrd2.cxx:478
 TGeoTrd2.cxx:479
 TGeoTrd2.cxx:480
 TGeoTrd2.cxx:481
 TGeoTrd2.cxx:482
 TGeoTrd2.cxx:483
 TGeoTrd2.cxx:484
 TGeoTrd2.cxx:485
 TGeoTrd2.cxx:486
 TGeoTrd2.cxx:487
 TGeoTrd2.cxx:488
 TGeoTrd2.cxx:489
 TGeoTrd2.cxx:490
 TGeoTrd2.cxx:491
 TGeoTrd2.cxx:492
 TGeoTrd2.cxx:493
 TGeoTrd2.cxx:494
 TGeoTrd2.cxx:495
 TGeoTrd2.cxx:496
 TGeoTrd2.cxx:497
 TGeoTrd2.cxx:498
 TGeoTrd2.cxx:499
 TGeoTrd2.cxx:500
 TGeoTrd2.cxx:501
 TGeoTrd2.cxx:502
 TGeoTrd2.cxx:503
 TGeoTrd2.cxx:504
 TGeoTrd2.cxx:505
 TGeoTrd2.cxx:506
 TGeoTrd2.cxx:507
 TGeoTrd2.cxx:508
 TGeoTrd2.cxx:509
 TGeoTrd2.cxx:510
 TGeoTrd2.cxx:511
 TGeoTrd2.cxx:512
 TGeoTrd2.cxx:513
 TGeoTrd2.cxx:514
 TGeoTrd2.cxx:515
 TGeoTrd2.cxx:516
 TGeoTrd2.cxx:517
 TGeoTrd2.cxx:518
 TGeoTrd2.cxx:519
 TGeoTrd2.cxx:520
 TGeoTrd2.cxx:521
 TGeoTrd2.cxx:522
 TGeoTrd2.cxx:523
 TGeoTrd2.cxx:524
 TGeoTrd2.cxx:525
 TGeoTrd2.cxx:526
 TGeoTrd2.cxx:527
 TGeoTrd2.cxx:528
 TGeoTrd2.cxx:529
 TGeoTrd2.cxx:530
 TGeoTrd2.cxx:531
 TGeoTrd2.cxx:532
 TGeoTrd2.cxx:533
 TGeoTrd2.cxx:534
 TGeoTrd2.cxx:535
 TGeoTrd2.cxx:536
 TGeoTrd2.cxx:537
 TGeoTrd2.cxx:538
 TGeoTrd2.cxx:539
 TGeoTrd2.cxx:540
 TGeoTrd2.cxx:541
 TGeoTrd2.cxx:542
 TGeoTrd2.cxx:543
 TGeoTrd2.cxx:544
 TGeoTrd2.cxx:545
 TGeoTrd2.cxx:546
 TGeoTrd2.cxx:547
 TGeoTrd2.cxx:548
 TGeoTrd2.cxx:549
 TGeoTrd2.cxx:550
 TGeoTrd2.cxx:551
 TGeoTrd2.cxx:552
 TGeoTrd2.cxx:553
 TGeoTrd2.cxx:554
 TGeoTrd2.cxx:555
 TGeoTrd2.cxx:556
 TGeoTrd2.cxx:557
 TGeoTrd2.cxx:558
 TGeoTrd2.cxx:559
 TGeoTrd2.cxx:560
 TGeoTrd2.cxx:561
 TGeoTrd2.cxx:562
 TGeoTrd2.cxx:563
 TGeoTrd2.cxx:564
 TGeoTrd2.cxx:565
 TGeoTrd2.cxx:566
 TGeoTrd2.cxx:567
 TGeoTrd2.cxx:568
 TGeoTrd2.cxx:569
 TGeoTrd2.cxx:570
 TGeoTrd2.cxx:571
 TGeoTrd2.cxx:572
 TGeoTrd2.cxx:573
 TGeoTrd2.cxx:574
 TGeoTrd2.cxx:575
 TGeoTrd2.cxx:576
 TGeoTrd2.cxx:577
 TGeoTrd2.cxx:578
 TGeoTrd2.cxx:579
 TGeoTrd2.cxx:580
 TGeoTrd2.cxx:581
 TGeoTrd2.cxx:582
 TGeoTrd2.cxx:583
 TGeoTrd2.cxx:584
 TGeoTrd2.cxx:585
 TGeoTrd2.cxx:586
 TGeoTrd2.cxx:587
 TGeoTrd2.cxx:588
 TGeoTrd2.cxx:589
 TGeoTrd2.cxx:590
 TGeoTrd2.cxx:591
 TGeoTrd2.cxx:592
 TGeoTrd2.cxx:593
 TGeoTrd2.cxx:594
 TGeoTrd2.cxx:595
 TGeoTrd2.cxx:596
 TGeoTrd2.cxx:597
 TGeoTrd2.cxx:598
 TGeoTrd2.cxx:599
 TGeoTrd2.cxx:600
 TGeoTrd2.cxx:601
 TGeoTrd2.cxx:602
 TGeoTrd2.cxx:603
 TGeoTrd2.cxx:604
 TGeoTrd2.cxx:605
 TGeoTrd2.cxx:606
 TGeoTrd2.cxx:607
 TGeoTrd2.cxx:608
 TGeoTrd2.cxx:609
 TGeoTrd2.cxx:610
 TGeoTrd2.cxx:611
 TGeoTrd2.cxx:612
 TGeoTrd2.cxx:613
 TGeoTrd2.cxx:614
 TGeoTrd2.cxx:615
 TGeoTrd2.cxx:616
 TGeoTrd2.cxx:617
 TGeoTrd2.cxx:618
 TGeoTrd2.cxx:619
 TGeoTrd2.cxx:620
 TGeoTrd2.cxx:621
 TGeoTrd2.cxx:622
 TGeoTrd2.cxx:623
 TGeoTrd2.cxx:624
 TGeoTrd2.cxx:625
 TGeoTrd2.cxx:626
 TGeoTrd2.cxx:627
 TGeoTrd2.cxx:628
 TGeoTrd2.cxx:629
 TGeoTrd2.cxx:630
 TGeoTrd2.cxx:631
 TGeoTrd2.cxx:632
 TGeoTrd2.cxx:633
 TGeoTrd2.cxx:634
 TGeoTrd2.cxx:635
 TGeoTrd2.cxx:636
 TGeoTrd2.cxx:637
 TGeoTrd2.cxx:638
 TGeoTrd2.cxx:639
 TGeoTrd2.cxx:640
 TGeoTrd2.cxx:641
 TGeoTrd2.cxx:642
 TGeoTrd2.cxx:643
 TGeoTrd2.cxx:644
 TGeoTrd2.cxx:645
 TGeoTrd2.cxx:646
 TGeoTrd2.cxx:647
 TGeoTrd2.cxx:648
 TGeoTrd2.cxx:649
 TGeoTrd2.cxx:650
 TGeoTrd2.cxx:651
 TGeoTrd2.cxx:652
 TGeoTrd2.cxx:653
 TGeoTrd2.cxx:654
 TGeoTrd2.cxx:655
 TGeoTrd2.cxx:656
 TGeoTrd2.cxx:657
 TGeoTrd2.cxx:658
 TGeoTrd2.cxx:659
 TGeoTrd2.cxx:660
 TGeoTrd2.cxx:661
 TGeoTrd2.cxx:662
 TGeoTrd2.cxx:663
 TGeoTrd2.cxx:664
 TGeoTrd2.cxx:665
 TGeoTrd2.cxx:666
 TGeoTrd2.cxx:667
 TGeoTrd2.cxx:668
 TGeoTrd2.cxx:669
 TGeoTrd2.cxx:670
 TGeoTrd2.cxx:671
 TGeoTrd2.cxx:672
 TGeoTrd2.cxx:673
 TGeoTrd2.cxx:674
 TGeoTrd2.cxx:675
 TGeoTrd2.cxx:676
 TGeoTrd2.cxx:677
 TGeoTrd2.cxx:678
 TGeoTrd2.cxx:679
 TGeoTrd2.cxx:680
 TGeoTrd2.cxx:681
 TGeoTrd2.cxx:682
 TGeoTrd2.cxx:683
 TGeoTrd2.cxx:684
 TGeoTrd2.cxx:685
 TGeoTrd2.cxx:686
 TGeoTrd2.cxx:687
 TGeoTrd2.cxx:688
 TGeoTrd2.cxx:689
 TGeoTrd2.cxx:690
 TGeoTrd2.cxx:691
 TGeoTrd2.cxx:692
 TGeoTrd2.cxx:693
 TGeoTrd2.cxx:694
 TGeoTrd2.cxx:695
 TGeoTrd2.cxx:696
 TGeoTrd2.cxx:697
 TGeoTrd2.cxx:698
 TGeoTrd2.cxx:699
 TGeoTrd2.cxx:700
 TGeoTrd2.cxx:701
 TGeoTrd2.cxx:702
 TGeoTrd2.cxx:703
 TGeoTrd2.cxx:704
 TGeoTrd2.cxx:705
 TGeoTrd2.cxx:706
 TGeoTrd2.cxx:707
 TGeoTrd2.cxx:708
 TGeoTrd2.cxx:709
 TGeoTrd2.cxx:710
 TGeoTrd2.cxx:711
 TGeoTrd2.cxx:712
 TGeoTrd2.cxx:713
 TGeoTrd2.cxx:714
 TGeoTrd2.cxx:715
 TGeoTrd2.cxx:716
 TGeoTrd2.cxx:717
 TGeoTrd2.cxx:718
 TGeoTrd2.cxx:719
 TGeoTrd2.cxx:720
 TGeoTrd2.cxx:721
 TGeoTrd2.cxx:722
 TGeoTrd2.cxx:723
 TGeoTrd2.cxx:724
 TGeoTrd2.cxx:725
 TGeoTrd2.cxx:726
 TGeoTrd2.cxx:727
 TGeoTrd2.cxx:728
 TGeoTrd2.cxx:729
 TGeoTrd2.cxx:730
 TGeoTrd2.cxx:731
 TGeoTrd2.cxx:732
 TGeoTrd2.cxx:733
 TGeoTrd2.cxx:734
 TGeoTrd2.cxx:735
 TGeoTrd2.cxx:736
 TGeoTrd2.cxx:737
 TGeoTrd2.cxx:738
 TGeoTrd2.cxx:739
 TGeoTrd2.cxx:740
 TGeoTrd2.cxx:741
 TGeoTrd2.cxx:742
 TGeoTrd2.cxx:743
 TGeoTrd2.cxx:744
 TGeoTrd2.cxx:745
 TGeoTrd2.cxx:746
 TGeoTrd2.cxx:747
 TGeoTrd2.cxx:748
 TGeoTrd2.cxx:749
 TGeoTrd2.cxx:750
 TGeoTrd2.cxx:751
 TGeoTrd2.cxx:752
 TGeoTrd2.cxx:753
 TGeoTrd2.cxx:754
 TGeoTrd2.cxx:755
 TGeoTrd2.cxx:756
 TGeoTrd2.cxx:757
 TGeoTrd2.cxx:758
 TGeoTrd2.cxx:759
 TGeoTrd2.cxx:760
 TGeoTrd2.cxx:761
 TGeoTrd2.cxx:762
 TGeoTrd2.cxx:763
 TGeoTrd2.cxx:764
 TGeoTrd2.cxx:765
 TGeoTrd2.cxx:766
 TGeoTrd2.cxx:767
 TGeoTrd2.cxx:768
 TGeoTrd2.cxx:769
 TGeoTrd2.cxx:770
 TGeoTrd2.cxx:771
 TGeoTrd2.cxx:772
 TGeoTrd2.cxx:773
 TGeoTrd2.cxx:774
 TGeoTrd2.cxx:775
 TGeoTrd2.cxx:776
 TGeoTrd2.cxx:777