library: libGeom
#include "TGeoTube.h"

TGeoTubeSeg


class description - header file - source file - inheritance tree (.pdf)

class TGeoTubeSeg : public TGeoTube

Inheritance Chart:
TObject
<-
TNamed
<-
TGeoShape
<-
TGeoBBox
<-
TGeoTube
<-
TGeoTubeSeg
<-
TGeoCtub

    public:
TGeoTubeSeg() TGeoTubeSeg(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2) TGeoTubeSeg(const char* name, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2) TGeoTubeSeg(Double_t* params) TGeoTubeSeg(const TGeoTubeSeg&) virtual ~TGeoTubeSeg() virtual Double_t Capacity() const static Double_t Capacity(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2) static TClass* Class() virtual void ComputeBBox() virtual void ComputeNormal(Double_t* point, Double_t* dir, Double_t* norm) static void ComputeNormalS(Double_t* point, Double_t* dir, Double_t* norm, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2) virtual Bool_t Contains(Double_t* point) const virtual Int_t DistancetoPrimitive(Int_t px, Int_t py) virtual Double_t DistFromInside(Double_t* point, Double_t* dir, Int_t iact = 1, Double_t step = TGeoShape::Big(), Double_t* safe = 0) const static Double_t DistFromInsideS(Double_t* point, Double_t* dir, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi) virtual Double_t DistFromOutside(Double_t* point, Double_t* dir, Int_t iact = 1, Double_t step = TGeoShape::Big(), Double_t* safe = 0) const static Double_t DistFromOutsideS(Double_t* point, Double_t* dir, Double_t rmin, Double_t rmax, Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi) virtual TGeoVolume* Divide(TGeoVolume* voldiv, const char* divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) virtual Double_t GetAxisRange(Int_t iaxis, Double_t& xlo, Double_t& xhi) const virtual void GetBoundingCylinder(Double_t* param) const virtual const TBuffer3D& GetBuffer3D(Int_t reqSections, Bool_t localFrame) const virtual Int_t GetByteCount() const virtual TGeoShape* GetMakeRuntimeShape(TGeoShape* mother, TGeoMatrix* mat) const virtual Int_t GetNmeshVertices() const Double_t GetPhi1() const Double_t GetPhi2() const virtual void InspectShape() const virtual TClass* IsA() const virtual TBuffer3D* MakeBuffer3D() const TGeoTubeSeg& operator=(const TGeoTubeSeg&) virtual Double_t Safety(Double_t* point, Bool_t in = kTRUE) const static Double_t SafetyS(Double_t* point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Int_t skipz = 0) virtual void SavePrimitive(ostream& out, Option_t* option = "") virtual void SetDimensions(Double_t* param) virtual void SetPoints(Double_t* points) const virtual void SetPoints(Float_t* points) const virtual void SetSegsAndPols(TBuffer3D& buff) const void SetTubsDimensions(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Sizeof3D() const virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
Double_t fPhi1 first phi limit Double_t fPhi2 second phi limit

Class Description

_____________________________________________________________________________
 TGeoTube - cylindrical tube class. It takes 3 parameters :
            inner radius, outer radius and half-length dz.

_____________________________________________________________________________
























_____________________________________________________________________________
 TGeoTubeSeg - a phi segment of a tube. Has 5 parameters :
            - the same 3 as a tube;
            - first phi limit (in degrees)
            - second phi limit

_____________________________________________________________________________












_____________________________________________________________________________
 TGeoCtub - a tube segment cut with 2 planes. Has 11 parameters :
            - the same 5 as a tube segment;
            - x, y, z components of the normal to the -dZ cut plane in
              point (0, 0, -dZ);
            - x, y, z components of the normal to the +dZ cut plane in
              point (0, 0, dZ);

_____________________________________________________________________________
TGeoTubeSeg()
 Default constructor
TGeoTubeSeg(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
 Default constructor specifying minimum and maximum radius
TGeoTubeSeg(const char *name, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
 Default constructor specifying minimum and maximum radius
TGeoTubeSeg(Double_t *param)
 Default constructor specifying minimum and maximum radius
 param[0] = Rmin
 param[1] = Rmax
 param[2] = dz
 param[3] = phi1
 param[4] = phi2
~TGeoTubeSeg()
 destructor
Double_t Capacity()
 Computes capacity of the shape in [length^3]
Double_t Capacity(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
 Computes capacity of the shape in [length^3]
void ComputeBBox()
 compute bounding box of the tube segment
void ComputeNormal(Double_t *point, Double_t *dir, Double_t *norm)
 Compute normal to closest surface from POINT.
void ComputeNormalS(Double_t *point, Double_t *dir, Double_t *norm, Double_t rmin, Double_t rmax, Double_t /*dz*/, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
{
 Compute normal to closest surface from POINT.
Double_t saf[2];
Double_t rsq = point[0]*point[0]+point[1]*point[1];
Double_t r = TMath::Sqrt(rsq);
saf[0] = (rmin>1E-10)?TMath::Abs(r-rmin):TGeoShape::Big();
saf[1] = TMath::Abs(rmax-r);
Int_t i = TMath::LocMin(2,saf);
if (TGeoShape::IsCloseToPhi(saf[i], point,c1,s1,c2,s2)) {
TGeoShape::NormalPhi(point,dir,norm,c1,s1,c2,s2);
return;
}
norm[2] = 0;
Double_t phi = TMath::ATan2(point[1], point[0]);
norm[0] = TMath::Cos(phi);
norm[1] = TMath::Sin(phi);
if (norm[0]*dir[0]+norm[1]*dir[1]<0) {
norm[0] = -norm[0];
norm[1] = -norm[1];
}
}

_____________________________________________________________________________
Bool_t TGeoTubeSeg::Contains(Double_t *point) const
{
 test if point is inside this tube segment
 first check if point is inside the tube
if (!TGeoTube::Contains(point)) return kFALSE;
return IsInPhiRange(point, fPhi1, fPhi2);
}

_____________________________________________________________________________
Int_t TGeoTubeSeg::DistancetoPrimitive(Int_t px, Int_t py)
{
 compute closest distance from point px,py to each corner
Int_t n = gGeoManager->GetNsegments()+1;
const Int_t numPoints = 4*n;
return ShapeDistancetoPrimitive(numPoints, px, py);
}

_____________________________________________________________________________
Double_t TGeoTubeSeg::DistFromInsideS(Double_t *point, Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz,
Double_t c1, Double_t s1, Double_t c2, Double_t s2, Double_t cm, Double_t sm, Double_t cdfi)
{
 Compute distance from inside point to surface of the tube segment (static)
 Boundary safe algorithm.
 Do Z
Double_t stube = TGeoTube::DistFromInsideS(point,dir,rmin,rmax,dz);
if (stube<=0) return 0.0;
Double_t sfmin = TGeoShape::Big();
Double_t rsq = point[0]*point[0]+point[1]*point[1];
Double_t r = TMath::Sqrt(rsq);
Double_t cpsi=point[0]*cm+point[1]*sm;
if (cpsi>r*cdfi+TGeoShape::Tolerance())  {
sfmin = TGeoShape::DistToPhiMin(point, dir, s1, c1, s2, c2, sm, cm);
return TMath::Min(stube,sfmin);
}
 Point on the phi boundary or outside
 which one: phi1 or phi2
Double_t ddotn, xi, yi;
if (TMath::Abs(point[1]-s1*r) < TMath::Abs(point[1]-s2*r)) {
ddotn = s1*dir[0]-c1*dir[1];
if (ddotn>=0) return 0.0;
ddotn = -s2*dir[0]+c2*dir[1];
if (ddotn<=0) return stube;
sfmin = s2*point[0]-c2*point[1];
if (sfmin<=0) return stube;
sfmin /= ddotn;
if (sfmin >= stube) return stube;
xi = point[0]+sfmin*dir[0];
yi = point[1]+sfmin*dir[1];
if (yi*cm-xi*sm<0) return stube;
return sfmin;
}
ddotn = -s2*dir[0]+c2*dir[1];
if (ddotn>=0) return 0.0;
ddotn = s1*dir[0]-c1*dir[1];
if (ddotn<=0) return stube;
sfmin = -s1*point[0]+c1*point[1];
if (sfmin<=0) return stube;
sfmin /= ddotn;
if (sfmin >= stube) return stube;
xi = point[0]+sfmin*dir[0];
yi = point[1]+sfmin*dir[1];
if (yi*cm-xi*sm>0) return stube;
return sfmin;
}

_____________________________________________________________________________
Double_t TGeoTubeSeg::DistFromInside(Double_t *point, Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
{
 Compute distance from inside point to surface of the tube segment
 Boundary safe algorithm.
if (iact<3 && safe) {
*safe = SafetyS(point, kTRUE, fRmin, fRmax, fDz, fPhi1, fPhi2);
if (iact==0) return TGeoShape::Big();
if ((iact==1) && (*safe>step)) return TGeoShape::Big();
}
Double_t phi1 = fPhi1*TMath::DegToRad();
Double_t phi2 = fPhi2*TMath::DegToRad();
Double_t c1 = TMath::Cos(phi1);
Double_t c2 = TMath::Cos(phi2);
Double_t s1 = TMath::Sin(phi1);
Double_t s2 = TMath::Sin(phi2);
Double_t phim = 0.5*(phi1+phi2);
Double_t cm = TMath::Cos(phim);
Double_t sm = TMath::Sin(phim);
Double_t dfi = 0.5*(phi2-phi1);
Double_t cdfi = TMath::Cos(dfi);

 compute distance to surface
return TGeoTubeSeg::DistFromInsideS(point,dir,fRmin,fRmax,fDz,c1,s1,c2,s2,cm,sm,cdfi);
}

_____________________________________________________________________________
Double_t TGeoTubeSeg::DistFromOutsideS(Double_t *point, Double_t *dir, Double_t rmin, Double_t rmax,
Double_t dz, Double_t c1, Double_t s1, Double_t c2, Double_t s2,
Double_t cm, Double_t sm, Double_t cdfi)
{
 Static method to compute distance to arbitrary tube segment from outside point
 Boundary safe algorithm.
Double_t r2, cpsi;
 check Z planes
Double_t xi, yi, zi;
zi = dz - TMath::Abs(point[2]);
Double_t rmaxsq = rmax*rmax;
Double_t rminsq = rmin*rmin;
Double_t s = TGeoShape::Big();
Double_t snxt=TGeoShape::Big();
Bool_t in = kFALSE;
Bool_t inz = (zi<0)?kFALSE:kTRUE;
if (!inz) {
if (point[2]*dir[2]>=0) return TGeoShape::Big();
s = -zi/TMath::Abs(dir[2]);
xi = point[0]+s*dir[0];
yi = point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
cpsi=(xi*cm+yi*sm)/TMath::Sqrt(r2);
if (cpsi>=cdfi) return s;
}
}

 check outer cyl. surface
Double_t rsq = point[0]*point[0]+point[1]*point[1];
Double_t r = TMath::Sqrt(rsq);
Double_t nsq=dir[0]*dir[0]+dir[1]*dir[1];
Double_t rdotn=point[0]*dir[0]+point[1]*dir[1];
Double_t b,d;
Bool_t inrmax = kFALSE;
Bool_t inrmin = kFALSE;
Bool_t inphi  = kFALSE;
if (rsq<=rmaxsq+TGeoShape::Tolerance()) inrmax = kTRUE;
if (rsq>=rminsq-TGeoShape::Tolerance()) inrmin = kTRUE;
cpsi=point[0]*cm+point[1]*sm;
if (cpsi>r*cdfi-TGeoShape::Tolerance())  inphi = kTRUE;
in = inz & inrmin & inrmax & inphi;
 If inside, we are most likely on a boundary within machine precision.
if (in) {
Bool_t checkout = kFALSE;
Double_t safphi = (cpsi-r*cdfi)*TMath::Sqrt(1.-cdfi*cdfi);
      Double_t sch, cch;
 check if on Z boundaries
if (zi<rmax-r) {
if ((rmin==0) || (zi<r-rmin)) {
if (zi<safphi) {
if (point[2]*dir[2]<0) return 0.0;
return TGeoShape::Big();
}
}
}
if ((rmaxsq-rsq) < (rsq-rminsq)) checkout = kTRUE;
 check if on Rmax boundary
if (checkout && (rmax-r<safphi)) {
if (rdotn>=0) return TGeoShape::Big();
return 0.0;
}
if (TMath::Abs(nsq)<TGeoShape::Tolerance()) return TGeoShape::Big();
 check if on phi boundary
if ((rmin==0) || (safphi<r-rmin)) {
 We may cross again a phi of rmin boundary
 check first if we are on phi1 or phi2
Double_t un;
if (TMath::Abs(point[1]-s1*r) < TMath::Abs(point[1]-s2*r)) {
un = dir[0]*s1-dir[1]*c1;
if (un < 0) return 0.0;
if (cdfi>=0) return TGeoShape::Big();
un = -dir[0]*s2+dir[1]*c2;
if (un<0) {
s = -point[0]*s2+point[1]*c2;
if (s>0) {
s /= (-un);
zi = point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi = point[0]+s*dir[0];
yi = point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
if ((yi*cm-xi*sm)>0) return s;
}
}
}
}
} else {
un = -dir[0]*s2+dir[1]*c2;
if (un < 0) return 0.0;
if (cdfi>=0) return TGeoShape::Big();
un = dir[0]*s1-dir[1]*c1;
if (un<0) {
s = point[0]*s1-point[1]*c1;
if (s>0) {
s /= (-un);
zi = point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi = point[0]+s*dir[0];
yi = point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
if ((yi*cm-xi*sm)<0) return s;
}
}
}
}
}
 We may also cross rmin, (+) solution
if (rdotn>=0) return TGeoShape::Big();
if (cdfi>=0) return TGeoShape::Big();
DistToTube(rsq, nsq, rdotn, rmin, b, d);
if (d>0) {
s=-b+d;
if (s>0) {
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
if ((xi*cm+yi*sm) >= rmin*cdfi) return s;
}
}
}
return TGeoShape::Big();
}
 we are on rmin boundary: we may cross again rmin or a phi facette
if (rdotn>=0) return 0.0;
DistToTube(rsq, nsq, rdotn, rmin, b, d);
if (d>0) {
s=-b+d;
if (s>0) {
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
 now check phi range
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((xi*cm+yi*sm) >= rmin*cdfi) return s;
 now we really have to check any phi crossing
Double_t un=-dir[0]*s1+dir[1]*c1;
if (un > 0) {
s=point[0]*s1-point[1]*c1;
if (s>=0) {
s /= un;
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
if ((yi*cm-xi*sm)<=0) {
if (s<snxt) snxt=s;
}
}
}
}
}
un=dir[0]*s2-dir[1]*c2;
if (un > 0) {
s=(point[1]*c2-point[0]*s2)/un;
if (s>=0 && s<snxt) {
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
if ((yi*cm-xi*sm)>=0) {
return s;
}
}
}
}
}
return snxt;
}
}
}
return TGeoShape::Big();
}
 only r>rmax has to be considered
if (TMath::Abs(nsq)<TGeoShape::Tolerance()) return TGeoShape::Big();
if (rsq>=rmax*rmax) {
if (rdotn>=0) return TGeoShape::Big();
TGeoTube::DistToTube(rsq, nsq, rdotn, rmax, b, d);
if (d>0) {
s=-b-d;
if (s>0) {
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
cpsi = xi*cm+yi*sm;
if (cpsi>=rmax*cdfi) return s;
}
}
}
}
 check inner cylinder
if (rmin>0) {
TGeoTube::DistToTube(rsq, nsq, rdotn, rmin, b, d);
if (d>0) {
s=-b+d;
if (s>0) {
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
cpsi = xi*cm+yi*sm;
if (cpsi>=rmin*cdfi) snxt=s;
}
}
}
}
 check phi planes
Double_t un=-dir[0]*s1+dir[1]*c1;
if (un > 0) {
s=point[0]*s1-point[1]*c1;
if (s>=0) {
s /= un;
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
if ((yi*cm-xi*sm)<=0) {
if (s<snxt) snxt=s;
}
}
}
}
}
un=dir[0]*s2-dir[1]*c2;
if (un > 0) {
s=point[1]*c2-point[0]*s2;
if (s>=0) {
s /= un;
zi=point[2]+s*dir[2];
if (TMath::Abs(zi)<=dz) {
xi=point[0]+s*dir[0];
yi=point[1]+s*dir[1];
r2=xi*xi+yi*yi;
if ((rminsq<=r2) && (r2<=rmaxsq)) {
if ((yi*cm-xi*sm)>=0) {
if (s<snxt) snxt=s;
}
}
}
}
}
return snxt;
}

_____________________________________________________________________________
Double_t TGeoTubeSeg::DistFromOutside(Double_t *point, Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
{
 compute distance from outside point to surface of the tube segment
 fist localize point w.r.t tube
if (iact<3 && safe) {
*safe = SafetyS(point, kFALSE, fRmin, fRmax, fDz, fPhi1, fPhi2);
if (iact==0) return TGeoShape::Big();
if ((iact==1) && (step<=*safe)) return TGeoShape::Big();
}
Double_t phi1 = fPhi1*TMath::DegToRad();
Double_t phi2 = fPhi2*TMath::DegToRad();
Double_t c1 = TMath::Cos(phi1);
Double_t s1 = TMath::Sin(phi1);
Double_t c2 = TMath::Cos(phi2);
Double_t s2 = TMath::Sin(phi2);
Double_t fio = 0.5*(phi1+phi2);
Double_t cm = TMath::Cos(fio);
Double_t sm = TMath::Sin(fio);
Double_t dfi = 0.5*(phi2-phi1);
Double_t cdfi = TMath::Cos(dfi);

 find distance to shape
return TGeoTubeSeg::DistFromOutsideS(point, dir, fRmin, fRmax, fDz, c1, s1, c2, s2, cm, sm, cdfi);
}

_____________________________________________________________________________
TGeoVolume *TGeoTubeSeg::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
Double_t start, Double_t step)
{
--- Divide this tube segment shape belonging to volume "voldiv" into ndiv volumes
 called divname, from start position with the given step. Returns pointer
 to created division cell volume in case of Z divisions. For radialdivision
 creates all volumes with different shapes and returns pointer to volume that
 was divided. In case a wrong division axis is supplied, returns pointer to
 volume that was divided.
Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi)
 Get range of shape for a given axis.
void GetBoundingCylinder(Double_t *param)
--- Fill vector param[4] with the bounding cylinder parameters. The order
 is the following : Rmin, Rmax, Phi1, Phi2
void InspectShape()
 print shape parameters
TBuffer3D * MakeBuffer3D()
 Creates a TBuffer3D describing *this* shape.
 Coordinates are in local reference frame.
void SetSegsAndPols(TBuffer3D &buff)
 Fill TBuffer3D structure for segments and polygons.
Double_t Safety(Double_t *point, Bool_t in)
 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 SafetyS(Double_t *point, Bool_t in, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Int_t skipz)
 Static method to compute the closest distance from given point to this shape.
void SavePrimitive(ostream &out, Option_t * /*option*/ /*= ""*/)
 Save a primitive as a C++ statement(s) on output stream "out".
void SetTubsDimensions(Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
 Set dimensions of the tube segment.
void SetDimensions(Double_t *param)
 Set dimensions of the tube segment starting from a list.
void SetPoints(Double_t *points)
 Create tube segment mesh points.
void SetPoints(Float_t *points)
 Create tube segment mesh points.
Int_t GetNmeshVertices()
 Return number of vertices of the mesh representation
void Sizeof3D()
/// fill size of this 3-D object
/    TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
/    if (!painter) return;

/    Int_t n = gGeoManager->GetNsegments()+1;
/    Int_t numPoints = n*4;
/    Int_t numSegs   = n*8;
/    Int_t numPolys  = n*4-2;

/    painter->AddSize3D(numPoints, numSegs, numPolys);
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame)
 Fills a static 3D buffer and returns a reference.
Bool_t Contains(Double_t *point)
Double_t DistFromInsideS(Double_t *point, Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Double_t DistFromInside(Double_t *point, Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0)
Double_t DistFromOutsideS(Double_t *point, Double_t *dir, Double_t rmin, Double_t rmax, Double_t dz)
Double_t DistFromOutside(Double_t *point, Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=0)
Int_t DistancetoPrimitive(Int_t px, Int_t py)
Int_t GetByteCount()
TGeoTubeSeg()
 constructors
Double_t GetPhi1()
Double_t GetPhi2()

Author: Andrei Gheata 24/10/01
Last update: root/geom:$Name: $:$Id: TGeoTube.cxx,v 1.68 2006/07/03 16:10:44 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.