TGeoVolume.cxx

Go to the documentation of this file.
00001 // @(#)root/geom:$Id: TGeoVolume.cxx 36632 2010-11-12 15:57:07Z agheata $
00002 // Author: Andrei Gheata   30/05/02
00003 // Divide(), CheckOverlaps() implemented by Mihaela Gheata
00004 
00005 /*************************************************************************
00006  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00007  * All rights reserved.                                                  *
00008  *                                                                       *
00009  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00010  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00011  *************************************************************************/
00012 
00013 //Begin_Html
00014 /*
00015 <img src="gif/t_volume.jpg">
00016 */
00017 //End_Html
00018 
00019 ////////////////////////////////////////////////////////////////////////////////
00020 //   TGeoVolume - the base class representing solids. 
00021 //
00022 //   Volumes are the basic objects used in building the geometrical hierarchy.
00023 // They represent unpositioned objects but store all information about the
00024 // placement of the other volumes they may contain. Therefore a volume can
00025 // be replicated several times in the geometry. In order to create a volume, one
00026 // has to put togeather a shape and a medium which are already defined. Volumes
00027 // have to be named by users at creation time. Every different name may represent a 
00028 // an unique volume object, but may also represent more general a family (class)
00029 // of volume objects having the same shape type and medium, but possibly
00030 // different shape parameters. It is the user's task to provide different names
00031 // for different volume families in order to avoid ambiguities at tracking time.
00032 // A generic family rather than a single volume is created only in two cases : 
00033 // when a generic shape is provided to the volume constructor or when a division
00034 // operation is applied. Each volume in the geometry stores an unique
00035 // ID corresponding to its family. In order to ease-up their creation, the manager
00036 // class is providing an API that allows making a shape and a volume in a single step.
00037 //
00038 //   Volumes are objects that can be visualized, therefore having visibility,
00039 // colour, line and fill attributes that can be defined or modified any time after
00040 // the volume creation. It is advisable however to define these properties just
00041 // after the first creation of a volume namespace, since in case of volume families
00042 // any new member created by the modeler inherits these properties. 
00043 //
00044 //    In order to provide navigation features, volumes have to be able to find
00045 // the proper container of any point defined in the local reference frame. This
00046 // can be the volume itself, one of its positioned daughter volumes or none if 
00047 // the point is actually outside. On the other hand, volumes have to provide also
00048 // other navigation methods such as finding the distances to its shape boundaries
00049 // or which daughter will be crossed first. The implementation of these features
00050 // is done at shape level, but the local mother-daughters management is handled
00051 // by volumes that builds additional optimisation structures upon geometry closure.
00052 // In order to have navigation features properly working one has to follow the
00053 // general rules for building a valid geometry (see TGeoManager class).
00054 //
00055 //   Now let's make a simple volume representing a copper wire. We suppose that
00056 // a medium is already created (see TGeoMedium class on how to create media). 
00057 // We will create a TUBE shape for our wire, having Rmin=0cm, Rmax=0.01cm
00058 // and a half-length dZ=1cm :
00059 //
00060 //   TGeoTube *tube = new TGeoTube("wire_tube", 0, 0.01, 1);
00061 //
00062 // One may ommit the name for the shape if no retreiving by name is further needed
00063 // during geometry building. The same shape can be shared by different volumes 
00064 // having different names and materials. Now let's make the volume for our wire.
00065 // The prototype for volumes constructor looks like :
00066 //
00067 //   TGeoVolume::TGeoVolume(const char *name, TGeoShape *shape, TGeoMedium *med)
00068 //
00069 // Since TGeoTube derives brom the base shape class, we can provide it to the volume
00070 // constructor :
00071 //
00072 //   TGeoVolume *wire_co = new TGeoVolume("WIRE_CO", tube, ptrCOPPER);
00073 //
00074 // Do not bother to delete neither the media, shapes or volumes that you have
00075 // created since all will be automatically cleaned on exit by the manager class.
00076 // If we would have taken a look inside TGeoManager::MakeTube() method, we would
00077 // have been able to create our wire with a single line :
00078 //
00079 //   TGeoVolume *wire_co = gGeoManager->MakeTube("WIRE_CO", ptrCOPPER, 0, 0.01, 1);
00080 //
00081 // The same applies for all primitive shapes, for which there can be found
00082 // corresponding MakeSHAPE() methods. Their usage is much more convenient unless 
00083 // a shape has to be shared between more volumes. Let's make now an aluminium wire 
00084 // having the same shape, supposing that we have created the copper wire with the 
00085 // line above :
00086 //
00087 //   TGeoVolume *wire_al = new TGeoVolume("WIRE_AL", wire_co->GetShape(), ptrAL);
00088 //
00089 // Now that we have learned how to create elementary volumes, let's see how we
00090 // can create a geometrical hierarchy.
00091 //
00092 //
00093 //   Positioning volumes
00094 // -----------------------
00095 //
00096 //   When creating a volume one does not specify if this will contain or not other
00097 // volumes. Adding daughters to a volume implies creating those and adding them 
00098 // one by one to the list of daughters. Since the volume has to know the position 
00099 // of all its daughters, we will have to supply at the same time a geometrical 
00100 // transformation with respect to its local reference frame for each of them.
00101 // The objects referencing a volume and a transformation are called NODES and
00102 // their creation is fully handled by the modeler. They represent the link 
00103 // elements in the hierarchy of volumes. Nodes are unique and distinct geometrical
00104 // objects ONLY from their container point of view. Since volumes can be replicated
00105 // in the geometry, the same node may be found on different branches.
00106 //
00107 //Begin_Html
00108 /*
00109 <img src="gif/t_example.jpg">
00110 */
00111 //End_Html
00112 //
00113 //   An important observation is that volume objects are owned by the TGeoManager
00114 // class. This stores a list of all volumes in the geometry, that is cleaned
00115 // upon destruction.
00116 //
00117 //   Let's consider positioning now our wire in the middle of a gas chamber. We 
00118 // need first to define the gas chamber :
00119 //
00120 //   TGeoVolume *chamber = gGeoManager->MakeTube("CHAMBER", ptrGAS, 0, 1, 1);
00121 // 
00122 // Now we can put the wire inside :
00123 //
00124 //   chamber->AddNode(wire_co, 1);
00125 //
00126 // If we inspect now the chamber volume in a browser, we will notice that it has 
00127 // one daughter. Of course the gas has some container also, but let's keep it like 
00128 // that for the sake of simplicity. The full prototype of AddNode() is :
00129 //
00130 //   TGeoVolume::AddNode(TGeoVolume *daughter, Int_t usernumber, 
00131 //                       TGeoMatrix *matrix=gGeoIdentity)
00132 //
00133 // Since we did not supplied the third argument, the wire will be positioned with
00134 // an identity transformation inside the chamber. One will notice that the inner
00135 // radii of the wire and chamber are both zero - therefore, aren't the two volumes
00136 // overlapping ? The answer is no, the modeler is even relaying on the fact that
00137 // any daughter is fully contained by its mother. On the other hand, neither of
00138 // the nodes positioned inside a volume should overlap with each other. We will
00139 // see that there are allowed some exceptions to those rules.
00140 //
00141 // Overlapping volumes
00142 // --------------------
00143 //
00144 //   Positioning volumes that does not overlap their neighbours nor extrude
00145 // their container is sometimes quite strong contrain. Some parts of the geometry
00146 // might overlap naturally, e.g. two crossing tubes. The modeller supports such
00147 // cases only if the overlapping nodes are declared by the user. In order to do
00148 // that, one should use TGeoVolume::AddNodeOverlap() instead of TGeoVolume::AddNode().
00149 //   When 2 or more positioned volumes are overlapping, not all of them have to
00150 // be declared so, but at least one. A point inside an overlapping region equally
00151 // belongs to all overlapping nodes, but the way these are defined can enforce
00152 // the modeler to give priorities.
00153 //   The general rule is that the deepest node in the hierarchy containing a point
00154 // have the highest priority. For the same geometry level, non-overlapping is
00155 // prioritized over overlapping. In order to illustrate this, we will consider 
00156 // few examples. We will designate non-overlapping nodes as ONLY and the others
00157 // MANY as in GEANT3, where this concept was introduced:
00158 //   1. The part of a MANY node B extruding its container A will never be "seen" 
00159 // during navigation, as if B was in fact the result of the intersection of A and B.
00160 //   2. If we have two nodes A (ONLY) and B (MANY) inside the same container, all
00161 // points in the overlapping region of A and B will be designated as belonging to A.
00162 //   3. If A an B in the above case were both MANY, points in the overlapping 
00163 // part will be designated to the one defined first. Both nodes must have the 
00164 // same medium.
00165 //   4. The silces of a divided MANY will be as well MANY.
00166 //
00167 // One needs to know that navigation inside geometry parts MANY nodes is much 
00168 // slower. Any overlapping part can be defined based on composite shapes - this
00169 // is always recommended. 
00170 
00171 //   Replicating volumes
00172 // -----------------------
00173 //
00174 //   What can we do if our chamber contains two identical wires instead of one ?
00175 // What if then we would need 1000 chambers in our detector ? Should we create
00176 // 2000 wires and 1000 chamber volumes ? No, we will just need to replicate the
00177 // ones that we have already created.
00178 //
00179 //   chamber->AddNode(wire_co, 1, new TGeoTranslation(-0.2,0,0));
00180 //   chamber->AddNode(wire_co, 2, new TGeoTranslation(0.2,0,0));
00181 //
00182 //   The 2 nodes that we have created inside chamber will both point to a wire_co
00183 // object, but will be completely distinct : WIRE_CO_1 and WIRE_CO_2. We will
00184 // want now to place symetrically 1000 chabmers on a pad, following a pattern
00185 // of 20 rows and 50 columns. One way to do this will be to replicate our chamber
00186 // by positioning it 1000 times in different positions of the pad. Unfortunatelly,
00187 // this is far from being the optimal way of doing what we want.
00188 // Imagine that we would like to find out which of the 1000 chambers is containing
00189 // a (x,y,z) point defined in the pad reference. You will never have to do that,
00190 // since the modeller will take care of it for you, but let's guess what it has 
00191 // to do. The most simple algorithm will just loop over all daughters, convert
00192 // the point from mother to local reference and check if the current chamber
00193 // contains the point or not. This might be efficient for pads with few chambers,
00194 // but definitely not for 1000. Fortunately the modeler is smarter than that and 
00195 // create for each volume some optimization structures called voxels (see Voxelization) 
00196 // to minimize the penalty having too many daughters, but if you have 100 pads like 
00197 // this in your geometry you will anyway loose a lot in your tracking performance.
00198 //
00199 //   The way out when volumes can be arranged acording to simple patterns is the
00200 // usage of divisions. We will describe them in detail later on. Let's think now
00201 // at a different situation : instead of 1000 chambers of the same type, we may
00202 // have several types of chambers. Let's say all chambers are cylindrical and have
00203 // a wire inside, but their dimensions are different. However, we would like all
00204 // to be represented by a single volume family, since they have the same properties.
00205 //
00206 //   Volume families
00207 // ------------------
00208 // A volume family is represented by the class TGeoVolumeMulti. It represents
00209 // a class of volumes having the same shape type and each member will be 
00210 // identified by the same name and volume ID. Any operation applied to a 
00211 // TGeoVolume equally affects all volumes in that family. The creation of a 
00212 // family is generally not a user task, but can be forced in particular cases:
00213 //
00214 //      TGeoManager::Volume(const char *vname, const char *shape, Int_t nmed);
00215 //
00216 // where VNAME is the family name, NMED is the medium number and SHAPE is the
00217 // shape type that can be:
00218 //   box    - for TGeoBBox
00219 //   trd1   - for TGeoTrd1
00220 //   trd2   - for TGeoTrd2
00221 //   trap   - for TGeoTrap
00222 //   gtra   - for TGeoGtra
00223 //   para   - for TGeoPara
00224 //   tube, tubs - for TGeoTube, TGeoTubeSeg
00225 //   cone, cons - for TGeoCone, TgeoCons
00226 //   eltu   - for TGeoEltu
00227 //   ctub   - for TGeoCtub
00228 //   pcon   - for TGeoPcon
00229 //   pgon   - for TGeoPgon
00230 //
00231 // Volumes are then added to a given family upon adding the generic name as node
00232 // inside other volume:
00233 //   TGeoVolume *box_family = gGeoManager->Volume("BOXES", "box", nmed);
00234 //   ...
00235 //   gGeoManager->Node("BOXES", Int_t copy_no, "mother_name", 
00236 //                     Double_t x, Double_t y, Double_t z, Int_t rot_index,
00237 //                     Bool_t is_only, Double_t *upar, Int_t npar);
00238 // here:
00239 //   BOXES   - name of the family of boxes
00240 //   copy_no - user node number for the created node
00241 //   mother_name - name of the volume to which we want to add the node
00242 //   x,y,z   - translation components
00243 //   rot_index   - indx of a rotation matrix in the list of matrices
00244 //   upar    - array of actual shape parameters
00245 //   npar    - number of parameters
00246 // The parameters order and number are the same as in the corresponding shape
00247 // constructors.
00248 //
00249 //   An other particular case where volume families are used is when we want
00250 // that a volume positioned inside a container to match one ore more container
00251 // limits. Suppose we want to position the same box inside 2 different volumes
00252 // and we want the Z size to match the one of each container:
00253 //
00254 //   TGeoVolume *container1 = gGeoManager->MakeBox("C1", imed, 10,10,30);
00255 //   TGeoVolume *container2 = gGeoManager->MakeBox("C2", imed, 10,10,20);
00256 //   TGeoVolume *pvol       = gGeoManager->MakeBox("PVOL", jmed, 3,3,-1);
00257 //   container1->AddNode(pvol, 1);
00258 //   container2->AddNode(pvol, 1);
00259 //
00260 //   Note that the third parameter of PVOL is negative, which does not make sense
00261 // as half-length on Z. This is interpreted as: when positioned, create a box
00262 // replacing all invalid parameters with the corresponding dimensions of the
00263 // container. This is also internally handled by the TGeoVolumeMulti class, which
00264 // does not need to be instanciated by users.
00265 //
00266 //   Dividing volumes
00267 // ------------------
00268 //
00269 //   Volumes can be divided according a pattern. The most simple division can
00270 // be done along one axis, that can be: X, Y, Z, Phi, Rxy or Rxyz. Let's take 
00271 // the most simple case: we would like to divide a box in N equal slices along X
00272 // coordinate, representing a new volume family. Supposing we already have created
00273 // the initial box, this can be done like:
00274 //
00275 //      TGeoVolume *slicex = box->Divide("SLICEX", 1, N);
00276 //
00277 // where SLICE is the name of the new family representing all slices and 1 is the
00278 // slicing axis. The meaning of the axis index is the following: for all volumes
00279 // having shapes like box, trd1, trd2, trap, gtra or para - 1,2,3 means X,Y,Z; for
00280 // tube, tubs, cone, cons - 1 means Rxy, 2 means phi and 3 means Z; for pcon and
00281 // pgon - 2 means phi and 3 means Z; for spheres 1 means R and 2 means phi.
00282 //   In fact, the division operation has the same effect as positioning volumes
00283 // in a given order inside the divided container - the advantage being that the 
00284 // navigation in such a structure is much faster. When a volume is divided, a
00285 // volume family corresponding to the slices is created. In case all slices can
00286 // be represented by a single shape, only one volume is added to the family and
00287 // positioned N times inside the divided volume, otherwise, each slice will be 
00288 // represented by a distinct volume in the family.
00289 //   Divisions can be also performed in a given range of one axis. For that, one
00290 // have to specify also the starting coordinate value and the step:
00291 //
00292 //      TGeoVolume *slicex = box->Divide("SLICEX", 1, N, start, step);
00293 //
00294 // A check is always done on the resulting division range : if not fitting into
00295 // the container limits, an error message is posted. If we will browse the divided
00296 // volume we will notice that it will contain N nodes starting with index 1 upto
00297 // N. The first one has the lower X limit at START position, while the last one
00298 // will have the upper X limit at START+N*STEP. The resulting slices cannot
00299 // be positioned inside an other volume (they are by default positioned inside the
00300 // divided one) but can be further divided and may contain other volumes:
00301 //
00302 //      TGeoVolume *slicey = slicex->Divide("SLICEY", 2, N1);
00303 //      slicey->AddNode(other_vol, index, some_matrix);
00304 //
00305 //   When doing that, we have to remember that SLICEY represents a family, therefore
00306 // all members of the family will be divided on Y and the other volume will be 
00307 // added as node inside all.
00308 //   In the example above all the resulting slices had the same shape as the
00309 // divided volume (box). This is not always the case. For instance, dividing a
00310 // volume with TUBE shape on PHI axis will create equal slices having TUBESEG 
00311 // shape. Other divisions can alsoo create slices having shapes with different
00312 // dimensins, e.g. the division of a TRD1 volume on Z. 
00313 //   When positioning volumes inside slices, one can do it using the generic
00314 // volume family (e.g. slicey). This should be done as if the coordinate system
00315 // of the generic slice was the same as the one of the divided volume. The generic
00316 // slice in case of PHI divisioned is centered with respect to X axis. If the
00317 // family contains slices of different sizes, ani volume positioned inside should 
00318 // fit into the smallest one.
00319 //    Examples for specific divisions according to shape types can be found inside
00320 // shape classes.
00321 // 
00322 //        TGeoVolume::Divide(N, Xmin, Xmax, "X"); 
00323 //
00324 //   The GEANT3 option MANY is supported by TGeoVolumeOverlap class. An overlapping
00325 // volume is in fact a virtual container that does not represent a physical object.
00326 // It contains a list of nodes that are not his daughters but that must be checked 
00327 // always before the container itself. This list must be defined by users and it 
00328 // is checked and resolved in a priority order. Note that the feature is non-standard
00329 // to geometrical modelers and it was introduced just to support conversions of 
00330 // GEANT3 geometries, therefore its extensive usage should be avoided.
00331 //
00332 //   The following picture represent how a simple geometry tree is built in
00333 // memory.
00334 
00335 #include "Riostream.h"
00336 #include "TString.h"
00337 #include "TBrowser.h"
00338 #include "TStyle.h"
00339 #include "TH2F.h"
00340 #include "TPad.h"
00341 #include "TROOT.h"
00342 #include "TClass.h"
00343 #include "TEnv.h"
00344 #include "TMap.h"
00345 #include "TFile.h"
00346 #include "TKey.h"
00347 
00348 #include "TGeoManager.h"
00349 #include "TGeoNode.h"
00350 #include "TGeoMatrix.h"
00351 #include "TVirtualGeoPainter.h"
00352 #include "TGeoVolume.h"
00353 #include "TGeoShapeAssembly.h"
00354 #include "TGeoScaledShape.h"
00355 #include "TGeoCompositeShape.h"
00356 #include "TGeoVoxelFinder.h"
00357 
00358 ClassImp(TGeoVolume)
00359 
00360 //_____________________________________________________________________________
00361 TGeoVolume::TGeoVolume()
00362 { 
00363 // dummy constructor
00364    fNodes    = 0;
00365    fShape    = 0;
00366    fFinder   = 0;
00367    fVoxels   = 0;
00368    fField    = 0;
00369    fMedium   = 0;
00370    fNumber   = 0;
00371    fNtotal   = 0;
00372    fOption   = "";
00373    fGeoManager = gGeoManager;
00374    TObject::ResetBit(kVolumeImportNodes);
00375 }
00376 
00377 //_____________________________________________________________________________
00378 TGeoVolume::TGeoVolume(const char *name, const TGeoShape *shape, const TGeoMedium *med)
00379            :TNamed(name, "")
00380 {
00381 // default constructor
00382    fName = fName.Strip();
00383    fNodes    = 0;
00384    fShape    = (TGeoShape*)shape;
00385    if (fShape) {
00386       if (fShape->TestShapeBit(TGeoShape::kGeoBad)) {
00387          Warning("Ctor", "volume %s has invalid shape", name);
00388       }
00389       if (!fShape->IsValid()) {
00390          Fatal("ctor", "Shape of volume %s invalid. Aborting!", fName.Data());
00391       }   
00392    }      
00393    fFinder   = 0;
00394    fVoxels   = 0;
00395    fField    = 0;
00396    fOption   = "";
00397    fMedium   = (TGeoMedium*)med;
00398    if (fMedium) {
00399       if (fMedium->GetMaterial()) fMedium->GetMaterial()->SetUsed();
00400    }   
00401    fNumber   = 0;
00402    fNtotal   = 0;
00403    fGeoManager = gGeoManager;
00404    if (fGeoManager) fNumber = fGeoManager->AddVolume(this);
00405    TObject::ResetBit(kVolumeImportNodes);
00406 }
00407 
00408 //_____________________________________________________________________________
00409 TGeoVolume::TGeoVolume(const TGeoVolume& gv) :
00410   TNamed(gv),
00411   TGeoAtt(gv),
00412   TAttLine(gv),
00413   TAttFill(gv),
00414   TAtt3D(gv),
00415   fNodes(gv.fNodes),
00416   fShape(gv.fShape),
00417   fMedium(gv.fMedium),
00418   fFinder(gv.fFinder),
00419   fVoxels(gv.fVoxels),
00420   fGeoManager(gv.fGeoManager),
00421   fField(gv.fField),
00422   fOption(gv.fOption),
00423   fNumber(gv.fNumber),
00424   fNtotal(gv.fNtotal)
00425 { 
00426    //copy constructor
00427 }
00428 
00429 //_____________________________________________________________________________
00430 TGeoVolume& TGeoVolume::operator=(const TGeoVolume& gv) 
00431 {
00432    //assignment operator
00433    if(this!=&gv) {
00434       TNamed::operator=(gv);
00435       TGeoAtt::operator=(gv);
00436       TAttLine::operator=(gv);
00437       TAttFill::operator=(gv);
00438       TAtt3D::operator=(gv);
00439       fNodes=gv.fNodes;
00440       fShape=gv.fShape;
00441       fMedium=gv.fMedium;
00442       fFinder=gv.fFinder;
00443       fVoxels=gv.fVoxels;
00444       fGeoManager=gv.fGeoManager;
00445       fField=gv.fField;
00446       fOption=gv.fOption;
00447       fNumber=gv.fNumber;
00448       fNtotal=gv.fNtotal;
00449    } 
00450    return *this;
00451 }
00452 
00453 //_____________________________________________________________________________
00454 TGeoVolume::~TGeoVolume()
00455 {
00456 // Destructor
00457    
00458    if (fNodes) { 
00459       if (!TObject::TestBit(kVolumeImportNodes)) {
00460          fNodes->Delete();
00461       }   
00462       delete fNodes;
00463    }
00464    if (fFinder && !TObject::TestBit(kVolumeImportNodes | kVolumeClone) ) delete fFinder;
00465    if (fVoxels) delete fVoxels;
00466 }
00467 
00468 //_____________________________________________________________________________
00469 void TGeoVolume::Browse(TBrowser *b)
00470 {
00471 // How to browse a volume
00472    if (!b) return;
00473 
00474 //   if (!GetNdaughters()) b->Add(this, GetName(), IsVisible());
00475    TGeoVolume *daughter;
00476    TString title;
00477    for (Int_t i=0; i<GetNdaughters(); i++) { 
00478       daughter = GetNode(i)->GetVolume();
00479       if(!strlen(daughter->GetTitle())) {
00480          if (daughter->IsAssembly()) title.TString::Format("Assembly with %d daughter(s)", 
00481                                                 daughter->GetNdaughters());
00482          else if (daughter->GetFinder()) {
00483             TString s1 = daughter->GetFinder()->ClassName();
00484             s1.ReplaceAll("TGeoPattern","");
00485             title.TString::Format("Volume having %s shape divided in %d %s slices",
00486                        daughter->GetShape()->ClassName(),daughter->GetNdaughters(), s1.Data()); 
00487                        
00488          } else title.TString::Format("Volume with %s shape having %d daughter(s)", 
00489                          daughter->GetShape()->ClassName(),daughter->GetNdaughters());
00490          daughter->SetTitle(title.Data());
00491       }   
00492       b->Add(daughter, daughter->GetName(), daughter->IsVisible());
00493 //      if (IsVisDaughters())
00494 //      b->AddCheckBox(daughter, daughter->IsVisible());
00495 //      else
00496 //         b->AddCheckBox(daughter, kFALSE);
00497    }
00498 }
00499 
00500 //_____________________________________________________________________________
00501 Double_t TGeoVolume::Capacity() const
00502 {
00503 // Computes the capacity of this [cm^3] as the capacity of its shape.
00504 // In case of assemblies, the capacity is computed as the sum of daughter's capacities.
00505    if (!IsAssembly()) return fShape->Capacity();
00506    Double_t capacity = 0.0;
00507    Int_t nd = GetNdaughters();
00508    Int_t i;
00509    for (i=0; i<nd; i++) capacity += GetNode(i)->GetVolume()->Capacity();
00510    return capacity;
00511 }   
00512 
00513 //_____________________________________________________________________________
00514 void TGeoVolume::CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
00515 {
00516 // Shoot nrays with random directions from starting point (startx, starty, startz)
00517 // in the reference frame of this volume. Track each ray until exiting geometry, then
00518 // shoot backwards from exiting point and compare boundary crossing points.
00519    TGeoVolume *old_vol = fGeoManager->GetTopVolume();
00520    if (old_vol!=this) fGeoManager->SetTopVolume((TGeoVolume*)this);
00521    else old_vol=0;
00522    fGeoManager->GetTopVolume()->Draw();
00523    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
00524    painter->CheckGeometry(nrays, startx, starty, startz);
00525 }         
00526 
00527 //_____________________________________________________________________________
00528 void TGeoVolume::CheckOverlaps(Double_t ovlp, Option_t *option) const
00529 {
00530 // Overlap checking tool. Check for illegal overlaps within a limit OVLP.
00531 // Use option="s[number]" to force overlap checking by sampling volume with
00532 // [number] points.
00533 // Ex: myVol->CheckOverlaps(0.01, "s10000000"); // shoot 10000000 points
00534 //     myVol->CheckOverlaps(0.01, "s"); // shoot the default value of 1e6 points
00535 
00536    if (!GetNdaughters() || fFinder) return;
00537    Bool_t sampling = kFALSE;
00538    TString opt(option);
00539    opt.ToLower();
00540    if (opt.Contains("s")) sampling = kTRUE;
00541    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
00542    if (!sampling) fGeoManager->SetNsegments(80);
00543    if (!fGeoManager->IsCheckingOverlaps()) {
00544       fGeoManager->ClearOverlaps();
00545 //      Info("CheckOverlaps", "=== Checking overlaps for volume %s ===\n", GetName());
00546    }   
00547    painter->CheckOverlaps(this, ovlp, option);
00548 //   if (sampling) return;
00549    if (!fGeoManager->IsCheckingOverlaps()) {
00550       fGeoManager->SortOverlaps();
00551       TObjArray *overlaps = fGeoManager->GetListOfOverlaps();
00552       Int_t novlps = overlaps->GetEntriesFast();
00553       TNamed *obj;
00554       TString name;
00555       for (Int_t i=0; i<novlps; i++) {
00556          obj = (TNamed*)overlaps->At(i);
00557          if (novlps<1000) name = TString::Format("ov%03d", i);
00558          else             name = TString::Format("ov%06d", i);
00559          obj->SetName(name);
00560       }   
00561       if (novlps) Info("CheckOverlaps", "Number of illegal overlaps/extrusions for volume %s: %d\n", GetName(), novlps);
00562    }   
00563 }
00564 
00565 //_____________________________________________________________________________
00566 void TGeoVolume::CleanAll()
00567 {
00568 // Clean data of the volume.
00569    ClearNodes();
00570    ClearShape();
00571 }
00572 
00573 //_____________________________________________________________________________
00574 void TGeoVolume::ClearShape()
00575 {
00576 // Clear the shape of this volume from the list held by the current manager.
00577    fGeoManager->ClearShape(fShape);
00578 }   
00579 
00580 //_____________________________________________________________________________
00581 void TGeoVolume::CheckShapes()
00582 {
00583 // check for negative parameters in shapes.
00584 // THIS METHOD LEAVES SOME GARBAGE NODES -> memory leak, to be fixed
00585 //   printf("---Checking daughters of volume %s\n", GetName());
00586    if (fShape->IsRunTimeShape()) {
00587       Error("CheckShapes", "volume %s has run-time shape", GetName());
00588       InspectShape();
00589       return;
00590    }   
00591    if (!fNodes) return;
00592    Int_t nd=fNodes->GetEntriesFast();
00593    TGeoNode *node = 0;
00594    TGeoNode *new_node;
00595    const TGeoShape *shape = 0;
00596    TGeoVolume *old_vol;
00597    for (Int_t i=0; i<nd; i++) {
00598       node=(TGeoNode*)fNodes->At(i);
00599       // check if node has name
00600       if (!strlen(node->GetName())) printf("Daughter %i of volume %s - NO NAME!!!\n",
00601                                            i, GetName());
00602       old_vol = node->GetVolume();
00603       shape = old_vol->GetShape();
00604       if (shape->IsRunTimeShape()) {
00605 //         printf("   Node %s/%s has shape with negative parameters. \n", 
00606 //                 GetName(), node->GetName());
00607 //         old_vol->InspectShape();
00608          // make a copy of the node
00609          new_node = node->MakeCopyNode();
00610          TGeoShape *new_shape = shape->GetMakeRuntimeShape(fShape, node->GetMatrix());
00611          if (!new_shape) {
00612             Error("CheckShapes","cannot resolve runtime shape for volume %s/%s\n",
00613                    GetName(),old_vol->GetName());
00614             continue;
00615          }         
00616          TGeoVolume *new_volume = old_vol->MakeCopyVolume(new_shape);
00617 //         printf(" new volume %s shape params :\n", new_volume->GetName());
00618 //         new_volume->InspectShape();
00619          new_node->SetVolume(new_volume);
00620          // decouple the old node and put the new one instead
00621          fNodes->AddAt(new_node, i);
00622 //         new_volume->CheckShapes();
00623       }
00624    }
00625 }     
00626 
00627 //_____________________________________________________________________________
00628 Int_t TGeoVolume::CountNodes(Int_t nlevels, Int_t option)
00629 {
00630 // Count total number of subnodes starting from this volume, nlevels down
00631 // option = 0 (default) - count only once per volume
00632 // option = 1           - count every time
00633 // option = 2           - count volumes on visible branches
00634 // option = 3           - return maximum level counted already with option = 0
00635    static Int_t maxlevel = 0;
00636    static Int_t nlev = 0;
00637    
00638    if (option<0 || option>3) option = 0;
00639    Int_t visopt = 0;
00640    Int_t nd = GetNdaughters();
00641    Bool_t last = (!nlevels || !nd)?kTRUE:kFALSE;
00642    switch (option) {
00643       case 0:
00644          if (fNtotal) return fNtotal;
00645       case 1:   
00646          fNtotal = 1;
00647          break;
00648       case 2:
00649          visopt = fGeoManager->GetVisOption();
00650          if (!IsVisDaughters()) last = kTRUE;
00651          switch (visopt) {
00652             case TVirtualGeoPainter::kGeoVisDefault:
00653                fNtotal = (IsVisible())?1:0;
00654                break;   
00655             case TVirtualGeoPainter::kGeoVisLeaves:
00656                fNtotal = (IsVisible() && last)?1:0;
00657          }
00658          if (!IsVisibleDaughters()) return fNtotal;
00659          break;
00660       case 3:
00661          return maxlevel;   
00662    }      
00663    if (last) return fNtotal;
00664    if (gGeoManager->GetTopVolume() == this) {
00665       maxlevel=0;
00666       nlev = 0;
00667    }   
00668    if (nlev>maxlevel) maxlevel = nlev;   
00669    TGeoNode *node;
00670    TGeoVolume *vol;
00671    nlev++;
00672    for (Int_t i=0; i<nd; i++) {
00673       node = GetNode(i);
00674       vol = node->GetVolume();
00675       fNtotal += vol->CountNodes(nlevels-1, option);
00676    }
00677    nlev--;
00678    return fNtotal;
00679 }
00680 
00681 //_____________________________________________________________________________
00682 Bool_t TGeoVolume::IsAllInvisible() const
00683 {
00684 // Return TRUE if volume and all daughters are invisible.
00685    if (IsVisible()) return kFALSE;
00686    Int_t nd = GetNdaughters();
00687    for (Int_t i=0; i<nd; i++) if (GetNode(i)->GetVolume()->IsVisible()) return kFALSE;
00688    return kTRUE;
00689 }   
00690 
00691 //_____________________________________________________________________________
00692 void TGeoVolume::InvisibleAll(Bool_t flag)
00693 {
00694 // Make volume and each of it daughters (in)visible.
00695    SetAttVisibility(!flag);
00696    Int_t nd = GetNdaughters();
00697    TObjArray *list = new TObjArray(nd+1);
00698    list->Add(this);
00699    TGeoVolume *vol;
00700    for (Int_t i=0; i<nd; i++) {
00701       vol = GetNode(i)->GetVolume();
00702       vol->SetAttVisibility(!flag);
00703       list->Add(vol);
00704    }
00705    TIter next(gROOT->GetListOfBrowsers());
00706    TBrowser *browser = 0;
00707    while ((browser=(TBrowser*)next())) {
00708       for (Int_t i=0; i<nd+1; i++) {
00709          vol = (TGeoVolume*)list->At(i);
00710          browser->CheckObjectItem(vol, !flag);
00711       }   
00712       browser->Refresh();
00713    }
00714    delete list;
00715    fGeoManager->SetVisOption(4);
00716 }   
00717 
00718 //_____________________________________________________________________________
00719 Bool_t TGeoVolume::IsFolder() const
00720 {
00721 // Return TRUE if volume contains nodes
00722 //   return (GetNdaughters()?kTRUE:kFALSE);
00723    return kTRUE;
00724 }
00725 
00726 //_____________________________________________________________________________
00727 Bool_t TGeoVolume::IsStyleDefault() const
00728 {
00729 // check if the visibility and attributes are the default ones
00730    if (!IsVisible()) return kFALSE;
00731    if (GetLineColor() != gStyle->GetLineColor()) return kFALSE;
00732    if (GetLineStyle() != gStyle->GetLineStyle()) return kFALSE;
00733    if (GetLineWidth() != gStyle->GetLineWidth()) return kFALSE;
00734    return kTRUE;
00735 }
00736 
00737 //_____________________________________________________________________________
00738 Bool_t TGeoVolume::IsTopVolume() const
00739 {
00740 // True if this is the top volume of the geometry
00741    if (fGeoManager->GetTopVolume() == this) return kTRUE;
00742    return kFALSE;
00743 }
00744 
00745 //_____________________________________________________________________________
00746 Bool_t TGeoVolume::IsRaytracing() const
00747 {
00748 // Check if the painter is currently ray-tracing the content of this volume.
00749    return TGeoAtt::IsVisRaytrace();
00750 }
00751 
00752 //_____________________________________________________________________________
00753 void TGeoVolume::InspectMaterial() const
00754 {
00755 // Inspect the material for this volume.
00756    fMedium->GetMaterial()->Print();
00757 }
00758 
00759 //_____________________________________________________________________________
00760 TGeoVolume *TGeoVolume::Import(const char *filename, const char *name, Option_t * /*option*/)
00761 {
00762 // Import a volume from a file.
00763    if (!gGeoManager) gGeoManager = new TGeoManager("geometry","");
00764    if (!filename) return 0;
00765    TGeoVolume *volume = 0;
00766    if (strstr(filename,".gdml")) {
00767    // import from a gdml file
00768    } else {
00769    // import from a root file
00770       TFile *old = gFile;
00771       TFile *f = TFile::Open(filename);
00772       if (!f || f->IsZombie()) {
00773          if (old) old->cd();
00774          printf("Error: TGeoVolume::Import : Cannot open file %s\n", filename);
00775          return 0;
00776       }
00777       if (name && strlen(name) > 0) {
00778          volume = (TGeoVolume*)f->Get(name);
00779       } else {
00780          TIter next(f->GetListOfKeys());
00781          TKey *key;
00782          while ((key = (TKey*)next())) {
00783             if (strcmp(key->GetClassName(),"TGeoVolume") != 0) continue;
00784             volume = (TGeoVolume*)key->ReadObj();
00785             break;
00786          }
00787       }
00788       if (old) old->cd();
00789       delete f;         
00790    }
00791    if (!volume) return NULL;
00792    volume->RegisterYourself();
00793    return volume;
00794 }
00795    
00796 //_____________________________________________________________________________
00797 Int_t TGeoVolume::Export(const char *filename, const char *name, Option_t *option)
00798 {
00799 // Export this volume to a file.
00800    //
00801    // -Case 1: root file or root/xml file
00802    //  if filename end with ".root". The key will be named name
00803    //  if filename end with ".xml" a root/xml file is produced.
00804    //
00805    // -Case 2: C++ script
00806    //  if filename end with ".C"
00807    //
00808    // -Case 3: gdml file
00809    //  if filename end with ".gdml"
00810    //  NOTE that to use this option, the PYTHONPATH must be defined like
00811    //      export PYTHONPATH=$ROOTSYS/lib:$ROOTSYS/gdml
00812    //
00813    TString sfile(filename);
00814    if (sfile.Contains(".C")) {
00815       //Save volume as a C++ script
00816       Info("Export","Exporting volume %s as C++ code", GetName());
00817       SaveAs(filename, "");
00818       return 1;
00819    }
00820    if (sfile.Contains(".gdml")) {
00821      //Save geometry as a gdml file
00822       Info("Export","Exporting %s as gdml code - not implemented yet", GetName());
00823       return 0;
00824    }   
00825    if (sfile.Contains(".root") || sfile.Contains(".xml")) {  
00826       //Save volume in a root file
00827       Info("Export","Exporting %s as root file.", GetName());
00828       TString opt(option);
00829       if (!opt.Length()) opt = "recreate";
00830       TFile *f = TFile::Open(filename,opt.Data());
00831       if (!f || f->IsZombie()) {
00832          Error("Export","Cannot open file");
00833          return 0;
00834       } 
00835       TString keyname(name);
00836       if (keyname.IsNull()) keyname = GetName();
00837       Int_t nbytes = Write(keyname);
00838       delete f;
00839       return nbytes;
00840    }
00841    return 0;
00842 }
00843 
00844 //_____________________________________________________________________________
00845 void TGeoVolume::cd(Int_t inode) const
00846 {
00847 // Actualize matrix of node indexed <inode>
00848    if (fFinder) fFinder->cd(inode-fFinder->GetDivIndex());
00849 }
00850 
00851 //_____________________________________________________________________________
00852 void TGeoVolume::AddNode(const TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t * /*option*/)
00853 {
00854 // Add a TGeoNode to the list of nodes. This is the usual method for adding
00855 // daughters inside the container volume.
00856    TGeoMatrix *matrix = mat;
00857    if (matrix==0) matrix = gGeoIdentity;
00858    else           matrix->RegisterYourself();
00859    if (!vol) {
00860       Error("AddNode", "Volume is NULL");
00861       return;
00862    }
00863    if (!vol->IsValid()) {
00864       Error("AddNode", "Won't add node with invalid shape");
00865       printf("### invalid volume was : %s\n", vol->GetName());
00866       return;
00867    }   
00868    if (!fNodes) fNodes = new TObjArray();   
00869 
00870    if (fFinder) {
00871       // volume already divided.
00872       Error("AddNode", "Cannot add node %s_%i into divided volume %s", vol->GetName(), copy_no, GetName());
00873       return;
00874    }
00875 
00876    TGeoNodeMatrix *node = 0;
00877    node = new TGeoNodeMatrix(vol, matrix);
00878    node->SetMotherVolume(this);
00879    fNodes->Add(node);
00880    TString name = TString::Format("%s_%d", vol->GetName(), copy_no);
00881    if (fNodes->FindObject(name))
00882       Warning("AddNode", "Volume %s : added node %s with same name", GetName(), name.Data());
00883    node->SetName(name);
00884    node->SetNumber(copy_no);
00885 }
00886 
00887 //_____________________________________________________________________________
00888 void TGeoVolume::AddNodeOffset(const TGeoVolume *vol, Int_t copy_no, Double_t offset, Option_t * /*option*/)
00889 {
00890 // Add a division node to the list of nodes. The method is called by
00891 // TGeoVolume::Divide() for creating the division nodes.
00892    if (!vol) {
00893       Error("AddNodeOffset", "invalid volume");
00894       return;
00895    }
00896    if (!vol->IsValid()) {
00897       Error("AddNode", "Won't add node with invalid shape");
00898       printf("### invalid volume was : %s\n", vol->GetName());
00899       return;
00900    }   
00901    if (!fNodes) fNodes = new TObjArray();
00902    TGeoNode *node = new TGeoNodeOffset(vol, copy_no, offset);
00903    node->SetMotherVolume(this);
00904    fNodes->Add(node);
00905    TString name = TString::Format("%s_%d", vol->GetName(), copy_no+1);
00906    node->SetName(name);
00907    node->SetNumber(copy_no+1);
00908 }
00909 
00910 //_____________________________________________________________________________
00911 void TGeoVolume::AddNodeOverlap(const TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option)
00912 {
00913 // Add a TGeoNode to the list of nodes. This is the usual method for adding
00914 // daughters inside the container volume.
00915    if (!vol) {
00916       Error("AddNodeOverlap", "Volume is NULL");
00917       return;
00918    }
00919    if (!vol->IsValid()) {
00920       Error("AddNodeOverlap", "Won't add node with invalid shape");
00921       printf("### invalid volume was : %s\n", vol->GetName());
00922       return;
00923    }
00924    if (vol->IsAssembly()) {
00925       Warning("AddNodeOverlap", "Declaring assembly %s as possibly overlapping inside %s not allowed. Using AddNode instead !",vol->GetName(),GetName());
00926       AddNode(vol, copy_no, mat, option);
00927       return;
00928    }   
00929    TGeoMatrix *matrix = mat;
00930    if (matrix==0) matrix = gGeoIdentity;
00931    else           matrix->RegisterYourself();
00932    if (!fNodes) fNodes = new TObjArray();   
00933 
00934    if (fFinder) {
00935       // volume already divided.
00936       Error("AddNodeOverlap", "Cannot add node %s_%i into divided volume %s", vol->GetName(), copy_no, GetName());
00937       return;
00938    }
00939 
00940    TGeoNodeMatrix *node = new TGeoNodeMatrix(vol, matrix);
00941    node->SetMotherVolume(this);
00942    fNodes->Add(node);
00943    TString name = TString::Format("%s_%d", vol->GetName(), copy_no);
00944    if (fNodes->FindObject(name))
00945       Warning("AddNode", "Volume %s : added node %s with same name", GetName(), name.Data());
00946    node->SetName(name);
00947    node->SetNumber(copy_no);
00948    node->SetOverlapping();
00949    if (vol->GetMedium() == fMedium)
00950    node->SetVirtual();
00951 }
00952 
00953 //_____________________________________________________________________________
00954 TGeoVolume *TGeoVolume::Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed, Option_t *option)
00955 {
00956 // Division a la G3. The volume will be divided along IAXIS (see shape classes), in NDIV
00957 // slices, from START with given STEP. The division volumes will have medium number NUMED.
00958 // If NUMED=0 they will get the medium number of the divided volume (this). If NDIV<=0,
00959 // all range of IAXIS will be divided and the resulting number of divisions will be centered on
00960 // IAXIS. If STEP<=0, the real STEP will be computed as the full range of IAXIS divided by NDIV.
00961 // Options (case insensitive):
00962 //  N  - divide all range in NDIV cells (same effect as STEP<=0) (GSDVN in G3)
00963 //  NX - divide range starting with START in NDIV cells          (GSDVN2 in G3)
00964 //  S  - divide all range with given STEP. NDIV is computed and divisions will be centered
00965 //         in full range (same effect as NDIV<=0)                (GSDVS, GSDVT in G3)
00966 //  SX - same as DVS, but from START position.                   (GSDVS2, GSDVT2 in G3)
00967 
00968    if (fFinder) {
00969    // volume already divided.
00970       Fatal("Divide","volume %s already divided", GetName());
00971       return 0;
00972    }
00973    TString opt(option);
00974    opt.ToLower();
00975    TString stype = fShape->ClassName();
00976    if (!fNodes) fNodes = new TObjArray();
00977    Double_t xlo, xhi, range;
00978    range = fShape->GetAxisRange(iaxis, xlo, xhi);
00979    // for phi divisions correct the range
00980    if (!strcmp(fShape->GetAxisName(iaxis), "PHI")) {
00981       if ((start-xlo)<-1E-3) start+=360.;
00982       if (TGeoShape::IsSameWithinTolerance(range,360)) {
00983          xlo = start;
00984          xhi = start+range;
00985       }   
00986    }   
00987    if (range <=0) {
00988       InspectShape();
00989       Fatal("Divide", "cannot divide volume %s (%s) on %s axis", GetName(), stype.Data(), fShape->GetAxisName(iaxis));
00990       return 0;
00991    }
00992    if (ndiv<=0 || opt.Contains("s")) {
00993       if (step<=0) {
00994          Fatal("Divide", "invalid division type for volume %s : ndiv=%i, step=%g", GetName(), ndiv, step);
00995          return 0;
00996       }   
00997       if (opt.Contains("x")) {
00998          if ((xlo-start)>1E-3 || (xhi-start)<-1E-3) {
00999             Fatal("Divide", "invalid START=%g for division on axis %s of volume %s. Range is (%g, %g)",
01000                   start, fShape->GetAxisName(iaxis), GetName(), xlo, xhi);
01001             return 0;
01002          }
01003          xlo = start;
01004          range = xhi-xlo;
01005       }            
01006       ndiv = Int_t((range+0.1*step)/step);
01007       Double_t ddx = range - ndiv*step;
01008       // always center the division in this case
01009       if (ddx>1E-3) Warning("Divide", "division of volume %s on %s axis (ndiv=%d) will be centered in the full range",
01010                             GetName(), fShape->GetAxisName(iaxis), ndiv);
01011       start = xlo + 0.5*ddx;
01012    }
01013    if (step<=0 || opt.Contains("n")) {
01014       if (opt.Contains("x")) {
01015          if ((xlo-start)>1E-3 || (xhi-start)<-1E-3) {
01016             Fatal("Divide", "invalid START=%g for division on axis %s of volume %s. Range is (%g, %g)",
01017                   start, fShape->GetAxisName(iaxis), GetName(), xlo, xhi);
01018             return 0;
01019          }
01020          xlo = start;
01021          range = xhi-xlo;
01022       }     
01023       step  = range/ndiv;
01024       start = xlo;
01025    }
01026    
01027    Double_t end = start+ndiv*step;
01028    if (((start-xlo)<-1E-3) || ((end-xhi)>1E-3)) {
01029       Fatal("Divide", "division of volume %s on axis %s exceed range (%g, %g)",
01030             GetName(), fShape->GetAxisName(iaxis), xlo, xhi);
01031       return 0;
01032    }         
01033    TGeoVolume *voldiv = fShape->Divide(this, divname, iaxis, ndiv, start, step);
01034    if (numed) {
01035       TGeoMedium *medium = fGeoManager->GetMedium(numed);
01036       if (!medium) {
01037          Fatal("Divide", "invalid medium number %d for division volume %s", numed, divname);
01038          return voldiv;
01039       }   
01040       voldiv->SetMedium(medium);
01041       if (medium->GetMaterial()) medium->GetMaterial()->SetUsed();
01042    }   
01043    return voldiv; 
01044 }
01045 
01046 //_____________________________________________________________________________
01047 Int_t TGeoVolume::DistancetoPrimitive(Int_t px, Int_t py)
01048 {
01049 // compute the closest distance of approach from point px,py to this volume
01050    if (gGeoManager != fGeoManager) gGeoManager = fGeoManager;
01051    TVirtualGeoPainter *painter = fGeoManager->GetPainter();
01052    Int_t dist = 9999;
01053    if (!painter) return dist;
01054    dist = painter->DistanceToPrimitiveVol(this, px, py);
01055    return dist;
01056 }
01057 
01058 //_____________________________________________________________________________
01059 void TGeoVolume::Draw(Option_t *option)
01060 {
01061 // draw top volume according to option
01062    if (gGeoManager != fGeoManager) gGeoManager = fGeoManager;
01063    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
01064    TGeoAtt::SetVisRaytrace(kFALSE);
01065    if (!IsVisContainers()) SetVisLeaves();
01066    if (option && strlen(option) > 0) {
01067       painter->DrawVolume(this, option); 
01068    } else {
01069       painter->DrawVolume(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
01070    }  
01071 }
01072 
01073 //_____________________________________________________________________________
01074 void TGeoVolume::DrawOnly(Option_t *option)
01075 {
01076 // draw only this volume
01077    if (IsAssembly()) {
01078       Info("DrawOnly", "Volume assemblies do not support this option.");
01079       return;
01080    }   
01081    if (gGeoManager != fGeoManager) gGeoManager = fGeoManager;
01082    SetVisOnly();
01083    TGeoAtt::SetVisRaytrace(kFALSE);
01084    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
01085    if (option && strlen(option) > 0) {
01086       painter->DrawVolume(this, option); 
01087    } else {
01088       painter->DrawVolume(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
01089    }  
01090 }
01091 
01092 //_____________________________________________________________________________
01093 Bool_t TGeoVolume::OptimizeVoxels()
01094 {
01095 // Perform an exensive sampling to find which type of voxelization is
01096 // most efficient.
01097    printf("Optimizing volume %s ...\n", GetName());
01098    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
01099    return painter->TestVoxels(this);   
01100 }
01101 
01102 //_____________________________________________________________________________
01103 void TGeoVolume::Paint(Option_t *option)
01104 {
01105 // paint volume
01106    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
01107    painter->SetTopVolume(this);
01108 //   painter->Paint(option);   
01109    if (option && strlen(option) > 0) {
01110       painter->Paint(option); 
01111    } else {
01112       painter->Paint(gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
01113    }  
01114 }
01115 
01116 //_____________________________________________________________________________
01117 void TGeoVolume::PrintVoxels() const
01118 {
01119 // Print the voxels for this volume.
01120    if (fVoxels) fVoxels->Print();
01121 }
01122 
01123 //_____________________________________________________________________________
01124 void TGeoVolume::ReplayCreation(const TGeoVolume *other)
01125 {
01126 // Recreate the content of the other volume without pointer copying. Voxels are 
01127 // ignored and supposed to be created in a later step via Voxelize.
01128    Int_t nd = other->GetNdaughters();
01129    if (!nd) return;
01130    TGeoPatternFinder *finder = other->GetFinder();
01131    if (finder) {
01132       Int_t iaxis = finder->GetDivAxis();
01133       Int_t ndiv = finder->GetNdiv();
01134       Double_t start = finder->GetStart();
01135       Double_t step = finder->GetStep();
01136       Int_t numed = other->GetNode(0)->GetVolume()->GetMedium()->GetId();
01137       TGeoVolume *voldiv = Divide(other->GetNode(0)->GetVolume()->GetName(), iaxis, ndiv, start, step, numed);
01138       voldiv->ReplayCreation(other->GetNode(0)->GetVolume());
01139       return;
01140    }   
01141    for (Int_t i=0; i<nd; i++) {
01142       TGeoNode *node = other->GetNode(i);
01143       if (node->IsOverlapping()) AddNodeOverlap(node->GetVolume(), node->GetNumber(), node->GetMatrix());
01144       else AddNode(node->GetVolume(), node->GetNumber(), node->GetMatrix());
01145    }
01146 }      
01147    
01148 //_____________________________________________________________________________
01149 void TGeoVolume::PrintNodes() const
01150 {
01151 // print nodes
01152    Int_t nd = GetNdaughters();
01153    for (Int_t i=0; i<nd; i++) {
01154       printf("%s\n", GetNode(i)->GetName());
01155       cd(i);
01156       GetNode(i)->GetMatrix()->Print();
01157    }   
01158 }
01159 //______________________________________________________________________________
01160 TH2F *TGeoVolume::LegoPlot(Int_t ntheta, Double_t themin, Double_t themax,
01161                             Int_t nphi,   Double_t phimin, Double_t phimax,
01162                             Double_t rmin, Double_t rmax, Option_t *option)
01163 {
01164 // Generate a lego plot fot the top volume, according to option.
01165    TVirtualGeoPainter *p = fGeoManager->GetGeomPainter();
01166    TGeoVolume *old_vol = fGeoManager->GetTopVolume();
01167    if (old_vol!=this) fGeoManager->SetTopVolume(this);
01168    else old_vol=0;
01169    TH2F *hist = p->LegoPlot(ntheta, themin, themax, nphi, phimin, phimax, rmin, rmax, option);   
01170    hist->Draw("lego1sph");
01171    return hist;
01172 }
01173 
01174 //_____________________________________________________________________________
01175 void TGeoVolume::RegisterYourself(Option_t *option)
01176 {
01177 // Register the volume and all materials/media/matrices/shapes to the manager.
01178    if (fGeoManager->GetListOfVolumes()->FindObject(this)) return;
01179    // Register volume
01180    fGeoManager->AddVolume(this);
01181    // Register shape
01182    if (!fGeoManager->GetListOfShapes()->FindObject(fShape)) {
01183       if (fShape->IsComposite()) {
01184          TGeoCompositeShape *comp = (TGeoCompositeShape*)fShape;
01185          comp->RegisterYourself();
01186       } else {
01187          fGeoManager->AddShape(fShape);   
01188       }
01189    }   
01190    // Register medium/material
01191    if (fMedium && !fGeoManager->GetListOfMedia()->FindObject(fMedium)) {
01192       fGeoManager->GetListOfMedia()->Add(fMedium);
01193       if (!fGeoManager->GetListOfMaterials()->FindObject(fMedium->GetMaterial()))
01194          fGeoManager->AddMaterial(fMedium->GetMaterial());
01195    }
01196    // Register matrices for nodes.
01197    TGeoMatrix *matrix;
01198    TGeoNode *node;
01199    Int_t nd = GetNdaughters();
01200    Int_t i;
01201    for (i=0; i<nd; i++) {
01202       node = GetNode(i);
01203       matrix = node->GetMatrix();
01204       if (!matrix->IsRegistered()) matrix->RegisterYourself();
01205       else if (!fGeoManager->GetListOfMatrices()->FindObject(matrix)) {
01206          fGeoManager->GetListOfMatrices()->Add(matrix);
01207       }
01208    }
01209    // Call RegisterYourself recursively
01210    for (i=0; i<nd; i++) GetNode(i)->GetVolume()->RegisterYourself(option);
01211 }      
01212       
01213 //_____________________________________________________________________________
01214 void TGeoVolume::RandomPoints(Int_t npoints, Option_t *option)
01215 {
01216 // Draw random points in the bounding box of this volume.
01217    if (gGeoManager != fGeoManager) gGeoManager = fGeoManager;
01218    TGeoVolume *old_vol = fGeoManager->GetTopVolume();
01219    if (old_vol!=this) fGeoManager->SetTopVolume(this);
01220    else old_vol=0;
01221    fGeoManager->RandomPoints(this, npoints, option);
01222    if (old_vol) fGeoManager->SetTopVolume(old_vol);
01223 }
01224 
01225 //_____________________________________________________________________________
01226 void TGeoVolume::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz)
01227 {
01228 // Random raytracing method.
01229    if (gGeoManager != fGeoManager) gGeoManager = fGeoManager;
01230    TGeoVolume *old_vol = fGeoManager->GetTopVolume();
01231    if (old_vol!=this) fGeoManager->SetTopVolume(this);
01232    else old_vol=0;
01233    fGeoManager->RandomRays(nrays, startx, starty, startz);
01234    if (old_vol) fGeoManager->SetTopVolume(old_vol);
01235 }
01236 
01237 //_____________________________________________________________________________
01238 void TGeoVolume::Raytrace(Bool_t flag)
01239 {
01240 // Draw this volume with current settings and perform raytracing in the pad.
01241    TGeoAtt::SetVisRaytrace(kFALSE);
01242    if (gGeoManager != fGeoManager) gGeoManager = fGeoManager;
01243    TVirtualGeoPainter *painter = fGeoManager->GetGeomPainter();
01244    Bool_t drawn = (painter->GetDrawnVolume()==this)?kTRUE:kFALSE;   
01245    if (!drawn) {
01246       painter->DrawVolume(this, "");
01247       TGeoAtt::SetVisRaytrace(flag);
01248       painter->ModifiedPad();
01249       return;
01250    }   
01251    TGeoAtt::SetVisRaytrace(flag);
01252    painter->ModifiedPad();
01253 }   
01254 
01255 //______________________________________________________________________________
01256 void TGeoVolume::SaveAs(const char *filename, Option_t *option) const
01257 {
01258 //  Save geometry having this as top volume as a C++ macro.
01259    if (!filename) return;
01260    ofstream out;
01261    out.open(filename, ios::out);
01262    if (out.bad()) {
01263       Error("SavePrimitive", "Bad file name: %s", filename);
01264       return;
01265    }
01266    if (fGeoManager->GetTopVolume() != this) fGeoManager->SetTopVolume((TGeoVolume*)this);
01267    
01268    TString fname(filename);
01269    Int_t ind = fname.Index(".");
01270    if (ind>0) fname.Remove(ind);
01271    out << "void "<<fname<<"() {" << endl;
01272    out << "   gSystem->Load(\"libGeom\");" << endl;
01273    ((TGeoVolume*)this)->SavePrimitive(out,option);
01274    out << "}" << endl;
01275 }   
01276 
01277 //______________________________________________________________________________
01278 void TGeoVolume::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
01279 {
01280    // Save a primitive as a C++ statement(s) on output stream "out".
01281    out.precision(6);
01282    out.setf(ios::fixed);
01283    Int_t i,icopy;
01284    Int_t nd = GetNdaughters();
01285    TGeoVolume *dvol;
01286    TGeoNode *dnode;
01287    TGeoMatrix *matrix;
01288 
01289    // check if we need to save shape/volume
01290    Bool_t mustDraw = kFALSE;
01291    if (fGeoManager->GetGeomPainter()->GetTopVolume()==this) mustDraw = kTRUE;
01292    if (!strlen(option)) {
01293       fGeoManager->SetAllIndex();
01294       out << "   new TGeoManager(\"" << fGeoManager->GetName() << "\", \"" << fGeoManager->GetTitle() << "\");" << endl << endl;
01295 //      if (mustDraw) out << "   Bool_t mustDraw = kTRUE;" << endl;
01296 //      else          out << "   Bool_t mustDraw = kFALSE;" << endl;
01297       out << "   Double_t dx,dy,dz;" << endl;
01298       out << "   Double_t dx1, dx2, dy1, dy2;" << endl;
01299       out << "   Double_t vert[20], par[20];" << endl;
01300       out << "   Double_t theta, phi, h1, bl1, tl1, alpha1, h2, bl2, tl2, alpha2;" << endl;
01301       out << "   Double_t twist;" << endl;
01302       out << "   Double_t origin[3];" << endl;
01303       out << "   Double_t rmin, rmax, rmin1, rmax1, rmin2, rmax2;" << endl;
01304       out << "   Double_t r, rlo, rhi;" << endl;
01305       out << "   Double_t phi1, phi2;" << endl;
01306       out << "   Double_t a,b;" << endl;
01307       out << "   Double_t point[3], norm[3];" << endl;
01308       out << "   Double_t rin, stin, rout, stout;" << endl;
01309       out << "   Double_t thx, phx, thy, phy, thz, phz;" << endl;
01310       out << "   Double_t alpha, theta1, theta2, phi1, phi2, dphi;" << endl;
01311       out << "   Double_t tr[3], rot[9];" << endl;
01312       out << "   Double_t z, density, radl, absl, w;" << endl;
01313       out << "   Double_t lx,ly,lz,tx,ty,tz;" << endl;
01314       out << "   Double_t xvert[50], yvert[50];" << endl;
01315       out << "   Double_t zsect,x0,y0,scale0;" << endl;
01316       out << "   Int_t nel, numed, nz, nedges, nvert;" << endl;
01317       out << "   TGeoBoolNode *pBoolNode = 0;" << endl << endl;
01318       // first save materials/media
01319       out << "   // MATERIALS, MIXTURES AND TRACKING MEDIA" << endl;
01320       SavePrimitive(out, "m");
01321       // then, save matrices
01322       out << endl << "   // TRANSFORMATION MATRICES" << endl;
01323       SavePrimitive(out, "x");
01324       // save this volume and shape
01325       SavePrimitive(out, "s");
01326       out << endl << "   // SET TOP VOLUME OF GEOMETRY" << endl;
01327       out << "   gGeoManager->SetTopVolume(" << GetPointerName() << ");" << endl;
01328       // save daughters
01329       out << endl << "   // SHAPES, VOLUMES AND GEOMETRICAL HIERARCHY" << endl;
01330       SavePrimitive(out, "d");
01331       out << endl << "   // CLOSE GEOMETRY" << endl;
01332       out << "   gGeoManager->CloseGeometry();" << endl;
01333       if (mustDraw) {
01334          if (!IsRaytracing()) out << "   gGeoManager->GetTopVolume()->Draw();" << endl;
01335          else                 out << "   gGeoManager->GetTopVolume()->Raytrace();" << endl;
01336       }
01337       return;
01338    }
01339    // check if we need to save shape/volume
01340    if (!strcmp(option, "s")) {
01341       // create the shape for this volume
01342       if (TestAttBit(TGeoAtt::kSavePrimitiveAtt)) return;
01343       if (!IsAssembly()) {
01344          fShape->SavePrimitive(out,option);      
01345          out << "   // Volume: " << GetName() << endl;
01346          out << "   " << GetPointerName() << " = new TGeoVolume(\"" << GetName() << "\"," << fShape->GetPointerName() << ", "<< fMedium->GetPointerName() << ");" << endl;
01347       } else {
01348          out << "   // Assembly: " << GetName() << endl;
01349          out << "   " << GetPointerName() << " = new TGeoVolumeAssembly(\"" << GetName() << "\"" << ");" << endl;
01350       }           
01351       if (fLineColor != 1) out << "   " << GetPointerName() << "->SetLineColor(" << fLineColor << ");" << endl;
01352       if (fLineWidth != 1) out << "   " << GetPointerName() << "->SetLineWidth(" << fLineWidth << ");" << endl;
01353       if (fLineStyle != 1) out << "   " << GetPointerName() << "->SetLineStyle(" << fLineStyle << ");" << endl;
01354       if (!IsVisible() && !IsAssembly()) out << "   " << GetPointerName() << "->SetVisibility(kFALSE);" << endl;
01355       if (!IsVisibleDaughters()) out << "   " << GetPointerName() << "->VisibleDaughters(kFALSE);" << endl;
01356       if (IsVisContainers()) out << "   " << GetPointerName() << "->SetVisContainers(kTRUE);" << endl;
01357       if (IsVisLeaves()) out << "   " << GetPointerName() << "->SetVisLeaves(kTRUE);" << endl;
01358       SetAttBit(TGeoAtt::kSavePrimitiveAtt);
01359    }   
01360    // check if we need to save the media
01361    if (!strcmp(option, "m")) {
01362       if (fMedium) fMedium->SavePrimitive(out,option);
01363       for (i=0; i<nd; i++) {
01364          dvol = GetNode(i)->GetVolume();
01365          dvol->SavePrimitive(out,option);
01366       }
01367       return;      
01368    }   
01369    // check if we need to save the matrices
01370    if (!strcmp(option, "x")) {
01371       if (fFinder) {
01372          dvol = GetNode(0)->GetVolume();
01373          dvol->SavePrimitive(out,option);
01374          return;
01375       }
01376       for (i=0; i<nd; i++) {
01377          dnode = GetNode(i);
01378          matrix = dnode->GetMatrix();
01379          if (!matrix->IsIdentity()) matrix->SavePrimitive(out,option);
01380          dnode->GetVolume()->SavePrimitive(out,option);
01381       }
01382       return;      
01383    } 
01384    // check if we need to save volume daughters
01385    if (!strcmp(option, "d")) {
01386       if (!nd) return;
01387       if (TestAttBit(TGeoAtt::kSaveNodesAtt)) return;
01388       SetAttBit(TGeoAtt::kSaveNodesAtt);     
01389       if (fFinder) {
01390          // volume divided: generate volume->Divide()
01391          dnode = GetNode(0);
01392          dvol = dnode->GetVolume();
01393          out << "   TGeoVolume *" << dvol->GetPointerName() << " = ";
01394          out << GetPointerName() << "->Divide(\"" << dvol->GetName() << "\", ";
01395          fFinder->SavePrimitive(out,option);
01396          if (fMedium != dvol->GetMedium()) {
01397             out << ", " << dvol->GetMedium()->GetId();
01398          }
01399          out << ");" << endl;   
01400          dvol->SavePrimitive(out,"d");   
01401          return;
01402       }
01403       for (i=0; i<nd; i++) {
01404          dnode = GetNode(i);
01405          dvol = dnode->GetVolume();
01406          dvol->SavePrimitive(out,"s");
01407          matrix = dnode->GetMatrix();
01408          icopy = dnode->GetNumber();
01409          // generate AddNode()
01410          out << "   " << GetPointerName() << "->AddNode";
01411          if (dnode->IsOverlapping()) out << "Overlap";
01412          out << "(" << dvol->GetPointerName() << ", " << icopy;
01413          if (!matrix->IsIdentity()) out << ", " << matrix->GetPointerName();
01414          out << ");" << endl;
01415       }
01416       // Recursive loop to daughters
01417       for (i=0; i<nd; i++) {
01418          dnode = GetNode(i);
01419          dvol = dnode->GetVolume();
01420          dvol->SavePrimitive(out,"d");
01421       } 
01422    }   
01423 }
01424 
01425 //_____________________________________________________________________________
01426 void TGeoVolume::UnmarkSaved()
01427 {
01428 // Reset SavePrimitive bits.
01429    ResetAttBit(TGeoAtt::kSavePrimitiveAtt);
01430    ResetAttBit(TGeoAtt::kSaveNodesAtt);
01431    if (fShape) fShape->ResetBit(TGeoShape::kGeoSavePrimitive);
01432 }   
01433 
01434 //_____________________________________________________________________________
01435 void TGeoVolume::ExecuteEvent(Int_t event, Int_t px, Int_t py)
01436 {
01437 // Execute mouse actions on this volume.
01438    TVirtualGeoPainter *painter = fGeoManager->GetPainter();
01439    if (!painter) return;
01440    painter->ExecuteVolumeEvent(this, event, px, py);
01441 }
01442 
01443 //_____________________________________________________________________________
01444 TGeoNode *TGeoVolume::FindNode(const char *name) const
01445 {
01446 // search a daughter inside the list of nodes
01447    return ((TGeoNode*)fNodes->FindObject(name));
01448 }
01449 
01450 //_____________________________________________________________________________
01451 Int_t TGeoVolume::GetNodeIndex(const TGeoNode *node, Int_t *check_list, Int_t ncheck) const
01452 {
01453 // Get the index of a daugther within check_list by providing the node pointer.
01454    TGeoNode *current = 0;
01455    for (Int_t i=0; i<ncheck; i++) {
01456       current = (TGeoNode*)fNodes->At(check_list[i]);
01457       if (current==node) return check_list[i];
01458    }
01459    return -1;
01460 }
01461 
01462 //_____________________________________________________________________________
01463 Int_t TGeoVolume::GetIndex(const TGeoNode *node) const
01464 {
01465 // get index number for a given daughter
01466    TGeoNode *current = 0;
01467    Int_t nd = GetNdaughters();
01468    if (!nd) return -1;
01469    for (Int_t i=0; i<nd; i++) {
01470       current = (TGeoNode*)fNodes->At(i);
01471       if (current==node) return i;
01472    }
01473    return -1;
01474 }
01475 
01476 //_____________________________________________________________________________
01477 char *TGeoVolume::GetObjectInfo(Int_t px, Int_t py) const
01478 {
01479 // Get volume info for the browser.
01480    TGeoVolume *vol = (TGeoVolume*)this;
01481    TVirtualGeoPainter *painter = fGeoManager->GetPainter();
01482    if (!painter) return 0;
01483    return (char*)painter->GetVolumeInfo(vol, px, py);
01484 }
01485 
01486 //_____________________________________________________________________________
01487 Bool_t TGeoVolume::GetOptimalVoxels() const
01488 {
01489 //--- Returns true if cylindrical voxelization is optimal.
01490    Int_t nd = GetNdaughters();
01491    if (!nd) return kFALSE;
01492    Int_t id;
01493    Int_t ncyl = 0;
01494    TGeoNode *node;
01495    for (id=0; id<nd; id++) {
01496       node = (TGeoNode*)fNodes->At(id);
01497       ncyl += node->GetOptimalVoxels();
01498    }
01499    if (ncyl>(nd/2)) return kTRUE;
01500    return kFALSE;
01501 }      
01502 
01503 //_____________________________________________________________________________
01504 char *TGeoVolume::GetPointerName() const
01505 {
01506 // Provide a pointer name containing uid.
01507    static TString name;
01508    name = TString::Format("p%s_%lx", GetName(), (ULong_t)this);
01509    return (char*)name.Data();
01510 }
01511 
01512 //_____________________________________________________________________________
01513 TGeoVoxelFinder *TGeoVolume::GetVoxels() const
01514 {
01515 // Getter for optimization structure.
01516    if (fVoxels && !fVoxels->IsInvalid()) return fVoxels;
01517    return NULL;
01518 }   
01519 
01520 //_____________________________________________________________________________
01521 void TGeoVolume::GrabFocus()
01522 {
01523 // Move perspective view focus to this volume
01524    TVirtualGeoPainter *painter = fGeoManager->GetPainter();
01525    if (painter) painter->GrabFocus();
01526 }   
01527 
01528 //_____________________________________________________________________________
01529 TGeoVolume *TGeoVolume::CloneVolume() const
01530 {
01531 // Clone this volume.
01532    // build a volume with same name, shape and medium
01533    TGeoVolume *vol = new TGeoVolume(GetName(), fShape, fMedium);
01534    Int_t i;
01535    // copy volume attributes
01536    vol->SetLineColor(GetLineColor());
01537    vol->SetLineStyle(GetLineStyle());
01538    vol->SetLineWidth(GetLineWidth());
01539    vol->SetFillColor(GetFillColor());
01540    vol->SetFillStyle(GetFillStyle());
01541    // copy other attributes
01542    Int_t nbits = 8*sizeof(UInt_t);
01543    for (i=0; i<nbits; i++) 
01544       vol->SetAttBit(1<<i, TGeoAtt::TestAttBit(1<<i));
01545    for (i=14; i<24; i++)
01546       vol->SetBit(1<<i, TestBit(1<<i));   
01547    
01548    // copy field
01549    vol->SetField(fField);
01550    // Set bits
01551    for (i=0; i<nbits; i++) 
01552       vol->SetBit(1<<i, TObject::TestBit(1<<i));
01553    vol->SetBit(kVolumeClone);   
01554    // copy nodes
01555 //   CloneNodesAndConnect(vol);
01556    vol->MakeCopyNodes(this);   
01557    // if volume is divided, copy finder
01558    vol->SetFinder(fFinder);
01559    // copy voxels
01560    TGeoVoxelFinder *voxels = 0;
01561    if (fVoxels) {
01562       voxels = new TGeoVoxelFinder(vol);
01563       vol->SetVoxelFinder(voxels);
01564    }   
01565    // copy option, uid
01566    vol->SetOption(fOption);
01567    vol->SetNumber(fNumber);
01568    vol->SetNtotal(fNtotal);
01569    return vol;
01570 }
01571 
01572 //_____________________________________________________________________________
01573 void TGeoVolume::CloneNodesAndConnect(TGeoVolume *newmother) const
01574 {
01575 // Clone the array of nodes.
01576    if (!fNodes) return;
01577    TGeoNode *node;
01578    Int_t nd = fNodes->GetEntriesFast();
01579    if (!nd) return;
01580    // create new list of nodes
01581    TObjArray *list = new TObjArray(nd);
01582    // attach it to new volume
01583    newmother->SetNodes(list);
01584 //   ((TObject*)newmother)->SetBit(kVolumeImportNodes);
01585    for (Int_t i=0; i<nd; i++) {
01586       //create copies of nodes and add them to list
01587       node = GetNode(i)->MakeCopyNode();
01588       node->SetMotherVolume(newmother);
01589       list->Add(node);
01590    }
01591 }
01592 
01593 //_____________________________________________________________________________
01594 void TGeoVolume::MakeCopyNodes(const TGeoVolume *other)
01595 {
01596 // make a new list of nodes and copy all nodes of other volume inside
01597    Int_t nd = other->GetNdaughters();
01598    if (!nd) return;
01599    if (fNodes) {
01600       if (!TObject::TestBit(kVolumeImportNodes)) fNodes->Delete();
01601       delete fNodes;   
01602    }   
01603    fNodes = new TObjArray();
01604    for (Int_t i=0; i<nd; i++) fNodes->Add(other->GetNode(i));
01605    TObject::SetBit(kVolumeImportNodes);
01606 }      
01607 
01608 //_____________________________________________________________________________
01609 TGeoVolume *TGeoVolume::MakeCopyVolume(TGeoShape *newshape)
01610 {
01611     // make a copy of this volume
01612    // build a volume with same name, shape and medium
01613    TGeoVolume *vol = new TGeoVolume(GetName(), newshape, fMedium);
01614    // copy volume attributes
01615    vol->SetVisibility(IsVisible());
01616    vol->SetLineColor(GetLineColor());
01617    vol->SetLineStyle(GetLineStyle());
01618    vol->SetLineWidth(GetLineWidth());
01619    vol->SetFillColor(GetFillColor());
01620    vol->SetFillStyle(GetFillStyle());
01621    // copy field
01622    vol->SetField(fField);
01623    // if divided, copy division object
01624    if (fFinder) {
01625 //       Error("MakeCopyVolume", "volume %s divided", GetName());
01626       vol->SetFinder(fFinder);
01627    }   
01628    CloneNodesAndConnect(vol);
01629 //   ((TObject*)vol)->SetBit(kVolumeImportNodes);
01630    ((TObject*)vol)->SetBit(kVolumeClone);
01631    return vol;       
01632 }    
01633 
01634 //_____________________________________________________________________________
01635 TGeoVolume *TGeoVolume::MakeReflectedVolume(const char *newname) const
01636 {
01637 // Make a copy of this volume which is reflected with respect to XY plane.
01638    static TMap map(100);
01639    if (!fGeoManager->IsClosed()) {
01640       Error("MakeReflectedVolume", "Geometry must be closed.");
01641       return NULL;
01642    }   
01643    TGeoVolume *vol = (TGeoVolume*)map.GetValue(this);
01644    if (vol) {
01645       if (strlen(newname)) vol->SetName(newname);
01646       return vol;
01647    }
01648 //   printf("Making reflection for volume: %s\n", GetName());   
01649    vol = CloneVolume();
01650    map.Add((TObject*)this, vol);
01651    if (strlen(newname)) vol->SetName(newname);
01652    delete vol->GetNodes();
01653    vol->SetNodes(NULL);
01654    vol->SetBit(kVolumeImportNodes, kFALSE);
01655    CloneNodesAndConnect(vol);
01656    // The volume is now properly cloned, but with the same shape.
01657    // Reflect the shape (if any) and connect it.
01658    if (fShape) {
01659       TGeoShape *reflected_shape = 
01660          TGeoScaledShape::MakeScaledShape("", fShape, new TGeoScale(1.,1.,-1.));
01661       vol->SetShape(reflected_shape);
01662    }   
01663    // Reflect the daughters.
01664    Int_t nd = vol->GetNdaughters();
01665    if (!nd) return vol;
01666    TGeoNodeMatrix *node;
01667    TGeoMatrix *local, *local_cloned;
01668    TGeoVolume *new_vol;
01669    if (!vol->GetFinder()) {
01670       for (Int_t i=0; i<nd; i++) {
01671          node = (TGeoNodeMatrix*)vol->GetNode(i);
01672          local = node->GetMatrix();
01673 //         printf("%s before\n", node->GetName());
01674 //         local->Print();
01675          Bool_t reflected = local->IsReflection();
01676          local_cloned = new TGeoCombiTrans(*local);
01677          local_cloned->RegisterYourself();
01678          node->SetMatrix(local_cloned);
01679          if (!reflected) {
01680          // We need to reflect only the translation and propagate to daughters.
01681             // H' = Sz * H * Sz
01682             local_cloned->ReflectZ(kTRUE);
01683             local_cloned->ReflectZ(kFALSE);
01684 //            printf("%s after\n", node->GetName());
01685 //            node->GetMatrix()->Print();
01686             new_vol = node->GetVolume()->MakeReflectedVolume();
01687             node->SetVolume(new_vol);
01688             continue;
01689          }
01690          // The next daughter is already reflected, so reflect on Z everything and stop
01691          local_cloned->ReflectZ(kTRUE); // rot + tr
01692 //         printf("%s already reflected... After:\n", node->GetName());
01693 //         node->GetMatrix()->Print();
01694       }
01695       if (vol->GetVoxels()) vol->GetVoxels()->Voxelize();
01696       return vol;
01697    }
01698    // Volume is divided, so we have to reflect the division.
01699 //   printf("   ... divided %s\n", fFinder->ClassName());
01700    TGeoPatternFinder *new_finder = fFinder->MakeCopy(kTRUE);
01701    new_finder->SetVolume(vol);
01702    vol->SetFinder(new_finder);
01703    TGeoNodeOffset *nodeoff;
01704    new_vol = 0;
01705    for (Int_t i=0; i<nd; i++) {
01706       nodeoff = (TGeoNodeOffset*)vol->GetNode(i);
01707       nodeoff->SetFinder(new_finder);
01708       new_vol = nodeoff->GetVolume()->MakeReflectedVolume();
01709       nodeoff->SetVolume(new_vol); 
01710    }   
01711    return vol;
01712 }
01713    
01714 //_____________________________________________________________________________
01715 void TGeoVolume::SetAsTopVolume()
01716 {
01717 // Set this volume as the TOP one (the whole geometry starts from here)
01718    fGeoManager->SetTopVolume(this);
01719 }
01720 
01721 //_____________________________________________________________________________
01722 void TGeoVolume::SetCurrentPoint(Double_t x, Double_t y, Double_t z)
01723 {
01724 // Set the current tracking point.
01725    fGeoManager->SetCurrentPoint(x,y,z);
01726 }
01727 
01728 //_____________________________________________________________________________
01729 void TGeoVolume::SetShape(const TGeoShape *shape)
01730 {
01731 // set the shape associated with this volume
01732    if (!shape) {
01733       Error("SetShape", "No shape");
01734       return;
01735    }
01736    fShape = (TGeoShape*)shape;  
01737 }
01738 
01739 //_____________________________________________________________________________
01740 void TGeoVolume::SortNodes()
01741 {
01742 // sort nodes by decreasing volume of the bounding box. ONLY nodes comes first,
01743 // then overlapping nodes and finally division nodes.
01744    if (!Valid()) {
01745       Error("SortNodes", "Bounding box not valid");
01746       return;
01747    }
01748    Int_t nd = GetNdaughters();
01749 //   printf("volume : %s, nd=%i\n", GetName(), nd);
01750    if (!nd) return;
01751    if (fFinder) return;
01752 //   printf("Nodes for %s\n", GetName());
01753    Int_t id = 0;
01754    TGeoNode *node = 0;
01755    TObjArray *nodes = new TObjArray(nd);
01756    Int_t inode = 0;
01757    // first put ONLY's
01758    for (id=0; id<nd; id++) {
01759       node = GetNode(id);
01760       if (node->InheritsFrom(TGeoNodeOffset::Class()) || node->IsOverlapping()) continue;
01761       nodes->Add(node);
01762 //      printf("inode %i ONLY\n", inode);
01763       inode++;
01764    }
01765    // second put overlapping nodes
01766    for (id=0; id<nd; id++) {
01767       node = GetNode(id);
01768       if (node->InheritsFrom(TGeoNodeOffset::Class()) || (!node->IsOverlapping())) continue;
01769       nodes->Add(node);
01770 //      printf("inode %i MANY\n", inode);
01771       inode++;
01772    }
01773    // third put the divided nodes
01774    if (fFinder) {
01775       fFinder->SetDivIndex(inode);
01776       for (id=0; id<nd; id++) {
01777          node = GetNode(id);
01778          if (!node->InheritsFrom(TGeoNodeOffset::Class())) continue;
01779          nodes->Add(node);
01780 //         printf("inode %i DIV\n", inode);
01781          inode++;
01782       }
01783    }
01784    if (inode != nd) printf(" volume %s : number of nodes does not match!!!\n", GetName());
01785    delete fNodes;
01786    fNodes = nodes;
01787 }
01788 
01789 //_____________________________________________________________________________
01790 void TGeoVolume::Streamer(TBuffer &R__b)
01791 {
01792    // Stream an object of class TGeoVolume.
01793    if (R__b.IsReading()) {
01794       R__b.ReadClassBuffer(TGeoVolume::Class(), this);
01795       if (fVoxels && fVoxels->IsInvalid()) Voxelize("");
01796    } else {
01797       if (!fVoxels) {
01798          R__b.WriteClassBuffer(TGeoVolume::Class(), this);
01799       } else {
01800          if (!fGeoManager->IsStreamingVoxels()) {
01801             TGeoVoxelFinder *voxels = fVoxels;
01802             fVoxels = 0;
01803             R__b.WriteClassBuffer(TGeoVolume::Class(), this);
01804             fVoxels = voxels;
01805          } else {
01806             R__b.WriteClassBuffer(TGeoVolume::Class(), this);
01807          }
01808       }
01809    }
01810 }
01811 
01812 //_____________________________________________________________________________
01813 void TGeoVolume::SetOption(const char * /*option*/)
01814 {
01815 // Set the current options (none implemented)
01816 }
01817 
01818 //_____________________________________________________________________________
01819 void TGeoVolume::SetLineColor(Color_t lcolor) 
01820 {
01821 // Set the line color.
01822    TAttLine::SetLineColor(lcolor);
01823 }   
01824 
01825 //_____________________________________________________________________________
01826 void TGeoVolume::SetLineStyle(Style_t lstyle) 
01827 {
01828 // Set the line style.
01829    TAttLine::SetLineStyle(lstyle);
01830 }   
01831 
01832 //_____________________________________________________________________________
01833 void TGeoVolume::SetLineWidth(Style_t lwidth) 
01834 {
01835 // Set the line width.
01836    TAttLine::SetLineWidth(lwidth);
01837 }   
01838 
01839 //_____________________________________________________________________________
01840 TGeoNode *TGeoVolume::GetNode(const char *name) const
01841 {
01842 // get the pointer to a daughter node
01843    if (!fNodes) return 0;
01844    TGeoNode *node = (TGeoNode *)fNodes->FindObject(name);
01845    return node;
01846 }
01847 
01848 //_____________________________________________________________________________
01849 Int_t TGeoVolume::GetByteCount() const
01850 {
01851 // get the total size in bytes for this volume
01852    Int_t count = 28+2+6+4+0;    // TNamed+TGeoAtt+TAttLine+TAttFill+TAtt3D
01853    count += strlen(GetName()) + strlen(GetTitle()); // name+title
01854    count += 4+4+4+4+4; // fShape + fMedium + fFinder + fField + fNodes
01855    count += 8 + strlen(fOption.Data()); // fOption
01856    if (fShape)  count += fShape->GetByteCount();
01857    if (fFinder) count += fFinder->GetByteCount();
01858    if (fNodes) {
01859       count += 32 + 4*fNodes->GetEntries(); // TObjArray
01860       TIter next(fNodes);
01861       TGeoNode *node;
01862       while ((node=(TGeoNode*)next())) count += node->GetByteCount();
01863    }
01864    return count;
01865 }
01866 
01867 //_____________________________________________________________________________
01868 void TGeoVolume::FindOverlaps() const
01869 {
01870 // loop all nodes marked as overlaps and find overlaping brothers
01871    if (!Valid()) {
01872       Error("FindOverlaps","Bounding box not valid");
01873       return;
01874    }   
01875    if (!fVoxels) return;
01876    Int_t nd = GetNdaughters();
01877    if (!nd) return;
01878    TGeoNode *node=0;
01879    Int_t inode = 0;
01880    for (inode=0; inode<nd; inode++) {
01881       node = GetNode(inode);
01882       if (!node->IsOverlapping()) continue;
01883       fVoxels->FindOverlaps(inode);
01884    }
01885 }
01886 
01887 //_____________________________________________________________________________
01888 void TGeoVolume::RemoveNode(TGeoNode *node) 
01889 {
01890 // Remove an existing daughter.
01891    if (!fNodes || !fNodes->GetEntriesFast()) return;
01892    if (!fNodes->Remove(node)) return;
01893    fNodes->Compress();
01894    if (fVoxels) fVoxels->SetNeedRebuild();
01895    if (IsAssembly()) fShape->ComputeBBox();
01896 }   
01897 
01898 //_____________________________________________________________________________
01899 TGeoNode *TGeoVolume::ReplaceNode(TGeoNode *nodeorig, TGeoShape *newshape, TGeoMatrix *newpos, TGeoMedium *newmed) 
01900 {
01901 // Replace an existing daughter with a new volume having the same name but
01902 // possibly a new shape, position or medium. Not allowed for positioned assemblies.
01903 // For division cells, the new shape/matrix are ignored.
01904    Int_t ind = GetIndex(nodeorig);
01905    if (ind < 0) return NULL;
01906    TGeoVolume *oldvol = nodeorig->GetVolume();
01907    if (oldvol->IsAssembly()) {
01908       Error("ReplaceNode", "Cannot replace node %s since it is an assembly", nodeorig->GetName());
01909       return NULL;
01910    }   
01911    TGeoShape  *shape = oldvol->GetShape();
01912    if (newshape && !nodeorig->IsOffset()) shape = newshape;
01913    TGeoMedium *med = oldvol->GetMedium();
01914    if (newmed) med = newmed;
01915    // Make a new volume
01916    TGeoVolume *vol = new TGeoVolume(oldvol->GetName(), shape, med);
01917    // copy volume attributes
01918    vol->SetVisibility(oldvol->IsVisible());
01919    vol->SetLineColor(oldvol->GetLineColor());
01920    vol->SetLineStyle(oldvol->GetLineStyle());
01921    vol->SetLineWidth(oldvol->GetLineWidth());
01922    vol->SetFillColor(oldvol->GetFillColor());
01923    vol->SetFillStyle(oldvol->GetFillStyle());
01924    // copy field
01925    vol->SetField(oldvol->GetField());
01926    // Make a copy of the node
01927    TGeoNode *newnode = nodeorig->MakeCopyNode();
01928    // Change the volume for the new node
01929    newnode->SetVolume(vol);
01930    // Replace the matrix
01931    if (newpos && !nodeorig->IsOffset()) {
01932       TGeoNodeMatrix *nodemat = (TGeoNodeMatrix*)newnode;
01933       nodemat->SetMatrix(newpos);
01934    }   
01935    // Replace nodeorig with new one
01936    fNodes->RemoveAt(ind);
01937    fNodes->AddAt(newnode, ind);   
01938    if (fVoxels) fVoxels->SetNeedRebuild();
01939    if (IsAssembly()) fShape->ComputeBBox();
01940    return newnode;
01941 }      
01942 
01943 //_____________________________________________________________________________
01944 void TGeoVolume::SelectVolume(Bool_t clear)
01945 {
01946 // Select this volume as matching an arbitrary criteria. The volume is added to
01947 // a static list and the flag TGeoVolume::kVolumeSelected is set. All flags need
01948 // to be reset at the end by calling the method with CLEAR=true. This will also clear 
01949 // the list.
01950    static TObjArray array(256);
01951    static Int_t len = 0;
01952    Int_t i;
01953    TObject *vol;
01954    if (clear) {
01955       for (i=0; i<len; i++) {
01956          vol = array.At(i);
01957          vol->ResetBit(TGeoVolume::kVolumeSelected);
01958       }
01959       array.Clear();
01960       len = 0;
01961       return;
01962    }
01963    SetBit(TGeoVolume::kVolumeSelected);
01964    array.AddAtAndExpand(this, len++);
01965 }      
01966 
01967 //_____________________________________________________________________________
01968 void TGeoVolume::SetVisibility(Bool_t vis)
01969 {
01970 // set visibility of this volume
01971    if (IsAssembly()) {
01972       Info("SetVisibility", "Volume %s: assemblies do not have visibility", GetName());
01973       return;
01974    }   
01975    TGeoAtt::SetVisibility(vis);
01976    if (fGeoManager->IsClosed()) SetVisTouched(kTRUE);
01977    fGeoManager->SetVisOption(4);
01978    TSeqCollection *brlist = gROOT->GetListOfBrowsers();
01979    TIter next(brlist);
01980    TBrowser *browser = 0;
01981    while ((browser=(TBrowser*)next())) {
01982       browser->CheckObjectItem(this, vis);
01983       browser->Refresh();
01984    }
01985 }   
01986 
01987 //_____________________________________________________________________________
01988 void TGeoVolume::SetVisContainers(Bool_t flag)
01989 {
01990 // Set visibility for containers.
01991    TGeoAtt::SetVisContainers(flag);
01992    if (fGeoManager && fGeoManager->IsClosed()) {
01993       if (flag) fGeoManager->SetVisOption(TVirtualGeoPainter::kGeoVisDefault);
01994       else      fGeoManager->SetVisOption(TVirtualGeoPainter::kGeoVisLeaves);
01995    }   
01996 }
01997    
01998 //_____________________________________________________________________________
01999 void TGeoVolume::SetVisLeaves(Bool_t flag)
02000 {
02001 // Set visibility for leaves.
02002    TGeoAtt::SetVisLeaves(flag);
02003    if (fGeoManager && fGeoManager->IsClosed()) {
02004       if (flag) fGeoManager->SetVisOption(TVirtualGeoPainter::kGeoVisLeaves);
02005       else      fGeoManager->SetVisOption(TVirtualGeoPainter::kGeoVisDefault);
02006    }   
02007 }
02008 
02009 //_____________________________________________________________________________
02010 void TGeoVolume::SetVisOnly(Bool_t flag)
02011 {
02012 // Set visibility for leaves.
02013    if (IsAssembly()) return;
02014    TGeoAtt::SetVisOnly(flag);
02015    if (fGeoManager && fGeoManager->IsClosed()) {
02016       if (flag) fGeoManager->SetVisOption(TVirtualGeoPainter::kGeoVisOnly);
02017       else      fGeoManager->SetVisOption(TVirtualGeoPainter::kGeoVisLeaves);
02018    }   
02019 }
02020 
02021 //_____________________________________________________________________________
02022 Bool_t TGeoVolume::Valid() const
02023 {
02024 // Check if the shape of this volume is valid.
02025    return fShape->IsValidBox();
02026 }
02027 
02028 //_____________________________________________________________________________
02029 Bool_t TGeoVolume::FindMatrixOfDaughterVolume(TGeoVolume *vol) const
02030 {
02031 // Find a daughter node having VOL as volume and fill TGeoManager::fHMatrix
02032 // with its global matrix.
02033    if (vol == this) return kTRUE;
02034    Int_t nd = GetNdaughters();
02035    if (!nd) return kFALSE;
02036    TGeoHMatrix *global = fGeoManager->GetHMatrix();
02037    TGeoNode *dnode;
02038    TGeoVolume *dvol;
02039    TGeoMatrix *local;
02040    Int_t i;
02041    for (i=0; i<nd; i++) {
02042       dnode = GetNode(i);
02043       dvol = dnode->GetVolume();
02044       if (dvol == vol) {
02045          local = dnode->GetMatrix();
02046          global->MultiplyLeft(local);
02047          return kTRUE;
02048       }
02049    }
02050    for (i=0; i<nd; i++) {
02051       dnode = GetNode(i);
02052       dvol = dnode->GetVolume();
02053       if (dvol->FindMatrixOfDaughterVolume(vol)) return kTRUE;
02054    }
02055    return kFALSE;
02056 }                    
02057 
02058 //_____________________________________________________________________________
02059 void TGeoVolume::VisibleDaughters(Bool_t vis)
02060 {
02061 // set visibility for daughters
02062    SetVisDaughters(vis);
02063    if (fGeoManager->IsClosed()) SetVisTouched(kTRUE);
02064    fGeoManager->SetVisOption(4);
02065 }
02066 
02067 //_____________________________________________________________________________
02068 void TGeoVolume::Voxelize(Option_t *option)
02069 {
02070 // build the voxels for this volume 
02071    if (!Valid()) {
02072       Error("Voxelize", "Bounding box not valid");
02073       return; 
02074    }   
02075    // do not voxelize divided volumes
02076    if (fFinder) return;
02077    // or final leaves
02078    Int_t nd = GetNdaughters();
02079    if (!nd) return;
02080    // If this is an assembly, re-compute bounding box
02081    if (IsAssembly()) fShape->ComputeBBox();
02082    // delete old voxelization if any
02083    if (fVoxels) {
02084       if (!TObject::TestBit(kVolumeClone)) delete fVoxels;
02085       fVoxels = 0;
02086    }   
02087    // Create the voxels structure
02088    fVoxels = new TGeoVoxelFinder(this);
02089    fVoxels->Voxelize(option);
02090    if (fVoxels) {
02091       if (fVoxels->IsInvalid()) {
02092          delete fVoxels;
02093          fVoxels = 0;
02094       }
02095    }      
02096 }
02097 
02098 //_____________________________________________________________________________
02099 Double_t TGeoVolume::Weight(Double_t precision, Option_t *option)
02100 {
02101 // Estimate the weight of a volume (in kg) with SIGMA(M)/M better than PRECISION.
02102 // Option can contain : v - verbose, a - analytical  (default)
02103    TGeoVolume *top = fGeoManager->GetTopVolume();
02104    if (top != this) fGeoManager->SetTopVolume(this);
02105    else top = 0;
02106    Double_t weight =  fGeoManager->Weight(precision, option);
02107    if (top) fGeoManager->SetTopVolume(top);
02108    return weight;
02109 }   
02110 
02111 //_____________________________________________________________________________
02112 Double_t TGeoVolume::WeightA() const
02113 {
02114 // Analytical computation of the weight.
02115    Double_t capacity = Capacity();
02116    Double_t weight = 0.0;
02117    Int_t i;
02118    Int_t nd = GetNdaughters();
02119    TGeoVolume *daughter;
02120    for (i=0; i<nd; i++) {
02121       daughter = GetNode(i)->GetVolume();
02122       weight += daughter->WeightA();
02123       capacity -= daughter->Capacity();
02124    }
02125    Double_t density = 0.0;
02126    if (!IsAssembly()) {
02127       if (fMedium) density = fMedium->GetMaterial()->GetDensity();
02128       if (density<0.01) density = 0.0; // do not weight gases
02129    }   
02130    weight += 0.001*capacity * density; //[kg]
02131    return weight;
02132 }
02133 
02134 ClassImp(TGeoVolumeMulti)
02135 
02136 
02137 //_____________________________________________________________________________
02138 TGeoVolumeMulti::TGeoVolumeMulti()
02139 { 
02140 // dummy constructor
02141    fVolumes   = 0;
02142    fDivision = 0;
02143    fNumed = 0;
02144    fNdiv = 0;
02145    fAxis = 0;
02146    fStart = 0;
02147    fStep = 0;
02148    fAttSet = kFALSE;
02149    TObject::SetBit(kVolumeMulti);
02150 }
02151 
02152 //_____________________________________________________________________________
02153 TGeoVolumeMulti::TGeoVolumeMulti(const char *name, TGeoMedium *med)
02154 {
02155 // default constructor
02156    fVolumes = new TObjArray();
02157    fDivision = 0;
02158    fNumed = 0;
02159    fNdiv = 0;
02160    fAxis = 0;
02161    fStart = 0;
02162    fStep = 0;
02163    fAttSet = kFALSE;
02164    TObject::SetBit(kVolumeMulti);
02165    SetName(name);
02166    SetMedium(med);
02167    fGeoManager->AddVolume(this);
02168 //   printf("--- volume multi %s created\n", name);
02169 }
02170 
02171 //_____________________________________________________________________________
02172 TGeoVolumeMulti::TGeoVolumeMulti(const TGeoVolumeMulti& vm) :
02173   TGeoVolume(vm),
02174   fVolumes(vm.fVolumes),
02175   fDivision(vm.fDivision),
02176   fNumed(vm.fNumed),
02177   fNdiv(vm.fNdiv),
02178   fAxis(vm.fAxis),
02179   fStart(vm.fStart),
02180   fStep(vm.fStep),
02181   fAttSet(vm.fAttSet)
02182 { 
02183    //copy constructor
02184 }
02185 
02186 //_____________________________________________________________________________
02187 TGeoVolumeMulti& TGeoVolumeMulti::operator=(const TGeoVolumeMulti& vm) 
02188 {
02189    //assignment operator
02190    if(this!=&vm) {
02191       TGeoVolume::operator=(vm);
02192       fVolumes=vm.fVolumes;
02193       fDivision=vm.fDivision;
02194       fNumed=vm.fNumed;
02195       fNdiv=vm.fNdiv;
02196       fAxis=vm.fAxis;
02197       fStart=vm.fStart;
02198       fStep=vm.fStep;
02199       fAttSet=vm.fAttSet;
02200    } 
02201    return *this;
02202 }
02203 
02204 //_____________________________________________________________________________
02205 TGeoVolumeMulti::~TGeoVolumeMulti()
02206 {
02207 // Destructor
02208    if (fVolumes) delete fVolumes;
02209 }
02210 
02211 //_____________________________________________________________________________
02212 void TGeoVolumeMulti::AddVolume(TGeoVolume *vol) 
02213 {
02214 // Add a volume with valid shape to the list of volumes. Copy all existing nodes
02215 // to this volume
02216    Int_t idx = fVolumes->GetEntriesFast();
02217    fVolumes->AddAtAndExpand(vol,idx);
02218    vol->SetUniqueID(idx+1);
02219    TGeoVolumeMulti *div;
02220    TGeoVolume *cell;
02221    if (fDivision) {
02222       div = (TGeoVolumeMulti*)vol->Divide(fDivision->GetName(), fAxis, fNdiv, fStart, fStep, fNumed, fOption.Data());
02223       for (Int_t i=0; i<div->GetNvolumes(); i++) {
02224          cell = div->GetVolume(i);
02225          fDivision->AddVolume(cell);
02226       }
02227    }      
02228    if (fNodes) {
02229       Int_t nd = fNodes->GetEntriesFast();
02230       for (Int_t id=0; id<nd; id++) {
02231          TGeoNode *node = (TGeoNode*)fNodes->At(id);
02232          Bool_t many = node->IsOverlapping();
02233          if (many) vol->AddNodeOverlap(node->GetVolume(), node->GetNumber(), node->GetMatrix());
02234          else      vol->AddNode(node->GetVolume(), node->GetNumber(), node->GetMatrix());
02235       }
02236    }      
02237 //      vol->MakeCopyNodes(this);
02238 }
02239    
02240 
02241 //_____________________________________________________________________________
02242 void TGeoVolumeMulti::AddNode(const TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option)
02243 {
02244 // Add a new node to the list of nodes. This is the usual method for adding
02245 // daughters inside the container volume.
02246    TGeoVolume::AddNode(vol, copy_no, mat, option);
02247    Int_t nvolumes = fVolumes->GetEntriesFast();
02248    TGeoVolume *volume = 0;
02249    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02250       volume = GetVolume(ivo);
02251       volume->SetLineColor(GetLineColor());
02252       volume->SetLineStyle(GetLineStyle());
02253       volume->SetLineWidth(GetLineWidth());
02254       volume->SetVisibility(IsVisible());
02255       volume->AddNode(vol, copy_no, mat, option); 
02256    }
02257 //   printf("--- vmulti %s : node %s added to %i components\n", GetName(), vol->GetName(), nvolumes);
02258 }
02259 
02260 //_____________________________________________________________________________
02261 void TGeoVolumeMulti::AddNodeOverlap(const TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option)
02262 {
02263 // Add a new node to the list of nodes, This node is possibly overlapping with other
02264 // daughters of the volume or extruding the volume.
02265    TGeoVolume::AddNodeOverlap(vol, copy_no, mat, option);
02266    Int_t nvolumes = fVolumes->GetEntriesFast();
02267    TGeoVolume *volume = 0;
02268    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02269       volume = GetVolume(ivo);
02270       volume->SetLineColor(GetLineColor());
02271       volume->SetLineStyle(GetLineStyle());
02272       volume->SetLineWidth(GetLineWidth());
02273       volume->SetVisibility(IsVisible());
02274       volume->AddNodeOverlap(vol, copy_no, mat, option); 
02275    }
02276 //   printf("--- vmulti %s : node ovlp %s added to %i components\n", GetName(), vol->GetName(), nvolumes);
02277 }
02278 
02279 
02280 //_____________________________________________________________________________
02281 TGeoVolume *TGeoVolumeMulti::Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed, const char *option)
02282 {
02283 // division of multiple volumes
02284    if (fDivision) {
02285       Error("Divide", "volume %s already divided", GetName());
02286       return 0;
02287    }   
02288    Int_t nvolumes = fVolumes->GetEntriesFast();
02289    TGeoMedium *medium = fMedium;
02290    if (numed) {
02291       medium = fGeoManager->GetMedium(numed);
02292       if (!medium) {
02293          Error("Divide", "Invalid medium number %d for division volume %s", numed, divname);
02294          medium = fMedium;
02295       }
02296    }      
02297    if (!nvolumes) {
02298       // this is a virtual volume
02299       fDivision = new TGeoVolumeMulti(divname, medium);
02300       fNumed = medium->GetId();
02301       fOption = option;
02302       fAxis = iaxis;
02303       fNdiv = ndiv;
02304       fStart = start;
02305       fStep = step;
02306       // nothing else to do at this stage
02307       return fDivision;
02308    }      
02309    TGeoVolume *vol = 0;
02310    fDivision = new TGeoVolumeMulti(divname, medium);
02311    if (medium) fNumed = medium->GetId();
02312    fOption = option;
02313    fAxis = iaxis;
02314    fNdiv = ndiv;
02315    fStart = start;
02316    fStep = step;
02317    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02318       vol = GetVolume(ivo);
02319       vol->SetLineColor(GetLineColor());
02320       vol->SetLineStyle(GetLineStyle());
02321       vol->SetLineWidth(GetLineWidth());
02322       vol->SetVisibility(IsVisible());
02323       fDivision->AddVolume(vol->Divide(divname,iaxis,ndiv,start,step, numed, option)); 
02324    }
02325 //   printf("--- volume multi %s (%i volumes) divided\n", GetName(), nvolumes);
02326    if (numed) fDivision->SetMedium(medium);
02327    return fDivision;
02328 }
02329 
02330 //_____________________________________________________________________________
02331 TGeoVolume *TGeoVolumeMulti::MakeCopyVolume(TGeoShape *newshape)
02332 {
02333    // Make a copy of this volume
02334    // build a volume with same name, shape and medium
02335    TGeoVolume *vol = new TGeoVolume(GetName(), newshape, fMedium);
02336    Int_t i=0;
02337    // copy volume attributes
02338    vol->SetVisibility(IsVisible());
02339    vol->SetLineColor(GetLineColor());
02340    vol->SetLineStyle(GetLineStyle());
02341    vol->SetLineWidth(GetLineWidth());
02342    vol->SetFillColor(GetFillColor());
02343    vol->SetFillStyle(GetFillStyle());
02344    // copy field
02345    vol->SetField(fField);
02346    // if divided, copy division object
02347 //    if (fFinder) {
02348 //       Error("MakeCopyVolume", "volume %s divided", GetName());
02349 //       vol->SetFinder(fFinder);
02350 //    }
02351    if (fDivision) {
02352       TGeoVolume *cell;
02353       TGeoVolumeMulti *div = (TGeoVolumeMulti*)vol->Divide(fDivision->GetName(), fAxis, fNdiv, fStart, fStep, fNumed, fOption.Data());
02354       for (i=0; i<div->GetNvolumes(); i++) {
02355          cell = div->GetVolume(i);
02356          fDivision->AddVolume(cell);
02357       }
02358    }      
02359                  
02360    if (!fNodes) return vol;
02361    TGeoNode *node;
02362    Int_t nd = fNodes->GetEntriesFast();
02363    if (!nd) return vol;
02364    // create new list of nodes
02365    TObjArray *list = new TObjArray();
02366    // attach it to new volume
02367    vol->SetNodes(list);
02368    ((TObject*)vol)->SetBit(kVolumeImportNodes);
02369    for (i=0; i<nd; i++) {
02370       //create copies of nodes and add them to list
02371       node = GetNode(i)->MakeCopyNode();
02372       node->SetMotherVolume(vol);
02373       list->Add(node);
02374    }
02375    return vol;       
02376 }    
02377 
02378 //_____________________________________________________________________________
02379 void TGeoVolumeMulti::SetLineColor(Color_t lcolor) 
02380 {
02381 // Set the line color for all components.
02382    TGeoVolume::SetLineColor(lcolor);
02383    Int_t nvolumes = fVolumes->GetEntriesFast();
02384    TGeoVolume *vol = 0;
02385    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02386       vol = GetVolume(ivo);
02387       vol->SetLineColor(lcolor); 
02388    }
02389 }
02390 
02391 //_____________________________________________________________________________
02392 void TGeoVolumeMulti::SetLineStyle(Style_t lstyle) 
02393 {
02394 // Set the line style for all components.
02395    TGeoVolume::SetLineStyle(lstyle); 
02396    Int_t nvolumes = fVolumes->GetEntriesFast();
02397    TGeoVolume *vol = 0;
02398    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02399       vol = GetVolume(ivo);
02400       vol->SetLineStyle(lstyle); 
02401    }
02402 }
02403 
02404 //_____________________________________________________________________________
02405 void TGeoVolumeMulti::SetLineWidth(Width_t lwidth) 
02406 {
02407 // Set the line width for all components.
02408    TGeoVolume::SetLineWidth(lwidth);
02409    Int_t nvolumes = fVolumes->GetEntriesFast();
02410    TGeoVolume *vol = 0;
02411    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02412       vol = GetVolume(ivo);
02413       vol->SetLineWidth(lwidth); 
02414    }
02415 }
02416 
02417 //_____________________________________________________________________________
02418 void TGeoVolumeMulti::SetMedium(TGeoMedium *med)
02419 {
02420 // Set medium for a multiple volume.
02421    TGeoVolume::SetMedium(med);
02422    Int_t nvolumes = fVolumes->GetEntriesFast();
02423    TGeoVolume *vol = 0;
02424    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02425       vol = GetVolume(ivo);
02426       vol->SetMedium(med); 
02427    }
02428 }   
02429 
02430 
02431 //_____________________________________________________________________________
02432 void TGeoVolumeMulti::SetVisibility(Bool_t vis) 
02433 {
02434 // Set visibility for all components.
02435    TGeoVolume::SetVisibility(vis); 
02436    Int_t nvolumes = fVolumes->GetEntriesFast();
02437    TGeoVolume *vol = 0;
02438    for (Int_t ivo=0; ivo<nvolumes; ivo++) {
02439       vol = GetVolume(ivo);
02440       vol->SetVisibility(vis); 
02441    }
02442 }
02443 
02444 ClassImp(TGeoVolumeAssembly)
02445 
02446 //_____________________________________________________________________________
02447 TGeoVolumeAssembly::TGeoVolumeAssembly()
02448                    :TGeoVolume()
02449 {
02450 // Default constructor
02451    fCurrent = -1;
02452    fNext = -1;
02453 }
02454 
02455 //_____________________________________________________________________________
02456 TGeoVolumeAssembly::TGeoVolumeAssembly(const char *name)
02457                    :TGeoVolume()
02458 {
02459 // Constructor. Just the name has to be provided. Assemblies does not have their own
02460 // shape or medium.
02461    fName = name;
02462    fName = fName.Strip();
02463    fCurrent = -1;
02464    fNext = -1;
02465    fShape = new TGeoShapeAssembly(this);
02466    if (fGeoManager) fNumber = fGeoManager->AddVolume(this);
02467 }
02468 
02469 //_____________________________________________________________________________
02470 TGeoVolumeAssembly::~TGeoVolumeAssembly()
02471 {
02472 // Destructor. The assembly is owner of its "shape".
02473    if (fShape) delete fShape;
02474 }   
02475 
02476 //_____________________________________________________________________________
02477 void TGeoVolumeAssembly::AddNode(const TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option)
02478 {
02479 // Add a component to the assembly. 
02480    TGeoVolume::AddNode(vol,copy_no,mat,option);
02481    ((TGeoShapeAssembly*)fShape)->RecomputeBoxLast();
02482 }   
02483 
02484 //_____________________________________________________________________________
02485 void TGeoVolumeAssembly::AddNodeOverlap(const TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat, Option_t *option)
02486 {
02487 // Add an overlapping node - not allowed for assemblies.
02488    Warning("AddNodeOverlap", "Declaring assembly %s as possibly overlapping inside %s not allowed. Using AddNode instead !",vol->GetName(),GetName());
02489    AddNode(vol, copy_no, mat, option);
02490 }   
02491 
02492 //_____________________________________________________________________________
02493 TGeoVolume *TGeoVolumeAssembly::CloneVolume() const
02494 {
02495 // Clone this volume.
02496    // build a volume with same name, shape and medium
02497    TGeoVolume *vol = new TGeoVolumeAssembly(GetName());
02498    Int_t i;
02499    // copy other attributes
02500    Int_t nbits = 8*sizeof(UInt_t);
02501    for (i=0; i<nbits; i++) 
02502       vol->SetAttBit(1<<i, TGeoAtt::TestAttBit(1<<i));
02503    for (i=14; i<24; i++)
02504       vol->SetBit(1<<i, TestBit(1<<i));   
02505    
02506    // copy field
02507    vol->SetField(fField);
02508    // Set bits
02509    for (i=0; i<nbits; i++) 
02510       vol->SetBit(1<<i, TObject::TestBit(1<<i));
02511    vol->SetBit(kVolumeClone);   
02512    // make copy nodes
02513    vol->MakeCopyNodes(this);
02514 //   CloneNodesAndConnect(vol);
02515    ((TGeoShapeAssembly*)vol->GetShape())->NeedsBBoxRecompute();
02516    // copy voxels
02517    TGeoVoxelFinder *voxels = 0;
02518    if (fVoxels) {
02519       voxels = new TGeoVoxelFinder(vol);
02520       vol->SetVoxelFinder(voxels);
02521    }   
02522    // copy option, uid
02523    vol->SetOption(fOption);
02524    vol->SetNumber(fNumber);
02525    vol->SetNtotal(fNtotal);
02526    return vol;
02527 }
02528 
02529 //_____________________________________________________________________________
02530 TGeoVolume *TGeoVolumeAssembly::Divide(const char *, Int_t, Int_t, Double_t, Double_t, Int_t, Option_t *)
02531 {
02532 // Division makes no sense for assemblies.
02533    Error("Divide","Assemblies cannot be divided");
02534    return 0;
02535 }
02536 
02537 //_____________________________________________________________________________
02538 TGeoVolumeAssembly *TGeoVolumeAssembly::MakeAssemblyFromVolume(TGeoVolume *volorig)
02539 {
02540 // Make a clone of volume VOL but which is an assembly.
02541    if (volorig->IsAssembly() || volorig->IsVolumeMulti()) return 0;
02542    Int_t nd = volorig->GetNdaughters();
02543    if (!nd) return 0;
02544    TGeoVolumeAssembly *vol = new TGeoVolumeAssembly(volorig->GetName());
02545    Int_t i;
02546    // copy other attributes
02547    Int_t nbits = 8*sizeof(UInt_t);
02548    for (i=0; i<nbits; i++) 
02549       vol->SetAttBit(1<<i, volorig->TestAttBit(1<<i));
02550    for (i=14; i<24; i++)
02551       vol->SetBit(1<<i, volorig->TestBit(1<<i));   
02552    
02553    // copy field
02554    vol->SetField(volorig->GetField());
02555    // Set bits
02556    for (i=0; i<nbits; i++) 
02557       vol->SetBit(1<<i, volorig->TestBit(1<<i));
02558    vol->SetBit(kVolumeClone);   
02559    // make copy nodes
02560    vol->MakeCopyNodes(volorig);
02561 //   volorig->CloneNodesAndConnect(vol);
02562    vol->GetShape()->ComputeBBox();
02563    // copy voxels
02564    TGeoVoxelFinder *voxels = 0;
02565    if (volorig->GetVoxels()) {
02566       voxels = new TGeoVoxelFinder(vol);
02567       vol->SetVoxelFinder(voxels);
02568    }   
02569    // copy option, uid
02570    vol->SetOption(volorig->GetOption());
02571    vol->SetNumber(volorig->GetNumber());
02572    vol->SetNtotal(volorig->GetNtotal());
02573    return vol;
02574 }   

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