TGeoShape.cxx

Go to the documentation of this file.
00001 // @(#)root/geom:$Id: TGeoShape.cxx 36535 2010-11-08 14:41:54Z agheata $
00002 // Author: Andrei Gheata   31/01/02
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //____________________________________________________________________________
00013 // TGeoShape - Base abstract class for all shapes. 
00014 //____________________________________________________________________________
00015 //
00016 //
00017 //   Shapes are geometrical objects that provide the basic modelling 
00018 // functionality. They provide the definition of the LOCAL frame of coordinates,
00019 // with respect to which they are defined. Any implementation of a shape deriving
00020 // from the base TGeoShape class has to provide methods for :
00021 //  - finding out if a point defined in their local frame is or not contained 
00022 // inside;
00023 //  - computing the distance from a local point to getting outside/entering the
00024 // shape, given a known direction;
00025 //  - computing the maximum distance in any direction from a local point that
00026 // does NOT result in a boundary crossing of the shape (safe distance); 
00027 //  - computing the cosines of the normal vector to the crossed shape surface,
00028 // given a starting local point and an ongoing direction.
00029 //   All the features above are globally managed by the modeller in order to
00030 // provide navigation functionality. In addition to those, shapes have also to
00031 // implement additional specific abstract methods :
00032 //  - computation of the minimal box bounding the shape, given that this box have
00033 // to be aligned with the local coordinates;
00034 //  - algorithms for dividing the shape along a given axis and producing resulting
00035 // divisions volumes.
00036 //
00037 //   The modeler currently provides a set of 16 basic shapes, which we will call
00038 // primitives. It also provides a special class allowing the creation of shapes
00039 // made as a result of boolean operations between primitives. These are called
00040 // composite shapes and the composition operation can be recursive (composition
00041 // of composites). This allows the creation of a quite large number of different
00042 // shape topologies and combinations.
00043 //
00044 //   Shapes are named objects and register themselves to the manager class at 
00045 // creation time. This is responsible for their final deletion. Shapes 
00046 // can be created without name if their retreival by name is no needed. Generally
00047 // shapes are objects that are usefull only at geometry creation stage. The pointer
00048 // to a shape is in fact needed only when referring to a given volume and it is 
00049 // always accessible at that level. A shape may be referenced by several volumes,
00050 // therefore its deletion is not possible once volumes were defined based on it.
00051 // 
00052 // 
00053 // 
00054 // Creating shapes
00055 //================
00056 //   Shape objects embeed only the minimum set of parameters that are fully
00057 // describing a valid physical shape. For instance, a tube is represented by
00058 // its half length, the minimum radius and the maximum radius. Shapes are used
00059 // togeather with media in order to create volumes, which in their turn 
00060 // are the main components of the geometrical tree. A specific shape can be created
00061 // stand-alone : 
00062 //
00063 //   TGeoBBox *box = new TGeoBBox("s_box", halfX, halfY, halfZ); // named
00064 //   TGeoTube *tub = new TGeoTube(rmin, rmax, halfZ);            // no name
00065 //   ...  (see each specific shape constructors)
00066 //
00067 //   Sometimes it is much easier to create a volume having a given shape in one 
00068 // step, since shapes are not direcly linked in the geometrical tree but volumes 
00069 // are :
00070 //
00071 //   TGeoVolume *vol_box = gGeoManager->MakeBox("BOX_VOL", "mat1", halfX, halfY, halfZ);
00072 //   TGeoVolume *vol_tub = gGeoManager->MakeTube("TUB_VOL", "mat2", rmin, rmax, halfZ);
00073 //   ...  (see MakeXXX() utilities in TGeoManager class)
00074 //
00075 //
00076 // Shape queries
00077 //===============
00078 // Note that global queries related to a geometry are handled by the manager class.
00079 // However, shape-related queries might be sometimes usefull.
00080 //
00081 // A) Bool_t TGeoShape::Contains(Double_t *point[3])
00082 //   - this method returns true if POINT is actually inside the shape. The point
00083 // has to be defined in the local shape reference. For instance, for a box having
00084 // DX, DY and DZ half-lengths a point will be considered inside if :
00085 //   | -DX <= point[0] <= DX
00086 //   | -DY <= point[1] <= DY
00087 //   | -DZ <= point[2] <= DZ
00088 //
00089 // B) Double_t TGeoShape::DistFromInside(Double_t *point[3], Double_t *dir[3],
00090 //                                  Int_t iact, Double_t step, Double_t *safe)
00091 //   - computes the distance to exiting a shape from a given point INSIDE, along
00092 // a given direction. The direction is given by its director cosines with respect
00093 // to the local shape coordinate system. This method provides additional
00094 // information according the value of IACT input parameter :
00095 //   IACT = 0     => compute only safe distance and fill it at the location 
00096 //                   given by SAFE
00097 //   IACT = 1     => a proposed STEP is supplied. The safe distance is computed
00098 //                   first. If this is bigger than STEP than the proposed step
00099 //                   is approved and returned by the method since it does not
00100 //                   cross the shape boundaries. Otherwise, the distance to
00101 //                   exiting the shape is computed and returned.
00102 //   IACT = 2     => compute both safe distance and distance to exiting, ignoring
00103 //                   the proposed step.
00104 //   IACT > 2     => compute only the distance to exiting, ignoring anything else.
00105 //
00106 // C) Double_t TGeoShape::DistFromOutside(Double_t *point[3], Double_t *dir[3],
00107 //                                  Int_t iact, Double_t step, Double_t *safe)
00108 //   - computes the distance to entering a shape from a given point OUTSIDE. Acts
00109 // in the same way as B).
00110 //
00111 // D) Double_t Safety(Double_t *point[3], Bool_t inside)
00112 //
00113 //   - compute maximum shift of a point in any direction that does not change its
00114 // INSIDE/OUTSIDE state (does not cross shape boundaries). The state of the point
00115 // have to be properly supplied.
00116 //
00117 // E) Double_t *Normal(Double_t *point[3], Double_t *dir[3], Bool_t inside)
00118 //
00119 //   - returns director cosines of normal to the crossed shape surface from a
00120 // given point towards a direction. One has to specify if the point is inside 
00121 // or outside shape. According to this, the normal will be outwards or inwards
00122 // shape respectively. Normal components are statically stored by shape class,
00123 // so it has to be copied after retreival in a different array. 
00124 //
00125 // Dividing shapes
00126 //=================
00127 //   Shapes can generally be divided along a given axis. Supported axis are
00128 // X, Y, Z, Rxy, Phi, Rxyz. A given shape cannot be divided however on any axis.
00129 // The general rule is that that divisions are possible on whatever axis that
00130 // produces still known shapes as slices. The division of shapes should not be
00131 // performed by TGeoShape::Divide() calls, but rather by TGeoVolume::Divide().
00132 // The algorithm for dividing a specific shape is known by the shape object, but
00133 // is always invoked in a generic way from the volume level. Details on how to
00134 // do that can be found in TGeoVolume class. One can see how all division options
00135 // are interpreted and which is their result inside specific shape classes.
00136 //_____________________________________________________________________________
00137 //
00138 //Begin_Html
00139 /*
00140 <img src="gif/t_shape.jpg">
00141 */
00142 //End_Html
00143 
00144 #include "TObjArray.h"
00145 #include "TEnv.h"
00146 #include "TError.h"
00147 
00148 #include "TGeoMatrix.h"
00149 #include "TGeoManager.h"
00150 #include "TGeoVolume.h"
00151 #include "TGeoShape.h"
00152 #include "TVirtualGeoPainter.h"
00153 #include "TBuffer3D.h"
00154 #include "TBuffer3DTypes.h"
00155 #include "TMath.h"
00156 
00157 ClassImp(TGeoShape)
00158 
00159 TGeoMatrix *TGeoShape::fgTransform = NULL;
00160 Double_t    TGeoShape::fgEpsMch = 2.220446049250313e-16;
00161 //_____________________________________________________________________________
00162 TGeoShape::TGeoShape()
00163 {
00164 // Default constructor
00165    fShapeBits = 0;
00166    fShapeId   = 0;
00167    if (!gGeoManager) {
00168       gGeoManager = new TGeoManager("Geometry", "default geometry");
00169       // gROOT->AddGeoManager(gGeoManager);
00170    }
00171 //   fShapeId = gGeoManager->GetListOfShapes()->GetSize();
00172 //   gGeoManager->AddShape(this);
00173 }
00174 
00175 //_____________________________________________________________________________
00176 TGeoShape::TGeoShape(const char *name)
00177           :TNamed(name, "")
00178 {
00179 // Default constructor
00180    fShapeBits = 0;
00181    fShapeId   = 0;
00182    if (!gGeoManager) {
00183       gGeoManager = new TGeoManager("Geometry", "default geometry");
00184       // gROOT->AddGeoManager(gGeoManager);
00185    }
00186    fShapeId = gGeoManager->GetListOfShapes()->GetSize();
00187    gGeoManager->AddShape(this);
00188 }
00189 
00190 //_____________________________________________________________________________
00191 TGeoShape::~TGeoShape()
00192 {
00193 // Destructor
00194    if (gGeoManager) gGeoManager->GetListOfShapes()->Remove(this);
00195 }
00196 
00197 //_____________________________________________________________________________
00198 Double_t TGeoShape::ComputeEpsMch()
00199 {
00200 // Compute machine round-off double precision error as the smallest number that
00201 // if added to 1.0 is different than 1.0.
00202    Double_t temp1 = 1.0; 
00203    Double_t temp2 = 1.0 + temp1;
00204    Double_t mchEps;
00205    while (temp2>1.0) {
00206       mchEps = temp1;
00207       temp1 /= 2;
00208       temp2 = 1.0 + temp1;
00209    }
00210    fgEpsMch = mchEps;
00211    return fgEpsMch;
00212 }   
00213    
00214 //_____________________________________________________________________________
00215 Double_t TGeoShape::EpsMch()
00216 {
00217    //static function returning the machine round-off error
00218    
00219    return fgEpsMch;
00220 }
00221    
00222 //_____________________________________________________________________________
00223 const char *TGeoShape::GetName() const
00224 {
00225 // Get the shape name.
00226    if (!strlen(fName)) {
00227       return ((TObject *)this)->ClassName();
00228    }
00229    return TNamed::GetName();
00230 }
00231 
00232 //_____________________________________________________________________________
00233 Int_t TGeoShape::ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
00234 {
00235 // Returns distance to shape primitive mesh.
00236    TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
00237    if (!painter) return 9999;
00238    return painter->ShapeDistancetoPrimitive(this, numpoints, px, py);
00239 }
00240 
00241 //_____________________________________________________________________________
00242 Bool_t TGeoShape::IsCloseToPhi(Double_t epsil, Double_t *point, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
00243 {
00244 // True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,s2
00245    Double_t saf1 = TGeoShape::Big();
00246    Double_t saf2 = TGeoShape::Big();
00247    if (point[0]*c1+point[1]*s1 >= 0) saf1 = TMath::Abs(-point[0]*s1 + point[1]*c1);
00248    if (point[0]*c2+point[1]*s2 >= 0) saf2 = TMath::Abs(point[0]*s2 - point[1]*c2);
00249    Double_t saf = TMath::Min(saf1,saf2);
00250    if (saf<epsil) return kTRUE;
00251    return kFALSE;
00252 }   
00253 
00254 //_____________________________________________________________________________
00255 Bool_t TGeoShape::IsInPhiRange(Double_t *point, Double_t phi1, Double_t phi2)
00256 {
00257 // Static method to check if a point is in the phi range (phi1, phi2) [degrees]
00258    Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
00259    while (phi<phi1) phi+=360.;
00260    Double_t ddp = phi-phi1;
00261    if (ddp>phi2-phi1) return kFALSE;
00262    return kTRUE;
00263 }   
00264 
00265 //_____________________________________________________________________________  
00266 Bool_t TGeoShape::IsCrossingSemiplane(Double_t *point, Double_t *dir, Double_t cphi, Double_t sphi, Double_t &snext, Double_t &rxy)
00267 {
00268 // Compute distance from POINT to semiplane defined by PHI angle along DIR. Computes
00269 // also radius at crossing point. This might be negative in case the crossing is
00270 // on the other side of the semiplane.
00271    snext = rxy = TGeoShape::Big();
00272    Double_t nx = -sphi;
00273    Double_t ny = cphi;
00274    Double_t rxy0 = point[0]*cphi+point[1]*sphi;
00275    Double_t rdotn = point[0]*nx + point[1]*ny;
00276    if (TMath::Abs(rdotn)<TGeoShape::Tolerance()) {
00277       snext = 0.0;
00278       rxy = rxy0;
00279       return kTRUE;
00280    }
00281    if (rdotn<0) {
00282       rdotn = -rdotn;
00283    } else {
00284       nx = -nx;
00285       ny = -ny;
00286    }
00287    Double_t ddotn = dir[0]*nx + dir[1]*ny;
00288    if (ddotn<=0) return kFALSE;
00289    snext = rdotn/ddotn;
00290    rxy = rxy0+snext*(dir[0]*cphi+dir[1]*sphi);
00291    if (rxy<0) return kFALSE;
00292    return kTRUE;
00293 }
00294 
00295 //_____________________________________________________________________________  
00296 Bool_t TGeoShape::IsSameWithinTolerance(Double_t a, Double_t b)
00297 {
00298 // Check if two numbers differ with less than a tolerance.
00299    if (TMath::Abs(a-b)<1.E-10) return kTRUE;
00300    return kFALSE;
00301 }   
00302 
00303 //_____________________________________________________________________________  
00304 Bool_t TGeoShape::IsSegCrossing(Double_t x1, Double_t y1, Double_t x2, Double_t y2,Double_t x3, Double_t y3,Double_t x4, Double_t y4)
00305 {
00306 // Check if segments (A,B) and (C,D) are crossing,
00307 // where: A(x1,y1), B(x2,y2), C(x3,y3), D(x4,y4)
00308    Double_t eps = TGeoShape::Tolerance();
00309    Bool_t stand1 = kFALSE;
00310    Double_t dx1 = x2-x1;
00311    Bool_t stand2 = kFALSE;
00312    Double_t dx2 = x4-x3;
00313    Double_t xm = 0.;
00314    Double_t ym = 0.;
00315    Double_t a1 = 0.;
00316    Double_t b1 = 0.;
00317    Double_t a2 = 0.;
00318    Double_t b2 = 0.;
00319    if (TMath::Abs(dx1) < eps) stand1 = kTRUE;
00320    if (TMath::Abs(dx2) < eps) stand2 = kTRUE;
00321    if (!stand1) {
00322       a1 = (x2*y1-x1*y2)/dx1;
00323       b1 = (y2-y1)/dx1;
00324    }   
00325    if (!stand2) {
00326       a2 = (x4*y3-x3*y4)/dx2;
00327       b2 = (y4-y3)/dx2;
00328    }   
00329    if (stand1 && stand2) {
00330       // Segments parallel and vertical
00331       if (TMath::Abs(x1-x3)<eps) {
00332          // Check if segments are overlapping
00333          if ((y3-y1)*(y3-y2)<-eps || (y4-y1)*(y4-y2)<-eps ||
00334              (y1-y3)*(y1-y4)<-eps || (y2-y3)*(y2-y4)<-eps) return kTRUE;
00335          return kFALSE;
00336       }
00337       // Different x values
00338       return kFALSE;
00339    }
00340    
00341    if (stand1) {
00342       // First segment vertical
00343       xm = x1;
00344       ym = a2+b2*xm;
00345    } else {
00346       if (stand2) {
00347          // Second segment vertical
00348          xm = x3;
00349          ym = a1+b1*xm;
00350       } else {
00351          // Normal crossing
00352          if (TMath::Abs(b1-b2)<eps) {
00353             // Parallel segments, are they aligned
00354             if (TMath::Abs(y3-(a1+b1*x3))>eps) return kFALSE;
00355             // Aligned segments, are they overlapping
00356             if ((x3-x1)*(x3-x2)<-eps || (x4-x1)*(x4-x2)<-eps ||
00357                 (x1-x3)*(x1-x4)<-eps || (x2-x3)*(x2-x4)<-eps) return kTRUE;
00358             return kFALSE;
00359          }
00360          xm = (a1-a2)/(b2-b1);
00361          ym = (a1*b2-a2*b1)/(b2-b1);
00362       }
00363    }
00364    // Check if crossing point is both between A,B and C,D
00365    Double_t check = (xm-x1)*(xm-x2)+(ym-y1)*(ym-y2);
00366    if (check > -eps) return kFALSE;
00367    check = (xm-x3)*(xm-x4)+(ym-y3)*(ym-y4);
00368    if (check > -eps) return kFALSE;
00369    return kTRUE;
00370 }         
00371 
00372 //_____________________________________________________________________________
00373 Double_t TGeoShape::DistToPhiMin(Double_t *point, Double_t *dir, Double_t s1, Double_t c1,
00374                                  Double_t s2, Double_t c2, Double_t sm, Double_t cm, Bool_t in)
00375 {
00376 // compute distance from point (inside phi) to both phi planes. Return minimum.
00377    Double_t sfi1=TGeoShape::Big();
00378    Double_t sfi2=TGeoShape::Big();
00379    Double_t s=0;
00380    Double_t un = dir[0]*s1-dir[1]*c1;
00381    if (!in) un=-un;
00382    if (un>0) {
00383       s=-point[0]*s1+point[1]*c1;
00384       if (!in) s=-s;
00385       if (s>=0) {
00386          s /= un;
00387          if (((point[0]+s*dir[0])*sm-(point[1]+s*dir[1])*cm)>=0) sfi1=s;
00388       }
00389    }
00390    un = -dir[0]*s2+dir[1]*c2;
00391    if (!in) un=-un;
00392    if (un>0) {
00393       s=point[0]*s2-point[1]*c2;
00394       if (!in) s=-s;
00395       if (s>=0) {
00396          s /= un;
00397          if ((-(point[0]+s*dir[0])*sm+(point[1]+s*dir[1])*cm)>=0) sfi2=s;
00398       }
00399    }
00400    return TMath::Min(sfi1, sfi2);
00401 }
00402 
00403 //_____________________________________________________________________________
00404 void TGeoShape::NormalPhi(Double_t *point, Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
00405 {
00406 // Static method to compute normal to phi planes.
00407    Double_t saf1 = TGeoShape::Big();
00408    Double_t saf2 = TGeoShape::Big();
00409    if (point[0]*c1+point[1]*s1 >= 0) saf1 = TMath::Abs(-point[0]*s1 + point[1]*c1);
00410    if (point[0]*c2+point[1]*s2 >= 0) saf2 =  TMath::Abs(point[0]*s2 - point[1]*c2);
00411    Double_t c,s;
00412    if (saf1<saf2) {
00413       c=c1;
00414       s=s1;
00415    } else {
00416       c=c2;
00417       s=s2;
00418    }
00419    norm[2] = 0;
00420    norm[0] = -s;
00421    norm[1] = c;
00422    if (dir[0]*norm[0]+dir[1]*norm[1] < 0) { 
00423       norm[0] = s;
00424       norm[1] = -c;
00425    }
00426 }           
00427  
00428 //_____________________________________________________________________________
00429 Double_t TGeoShape::SafetyPhi(Double_t *point, Bool_t in, Double_t phi1, Double_t phi2)
00430 {
00431 // Static method to compute safety w.r.t a phi corner defined by cosines/sines
00432 // of the angles phi1, phi2.
00433    Bool_t inphi = TGeoShape::IsInPhiRange(point, phi1, phi2);
00434    if (inphi && !in) return -TGeoShape::Big(); 
00435    phi1 *= TMath::DegToRad();
00436    phi2 *= TMath::DegToRad();  
00437    Double_t c1 = TMath::Cos(phi1);
00438    Double_t s1 = TMath::Sin(phi1);
00439    Double_t c2 = TMath::Cos(phi2);
00440    Double_t s2 = TMath::Sin(phi2);
00441    Double_t rsq = point[0]*point[0]+point[1]*point[1];
00442    Double_t rproj = point[0]*c1+point[1]*s1;
00443    Double_t safsq = rsq-rproj*rproj;
00444    if (safsq<0) return 0.;
00445    Double_t saf1 = (rproj<0)?TGeoShape::Big():TMath::Sqrt(safsq);
00446    rproj = point[0]*c2+point[1]*s2;
00447    safsq = rsq-rproj*rproj;
00448    if (safsq<0) return 0.;   
00449    Double_t saf2 = (rproj<0)?TGeoShape::Big():TMath::Sqrt(safsq);
00450    Double_t safe = TMath::Min(saf1, saf2); // >0
00451    if (safe>1E10) {
00452       if (in) return TGeoShape::Big();
00453       return -TGeoShape::Big();
00454    }
00455    return safe;   
00456 }        
00457 
00458 //_____________________________________________________________________________
00459 void TGeoShape::SetShapeBit(UInt_t f, Bool_t set)
00460 {
00461 // Equivalent of TObject::SetBit.
00462    if (set) {
00463       SetShapeBit(f);
00464    } else {
00465       ResetShapeBit(f);
00466    }
00467 }
00468 
00469 //_____________________________________________________________________________
00470 TGeoMatrix *TGeoShape::GetTransform()
00471 {
00472 // Returns current transformation matrix that applies to shape.
00473    return fgTransform;
00474 }   
00475 
00476 //_____________________________________________________________________________
00477 void TGeoShape::SetTransform(TGeoMatrix *matrix)
00478 {
00479 // Set current transformation matrix that applies to shape.
00480    fgTransform = matrix;
00481 }   
00482 
00483 //_____________________________________________________________________________
00484 void TGeoShape::TransformPoints(Double_t *points, UInt_t NbPnts) const
00485 {
00486    // Tranform a set of points (LocalToMaster)
00487    UInt_t i,j;
00488    Double_t dlocal[3];
00489    Double_t dmaster[3];
00490    if (fgTransform) {
00491       for (j = 0; j < NbPnts; j++) {
00492          i = 3*j;
00493          fgTransform->LocalToMaster(&points[i], dmaster);
00494          points[i]   = dmaster[0];
00495          points[i+1] = dmaster[1];
00496          points[i+2] = dmaster[2];
00497       }
00498       return;   
00499    }   
00500    if (!gGeoManager) return;
00501    Bool_t bomb = (gGeoManager->GetBombMode()==0)?kFALSE:kTRUE;
00502 
00503    for (j = 0; j < NbPnts; j++) {
00504       i = 3*j;
00505       dlocal[0] = points[3*j];
00506       dlocal[1] = points[3*j+1];
00507       dlocal[2] = points[3*j+2];
00508       if (gGeoManager->IsMatrixTransform()) {
00509          TGeoHMatrix *glmat = gGeoManager->GetGLMatrix();
00510          if (bomb) glmat->LocalToMasterBomb(&points[i], dmaster);
00511          else      glmat->LocalToMaster(&points[i], dmaster);
00512       } else {
00513          if (bomb) gGeoManager->LocalToMasterBomb(&points[i], dmaster);
00514          else      gGeoManager->LocalToMaster(&points[i],dmaster);
00515       }
00516       points[i]   = dmaster[0];
00517       points[i+1] = dmaster[1];
00518       points[i+2] = dmaster[2];
00519    }
00520 }
00521 
00522 //_____________________________________________________________________________
00523 void TGeoShape::FillBuffer3D(TBuffer3D & buffer, Int_t reqSections, Bool_t localFrame) const
00524 {
00525    // Fill the supplied buffer, with sections in desired frame
00526    // See TBuffer3D.h for explanation of sections, frame etc.
00527   
00528    // Catch this common potential error here
00529    // We have to set kRawSize (unless already done) to allocate buffer space 
00530    // before kRaw can be filled
00531    if (reqSections & TBuffer3D::kRaw) {
00532       if (!(reqSections & TBuffer3D::kRawSizes) && !buffer.SectionsValid(TBuffer3D::kRawSizes)) {
00533          R__ASSERT(kFALSE);
00534       }
00535    }
00536 
00537    if (reqSections & TBuffer3D::kCore) {
00538       // If writing core section all others will be invalid
00539       buffer.ClearSectionsValid();
00540  
00541       // Check/grab some objects we need
00542       if (!gGeoManager) { 
00543          R__ASSERT(kFALSE); 
00544          return; 
00545       }
00546       const TGeoVolume * paintVolume = gGeoManager->GetPaintVolume();
00547       if (!paintVolume) paintVolume = gGeoManager->GetTopVolume();
00548       if (!paintVolume) { 
00549          buffer.fID = const_cast<TGeoShape *>(this);
00550          buffer.fColor = 0;
00551          buffer.fTransparency = 0;
00552 //         R__ASSERT(kFALSE); 
00553 //         return; 
00554       } else {
00555          buffer.fID = const_cast<TGeoVolume *>(paintVolume);
00556          buffer.fColor = paintVolume->GetLineColor();
00557 
00558          buffer.fTransparency = paintVolume->GetTransparency();
00559          Double_t visdensity = gGeoManager->GetVisDensity();
00560          if (visdensity>0 && paintVolume->GetMedium()) {
00561             if (paintVolume->GetMaterial()->GetDensity() < visdensity) {
00562                buffer.fTransparency = 90;
00563             }
00564          }   
00565       }
00566 
00567       buffer.fLocalFrame = localFrame;
00568       Bool_t r1,r2=kFALSE;
00569       r1 = gGeoManager->IsMatrixReflection();
00570       if (paintVolume && paintVolume->GetShape()) {
00571          if (paintVolume->GetShape()->IsReflected()) {
00572          // Temporary trick to deal with reflected shapes.
00573          // Still lighting gets wrong...
00574             if (buffer.Type() < TBuffer3DTypes::kTube) r2 = kTRUE;
00575          }
00576       }      
00577       buffer.fReflection = ((r1&(!r2))|(r2&!(r1)));
00578       
00579       // Set up local -> master translation matrix
00580       if (localFrame) {
00581          TGeoMatrix * localMasterMat = 0;
00582          if (TGeoShape::GetTransform()) {
00583             localMasterMat = TGeoShape::GetTransform();
00584          } else {   
00585             localMasterMat = gGeoManager->GetCurrentMatrix();
00586 
00587          // For overlap drawing the correct matrix needs to obtained in 
00588          // from GetGLMatrix() - this should not be applied in the case
00589          // of composite shapes
00590             if (gGeoManager->IsMatrixTransform() && !IsComposite()) {
00591                localMasterMat = gGeoManager->GetGLMatrix();
00592             }
00593          }
00594          if (!localMasterMat) { 
00595             R__ASSERT(kFALSE); 
00596             return; 
00597          }
00598          localMasterMat->GetHomogenousMatrix(buffer.fLocalMaster);
00599       } else {
00600          buffer.SetLocalMasterIdentity();
00601       }
00602 
00603       buffer.SetSectionsValid(TBuffer3D::kCore);
00604    }
00605 }
00606 
00607 //_____________________________________________________________________________
00608 Int_t TGeoShape::GetBasicColor() const
00609 {
00610 // Get the basic color (0-7).
00611    Int_t basicColor = 0; // TODO: Check on sensible fallback
00612    if (gGeoManager) {
00613       const TGeoVolume * volume = gGeoManager->GetPaintVolume();
00614       if (volume) {
00615             basicColor = ((volume->GetLineColor() %8) -1) * 4;
00616             if (basicColor < 0) basicColor = 0;
00617       }
00618    }
00619    return basicColor;
00620 }
00621 
00622 //_____________________________________________________________________________
00623 const TBuffer3D &TGeoShape::GetBuffer3D(Int_t /*reqSections*/, Bool_t /*localFrame*/) const
00624 {
00625    // Stub implementation to avoid forcing implementation at this stage
00626    static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
00627    Warning("GetBuffer3D", "this must be implemented for shapes in a TGeoPainter hierarchy. This will be come a pure virtual fn eventually.");
00628    return buffer;
00629 }
00630 
00631 //_____________________________________________________________________________
00632 const char *TGeoShape::GetPointerName() const
00633 {
00634 // Provide a pointer name containing uid.
00635    static TString name;
00636    Int_t uid = GetUniqueID();
00637    if (uid) name = TString::Format("p%s_%d", GetName(),uid);
00638    else     name = TString::Format("p%s", GetName());
00639    return name.Data();
00640 }
00641 
00642 //_____________________________________________________________________________
00643 void TGeoShape::ExecuteEvent(Int_t event, Int_t px, Int_t py)
00644 {
00645 // Execute mouse actions on this shape.
00646    if (!gGeoManager) return;
00647    TVirtualGeoPainter *painter = gGeoManager->GetPainter();
00648    painter->ExecuteShapeEvent(this, event, px, py);
00649 }
00650 
00651 //_____________________________________________________________________________
00652 void TGeoShape::Draw(Option_t *option)
00653 {
00654 // Draw this shape.
00655    TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
00656    if (option && strlen(option) > 0) {
00657       painter->DrawShape(this, option); 
00658    } else {
00659       painter->DrawShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
00660    }  
00661 }
00662 
00663 //_____________________________________________________________________________
00664 void TGeoShape::Paint(Option_t *option)
00665 {
00666 // Paint this shape.
00667    TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
00668    if (option && strlen(option) > 0) {
00669       painter->PaintShape(this, option); 
00670    } else {
00671       painter->PaintShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
00672    }  
00673 }

Generated on Tue Jul 5 14:12:57 2011 for ROOT_528-00b_version by  doxygen 1.5.1