ROOT logo
// @(#):$Id: TGeoBoolNode.cxx 27731 2009-03-09 17:40:56Z brun $
// Author: Andrei Gheata   30/05/02
// TGeoBoolNode::Contains and parser implemented by Mihaela Gheata

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

#include "Riostream.h"

#include "TGeoCompositeShape.h"
#include "TGeoMatrix.h"
#include "TGeoManager.h"

#include "TGeoBoolNode.h"

#include "TVirtualPad.h"
#include "TVirtualViewer3D.h"
#include "TBuffer3D.h"
#include "TBuffer3DTypes.h"
#include "TMath.h"

//_____________________________________________________________________________
//  TGeoBoolNode - base class for Boolean operations between two shapes.
//===============
// A Boolean node describes a Boolean operation between 'left' and 'right' 
// shapes positioned with respect to an ARBITRARY reference frame. The boolean
// node is referenced by a mother composite shape and its shape components may
// be primitive but also composite shapes. The later situation leads to a binary
// tree hierarchy. When the parent composite shape is used to create a volume,
// the reference frame of the volume is chosen to match the frame in which 
// node shape components were defined.
//
// The positioned shape components may or may not be disjoint. The specific 
// implementations for Boolean nodes are:
//
//    TGeoUnion - representing the Boolean  union of two positioned shapes
//
//    TGeoSubtraction - representing the Boolean subtraction of two positioned 
//                shapes
// 
//    TGeoIntersection - representing the Boolean intersection of two positioned
//                shapes
//_____________________________________________________________________________

ClassImp(TGeoBoolNode)

//-----------------------------------------------------------------------------
TGeoBoolNode::TGeoBoolNode()
{
// Default constructor
   fLeft     = 0;
   fRight    = 0;
   fLeftMat  = 0;
   fRightMat = 0;
   fSelected = 0;
}
//-----------------------------------------------------------------------------
TGeoBoolNode::TGeoBoolNode(const char *expr1, const char *expr2)
{
// Constructor called by TGeoCompositeShape providing 2 subexpressions for the 2 branches.
   fLeft     = 0;
   fRight    = 0;
   fLeftMat  = 0;
   fRightMat = 0;
   fSelected = 0;
   if (!MakeBranch(expr1, kTRUE)) {
      return;
   }
   if (!MakeBranch(expr2, kFALSE)) {
      return;
   }
}

//-----------------------------------------------------------------------------
TGeoBoolNode::TGeoBoolNode(TGeoShape *left, TGeoShape *right, TGeoMatrix *lmat, TGeoMatrix *rmat)
{
// Constructor providing left and right shapes and matrices (in the Boolean operation).
   fSelected = 0;
   fLeft = left;
   fRight = right;
   fLeftMat = lmat;
   if (!fLeftMat) fLeftMat = gGeoIdentity;
   else fLeftMat->RegisterYourself();
   fRightMat = rmat;
   if (!fRightMat) fRightMat = gGeoIdentity;
   else fRightMat->RegisterYourself();
   if (!fLeft) {
      Error("ctor", "left shape is NULL");
      return;
   }   
   if (!fRight) {
      Error("ctor", "right shape is NULL");
      return;
   }   
}

//-----------------------------------------------------------------------------
TGeoBoolNode::~TGeoBoolNode()
{
// Destructor.
// --- deletion of components handled by TGeoManager class.
}
//-----------------------------------------------------------------------------
Bool_t TGeoBoolNode::MakeBranch(const char *expr, Bool_t left)
{
// Expands the boolean expression either on left or right branch, creating
// component elements (composite shapes and boolean nodes). Returns true on success.
   TString sleft, sright, stransf;
   Int_t boolop = TGeoManager::Parse(expr, sleft, sright, stransf);
   if (boolop<0) {
      Error("MakeBranch", "invalid expresion");
      return kFALSE;
   }
   TGeoShape *shape = 0;
   TGeoMatrix *mat;
   TString newshape;

   if (stransf.Length() == 0) {
      mat = gGeoIdentity;
   } else {   
      mat = (TGeoMatrix*)gGeoManager->GetListOfMatrices()->FindObject(stransf.Data());    
   }
   if (!mat) {
      Error("MakeBranch", "transformation %s not found", stransf.Data());
      return kFALSE;
   }
   switch (boolop) {
      case 0:
         // elementary shape
         shape = (TGeoShape*)gGeoManager->GetListOfShapes()->FindObject(sleft.Data()); 
         if (!shape) {
            Error("MakeBranch", "shape %s not found", sleft.Data());
            return kFALSE;
         }
         break;
      case 1:
         // composite shape - union
         newshape = sleft;
         newshape += "+";
         newshape += sright;
         shape = new TGeoCompositeShape(newshape.Data());
         break;
      case 2:
         // composite shape - difference
         newshape = sleft;
         newshape += "-";
         newshape += sright;
         shape = new TGeoCompositeShape(newshape.Data());
         break;
      case 3:
         // composite shape - intersection
         newshape = sleft;
         newshape += "*";
         newshape += sright;
         shape = new TGeoCompositeShape(newshape.Data());
         break;
   }      
   if (boolop && !shape->IsValid()) {
      Error("MakeBranch", "Shape %s not valid", newshape.Data());
      delete shape;
      return kFALSE;
   }      
   if (left) {
      fLeft = shape;
      fLeftMat = mat;
   } else {
      fRight = shape;
      fRightMat = mat;
   }
   return kTRUE;                  
}
//-----------------------------------------------------------------------------
void TGeoBoolNode::Paint(Option_t * option)
{
// Special schema for feeding the 3D buffers to the painter client.
   TVirtualViewer3D * viewer = gPad->GetViewer3D();
   if (!viewer) return;

   // Components of composite shape hierarchies for local frame viewers are painted 
   // in coordinate frame of the top level composite shape. So we force 
   // conversion to this.  See TGeoPainter::PaintNode for loading of GLMatrix.
   Bool_t localFrame = kFALSE; //viewer->PreferLocalFrame();

   TGeoHMatrix *glmat = (TGeoHMatrix*)TGeoShape::GetTransform();
   TGeoHMatrix mat;
   mat = glmat; // keep a copy

   // Now perform fetch and add of the two components buffers.
   // Note we assume that composite shapes are always completely added
   // so don't bother to get addDaughters flag from viewer->AddObject()

   // Setup matrix and fetch/add the left component buffer
   glmat->Multiply(fLeftMat);
   //fLeft->Paint(option);
   if (TGeoCompositeShape *left = dynamic_cast<TGeoCompositeShape *>(fLeft)) left->PaintComposite(option);
   else {
      const TBuffer3D & leftBuffer = fLeft->GetBuffer3D(TBuffer3D::kAll, localFrame);
      viewer->AddObject(leftBuffer);
   }

   // Setup matrix and fetch/add the right component buffer
   *glmat = &mat;
   glmat->Multiply(fRightMat);
   //fRight->Paint(option);
   if (TGeoCompositeShape *right = dynamic_cast<TGeoCompositeShape *>(fRight)) right->PaintComposite(option);
   else {
      const TBuffer3D & rightBuffer = fRight->GetBuffer3D(TBuffer3D::kAll, localFrame);
      viewer->AddObject(rightBuffer);
   }

   *glmat = &mat;   
}
//_____________________________________________________________________________
void TGeoBoolNode::RegisterMatrices()
{
// Register all matrices of the boolean node and descendents.
   if (!fLeftMat->IsIdentity()) fLeftMat->RegisterYourself();   
   if (!fRightMat->IsIdentity()) fRightMat->RegisterYourself();   
   if (fLeft->IsComposite()) ((TGeoCompositeShape*)fLeft)->GetBoolNode()->RegisterMatrices();
   if (fRight->IsComposite()) ((TGeoCompositeShape*)fRight)->GetBoolNode()->RegisterMatrices();
}
//_____________________________________________________________________________
void TGeoBoolNode::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
// Save a primitive as a C++ statement(s) on output stream "out".
   fLeft->SavePrimitive(out,option);
   fRight->SavePrimitive(out,option);
   if (!fLeftMat->IsIdentity()) {
      fLeftMat->RegisterYourself();
      fLeftMat->SavePrimitive(out,option);
   }      
   if (!fRightMat->IsIdentity()) {
      fRightMat->RegisterYourself();
      fRightMat->SavePrimitive(out,option);
   }      
}
//-----------------------------------------------------------------------------
void TGeoBoolNode::Sizeof3D() const
{
// Register size of this 3D object
   fLeft->Sizeof3D();
   fRight->Sizeof3D();
}
ClassImp(TGeoUnion)

//-----------------------------------------------------------------------------
void TGeoUnion::Paint(Option_t *option)
{
// Paint method.
   TVirtualViewer3D *viewer = gPad->GetViewer3D();

   if (!viewer) {
      Error("Paint", "gPad->GetViewer3D() returned 0, cannot work with composite!\n");
      return;
   }

   viewer->AddCompositeOp(TBuffer3D::kCSUnion);

   TGeoBoolNode::Paint(option);
}

//-----------------------------------------------------------------------------
TGeoUnion::TGeoUnion()
{
// Default constructor
}
//-----------------------------------------------------------------------------
TGeoUnion::TGeoUnion(const char *expr1, const char *expr2)
          :TGeoBoolNode(expr1, expr2)
{
// Constructor
}

//-----------------------------------------------------------------------------
TGeoUnion::TGeoUnion(TGeoShape *left, TGeoShape *right, TGeoMatrix *lmat, TGeoMatrix *rmat)
          :TGeoBoolNode(left,right,lmat,rmat)
{
// Constructor providing pointers to components
   if (left->TestShapeBit(TGeoShape::kGeoHalfSpace) || right->TestShapeBit(TGeoShape::kGeoHalfSpace)) {
      Fatal("TGeoUnion", "Unions with a half-space (%s + %s) not allowed", left->GetName(), right->GetName());
   }
}

//-----------------------------------------------------------------------------
TGeoUnion::~TGeoUnion()
{
// Destructor
// --- deletion of components handled by TGeoManager class.
}
//-----------------------------------------------------------------------------
void TGeoUnion::ComputeBBox(Double_t &dx, Double_t &dy, Double_t &dz, Double_t *origin)
{
// Compute bounding box corresponding to a union of two shapes.
   if (((TGeoBBox*)fLeft)->IsNullBox()) fLeft->ComputeBBox();
   if (((TGeoBBox*)fRight)->IsNullBox()) fRight->ComputeBBox();
   Double_t vert[48];
   Double_t pt[3];
   Int_t i;
   Double_t xmin, xmax, ymin, ymax, zmin, zmax;
   xmin = ymin = zmin = TGeoShape::Big();
   xmax = ymax = zmax = -TGeoShape::Big();
   ((TGeoBBox*)fLeft)->SetBoxPoints(&vert[0]);
   ((TGeoBBox*)fRight)->SetBoxPoints(&vert[24]);
   for (i=0; i<8; i++) {
      fLeftMat->LocalToMaster(&vert[3*i], &pt[0]);
      if (pt[0]<xmin) xmin=pt[0];
      if (pt[0]>xmax) xmax=pt[0];
      if (pt[1]<ymin) ymin=pt[1];
      if (pt[1]>ymax) ymax=pt[1];
      if (pt[2]<zmin) zmin=pt[2];
      if (pt[2]>zmax) zmax=pt[2];
   }   
   for (i=8; i<16; i++) {
      fRightMat->LocalToMaster(&vert[3*i], &pt[0]);
      if (pt[0]<xmin) xmin=pt[0];
      if (pt[0]>xmax) xmax=pt[0];
      if (pt[1]<ymin) ymin=pt[1];
      if (pt[1]>ymax) ymax=pt[1];
      if (pt[2]<zmin) zmin=pt[2];
      if (pt[2]>zmax) zmax=pt[2];
   }   
   dx = 0.5*(xmax-xmin);
   origin[0] = 0.5*(xmin+xmax);
   dy = 0.5*(ymax-ymin);
   origin[1] = 0.5*(ymin+ymax);
   dz = 0.5*(zmax-zmin);
   origin[2] = 0.5*(zmin+zmax);
}   
//-----------------------------------------------------------------------------
Bool_t TGeoUnion::Contains(Double_t *point) const
{
// Find if a union of two shapes contains a given point
   Double_t local[3];
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   fLeftMat->MasterToLocal(point, &local[0]);
   Bool_t inside = fLeft->Contains(&local[0]);
   if (inside) {
      node->SetSelected(1);
      return kTRUE;
   }   
   fRightMat->MasterToLocal(point, &local[0]);
   inside = fRight->Contains(&local[0]);
   if (inside) node->SetSelected(2);
   return inside;
}

//_____________________________________________________________________________
void TGeoUnion::ComputeNormal(Double_t *point, Double_t *dir, Double_t *norm)
{
// Normal computation in POINT. The orientation is chosen so that DIR.dot.NORM>0.
   norm[0] = norm[1] = 0.;
   norm[2] = 1.;
   Double_t local[3];
   Double_t ldir[3], lnorm[3];
   if (fSelected == 1) {
      fLeftMat->MasterToLocal(point, local);
      fLeftMat->MasterToLocalVect(dir, ldir);
      fLeft->ComputeNormal(local,ldir,lnorm);
      fLeftMat->LocalToMasterVect(lnorm, norm);
      return;
   }
   if (fSelected == 2) {
      fRightMat->MasterToLocal(point, local);
      fRightMat->MasterToLocalVect(dir, ldir);
      fRight->ComputeNormal(local,ldir,lnorm);
      fRightMat->LocalToMasterVect(lnorm, norm);
      return;
   }            
   fLeftMat->MasterToLocal(point, local);
   if (fLeft->Contains(local)) {
      fLeftMat->MasterToLocalVect(dir, ldir);
      fLeft->ComputeNormal(local,ldir,lnorm);
      fLeftMat->LocalToMasterVect(lnorm, norm);
      return;
   }   
   fRightMat->MasterToLocal(point, local);
   if (fRight->Contains(local)) {
      fRightMat->MasterToLocalVect(dir, ldir);
      fRight->ComputeNormal(local,ldir,lnorm);
      fRightMat->LocalToMasterVect(lnorm, norm);
      return;
   }   
   // Propagate forward/backward to see which of the components is intersected first
   local[0] = point[0] + 1E-5*dir[0];
   local[1] = point[1] + 1E-5*dir[1];
   local[2] = point[2] + 1E-5*dir[2];

   if (!Contains(local)) {
      local[0] = point[0] - 1E-5*dir[0];
      local[1] = point[1] - 1E-5*dir[1];
      local[2] = point[2] - 1E-5*dir[2];
      if (!Contains(local)) return;
   }
   ComputeNormal(local,dir,norm);   
}

//-----------------------------------------------------------------------------
Int_t TGeoUnion::DistanceToPrimitive(Int_t /*px*/, Int_t /*py*/)
{
// Compute minimum distance to shape vertices.
   return 9999;
}
//-----------------------------------------------------------------------------
Double_t TGeoUnion::DistFromInside(Double_t *point, Double_t *dir, Int_t iact,
                              Double_t step, Double_t *safe) const
{
// Computes distance from a given point inside the shape to its boundary.
   if (iact<3 && safe) {
      // compute safe distance
      *safe = Safety(point,kTRUE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }

   Double_t local[3], master[3], ldir[3], rdir[3], pushed[3];
   memcpy(master, point, 3*sizeof(Double_t));
   Int_t i;
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   Double_t d1=0., d2=0., snxt=0.;
   fLeftMat->MasterToLocalVect(dir, ldir);
   fRightMat->MasterToLocalVect(dir, rdir);
   fLeftMat->MasterToLocal(point, local);
   Bool_t inside1 = fLeft->Contains(local);
   if (inside1) d1 = fLeft->DistFromInside(local, ldir, 3);
   fRightMat->MasterToLocal(point, local);
   Bool_t inside2 = fRight->Contains(local);
   if (inside2) d2 = fRight->DistFromInside(local, rdir, 3);

   while (inside1 || inside2) {
      if (inside1 && inside2) {
         if (d1<d2) {      
            snxt += d1;
            node->SetSelected(1);
            // propagate to exit of left shape
            inside1 = kFALSE;
            for (i=0; i<3; i++) master[i] += d1*dir[i];
            // check if propagated point is in right shape        
            fRightMat->MasterToLocal(master, local);
            inside2 = fRight->Contains(local);
            if (!inside2) return snxt;
            d2 = fRight->DistFromInside(local, rdir, 3);
            if (d2 < TGeoShape::Tolerance()) return snxt;
         } else {
            snxt += d2;
            node->SetSelected(2);
            // propagate to exit of right shape
            inside2 = kFALSE;
            for (i=0; i<3; i++) master[i] += d2*dir[i];
            // check if propagated point is in left shape        
            fLeftMat->MasterToLocal(master, local);
            inside1 = fLeft->Contains(local);
            if (!inside1) return snxt;
            d1 = fLeft->DistFromInside(local, ldir, 3);
            if (d1 < TGeoShape::Tolerance()) return snxt;
         }
      } 
      if (inside1) {
         snxt += d1;
         node->SetSelected(1);
         // propagate to exit of left shape
         inside1 = kFALSE;
         for (i=0; i<3; i++) {
            master[i] += d1*dir[i];
            pushed[i] = master[i]+(1.+d1)*TGeoShape::Tolerance()*dir[i];
         }   
         // check if propagated point is in right shape        
         fRightMat->MasterToLocal(pushed, local);
         inside2 = fRight->Contains(local);
         if (!inside2) return snxt;
         d2 = fRight->DistFromInside(local, rdir, 3);
         if (d2 < TGeoShape::Tolerance()) return snxt;
         d2 += (1.+d1)*TGeoShape::Tolerance();
      }   
      if (inside2) {
         snxt += d2;
         node->SetSelected(2);
         // propagate to exit of right shape
         inside2 = kFALSE;
         for (i=0; i<3; i++) {
            master[i] += d2*dir[i];
            pushed[i] = master[i]+(1.+d2)*TGeoShape::Tolerance()*dir[i];
         }   
         // check if propagated point is in left shape        
         fLeftMat->MasterToLocal(pushed, local);
         inside1 = fLeft->Contains(local);
         if (!inside1) return snxt;
         d1 = fLeft->DistFromInside(local, ldir, 3);
         if (d1 < TGeoShape::Tolerance()) return snxt;
         d1 += (1.+d2)*TGeoShape::Tolerance();
      }
   }      
   return snxt;
}
//-----------------------------------------------------------------------------
Double_t TGeoUnion::DistFromOutside(Double_t *point, Double_t *dir, Int_t iact,
                              Double_t step, Double_t *safe) const
{
// Compute distance from a given outside point to the shape.
   if (iact<3 && safe) {
      // compute safe distance
      *safe = Safety(point,kFALSE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   Double_t local[3], ldir[3], rdir[3];
   Double_t d1, d2, snxt;
   fLeftMat->MasterToLocal(point, &local[0]);
   fLeftMat->MasterToLocalVect(dir, &ldir[0]);
   fRightMat->MasterToLocalVect(dir, &rdir[0]);
   d1 = fLeft->DistFromOutside(&local[0], &ldir[0], iact, step, safe);
   fRightMat->MasterToLocal(point, &local[0]);
   d2 = fRight->DistFromOutside(&local[0], &rdir[0], iact, step, safe);
   if (d1<d2) {
      snxt = d1;
      node->SetSelected(1);
   } else {
      snxt = d2;
      node->SetSelected(2);
   }      
   return snxt;
}
//-----------------------------------------------------------------------------
Int_t TGeoUnion::GetNpoints() const
{
// Returns number of vertices for the composite shape described by this union.
   return 0;
}
//-----------------------------------------------------------------------------
void TGeoUnion::SetPoints(Double_t * /*points*/) const
{
// Fill buffer with shape vertices.
}

//-----------------------------------------------------------------------------
Double_t TGeoUnion::Safety(Double_t *point, Bool_t in) const
{
// Compute safety distance for a union node;
   Double_t local1[3], local2[3];
   fLeftMat->MasterToLocal(point,local1);
   Bool_t in1 = fLeft->Contains(local1);
   fRightMat->MasterToLocal(point,local2);
   Bool_t in2 = fRight->Contains(local2);
   Bool_t intrue = in1 | in2;
   if (intrue^in) return 0.0;
   Double_t saf1 = fLeft->Safety(local1, in1);
   Double_t saf2 = fRight->Safety(local2, in2);
   if (in1 && in2) return TMath::Min(saf1, saf2);
   if (in1)        return saf1;
   if (in2)        return saf2;
   return TMath::Min(saf1,saf2);
}   

//_____________________________________________________________________________
void TGeoUnion::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
// Save a primitive as a C++ statement(s) on output stream "out".
   TGeoBoolNode::SavePrimitive(out,option);
   out << "   pBoolNode = new TGeoUnion(";
   out << fLeft->GetPointerName() << ",";
   out << fRight->GetPointerName() << ",";
   if (!fLeftMat->IsIdentity()) out << fLeftMat->GetPointerName() << ",";
   else                         out << "0,";
   if (!fRightMat->IsIdentity()) out << fRightMat->GetPointerName() << ");" << endl;
   else                         out << "0);" << endl;
}   

//-----------------------------------------------------------------------------
void TGeoUnion::SetPoints(Float_t * /*points*/) const
{
// Fill buffer with shape vertices.
}
//-----------------------------------------------------------------------------
void TGeoUnion::Sizeof3D() const
{
// Register 3D size of this shape.
   TGeoBoolNode::Sizeof3D();
}
   

ClassImp(TGeoSubtraction)

//-----------------------------------------------------------------------------
void TGeoSubtraction::Paint(Option_t *option)
{
// Paint method.
   TVirtualViewer3D *viewer = gPad->GetViewer3D();

   if (!viewer) {
      Error("Paint", "gPad->GetViewer3D() returned 0, cannot work with composite!\n");
      return;
   }

   viewer->AddCompositeOp(TBuffer3D::kCSDifference);

   TGeoBoolNode::Paint(option);
}

//-----------------------------------------------------------------------------
TGeoSubtraction::TGeoSubtraction()
{
// Default constructor
}

//-----------------------------------------------------------------------------
TGeoSubtraction::TGeoSubtraction(const char *expr1, const char *expr2)
          :TGeoBoolNode(expr1, expr2)
{
// Constructor
}

//-----------------------------------------------------------------------------
TGeoSubtraction::TGeoSubtraction(TGeoShape *left, TGeoShape *right, TGeoMatrix *lmat, TGeoMatrix *rmat)
                :TGeoBoolNode(left,right,lmat,rmat)
{
// Constructor providing pointers to components
   if (left->TestShapeBit(TGeoShape::kGeoHalfSpace)) {
      Fatal("TGeoSubstraction", "Substractions from a half-space (%s) not allowed", left->GetName());
   }
}

//-----------------------------------------------------------------------------
TGeoSubtraction::~TGeoSubtraction()
{
// Destructor
// --- deletion of components handled by TGeoManager class.
}

//-----------------------------------------------------------------------------
void TGeoSubtraction::ComputeBBox(Double_t &dx, Double_t &dy, Double_t &dz, Double_t *origin)
{
// Compute bounding box corresponding to a subtraction of two shapes.
   TGeoBBox *box = (TGeoBBox*)fLeft;
   if (box->IsNullBox()) fLeft->ComputeBBox();
   Double_t vert[24];
   Double_t pt[3];
   Int_t i;
   Double_t xmin, xmax, ymin, ymax, zmin, zmax;
   xmin = ymin = zmin = TGeoShape::Big();
   xmax = ymax = zmax = -TGeoShape::Big();
   box->SetBoxPoints(&vert[0]);
   for (i=0; i<8; i++) {
      fLeftMat->LocalToMaster(&vert[3*i], &pt[0]);
      if (pt[0]<xmin) xmin=pt[0];
      if (pt[0]>xmax) xmax=pt[0];
      if (pt[1]<ymin) ymin=pt[1];
      if (pt[1]>ymax) ymax=pt[1];
      if (pt[2]<zmin) zmin=pt[2];
      if (pt[2]>zmax) zmax=pt[2];
   }   
   dx = 0.5*(xmax-xmin);
   origin[0] = 0.5*(xmin+xmax);
   dy = 0.5*(ymax-ymin);
   origin[1] = 0.5*(ymin+ymax);
   dz = 0.5*(zmax-zmin);
   origin[2] = 0.5*(zmin+zmax);
}   

//_____________________________________________________________________________
void TGeoSubtraction::ComputeNormal(Double_t *point, Double_t *dir, Double_t *norm)
{
// Normal computation in POINT. The orientation is chosen so that DIR.dot.NORM>0.
   norm[0] = norm[1] = 0.;
   norm[2] = 1.;
   Double_t local[3], ldir[3], lnorm[3];
   if (fSelected == 1) {
      fLeftMat->MasterToLocal(point, local);
      fLeftMat->MasterToLocalVect(dir, ldir);
      fLeft->ComputeNormal(local,ldir,lnorm);
      fLeftMat->LocalToMasterVect(lnorm, norm);
      return;
   }
   if (fSelected == 2) {
      fRightMat->MasterToLocal(point, local);
      fRightMat->MasterToLocalVect(dir, ldir);
      fRight->ComputeNormal(local,ldir,lnorm);
      fRightMat->LocalToMasterVect(lnorm, norm);
      return;
   }
   fRightMat->MasterToLocal(point,local);
   if (fRight->Contains(local)) {
      fRightMat->MasterToLocalVect(dir,ldir);
      fRight->ComputeNormal(local,ldir, lnorm);
      fRightMat->LocalToMasterVect(lnorm,norm);
      return;
   }   
   fLeftMat->MasterToLocal(point,local);
   if (!fLeft->Contains(local)) {
      fLeftMat->MasterToLocalVect(dir,ldir);
      fLeft->ComputeNormal(local,ldir, lnorm);
      fLeftMat->LocalToMasterVect(lnorm,norm);
      return;
   }
   // point is inside left shape, but not inside the right
   local[0] = point[0]+1E-5*dir[0];
   local[1] = point[1]+1E-5*dir[1];
   local[2] = point[2]+1E-5*dir[2];
   if (Contains(local)) {
      local[0] = point[0]-1E-5*dir[0];
      local[1] = point[1]-1E-5*dir[1];
      local[2] = point[2]-1E-5*dir[2];
      if (Contains(local)) return;
   }  
   ComputeNormal(local,dir,norm);
}

//-----------------------------------------------------------------------------
Bool_t TGeoSubtraction::Contains(Double_t *point) const
{
// Find if a subtraction of two shapes contains a given point
   Double_t local[3];
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   fLeftMat->MasterToLocal(point, &local[0]);
   Bool_t inside = fLeft->Contains(&local[0]);
   if (inside) node->SetSelected(1);
   else return kFALSE;
   fRightMat->MasterToLocal(point, &local[0]);
   inside = !fRight->Contains(&local[0]);
   if (!inside) node->SetSelected(2);
   return inside;
}
//-----------------------------------------------------------------------------
Int_t TGeoSubtraction::DistanceToPrimitive(Int_t /*px*/, Int_t /*py*/)
{
// Compute minimum distance to shape vertices
   return 9999;
}
//-----------------------------------------------------------------------------
Double_t TGeoSubtraction::DistFromInside(Double_t *point, Double_t *dir, Int_t iact,
                              Double_t step, Double_t *safe) const
{
// Compute distance from a given point inside to the shape boundary.
   if (iact<3 && safe) {
      // compute safe distance
      *safe = Safety(point,kTRUE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   Double_t local[3], ldir[3], rdir[3];
   Double_t d1, d2, snxt=0.;
   fLeftMat->MasterToLocal(point, &local[0]);
   fLeftMat->MasterToLocalVect(dir, &ldir[0]);
   fRightMat->MasterToLocalVect(dir, &rdir[0]);
   d1 = fLeft->DistFromInside(&local[0], &ldir[0], iact, step, safe);
   fRightMat->MasterToLocal(point, &local[0]);
   d2 = fRight->DistFromOutside(&local[0], &rdir[0], iact, step, safe);
   if (d1<d2) {
      snxt = d1;
      node->SetSelected(1);
   } else {
      snxt = d2;
      node->SetSelected(2);
   }      
   return snxt;
}   
//-----------------------------------------------------------------------------
Double_t TGeoSubtraction::DistFromOutside(Double_t *point, Double_t *dir, Int_t iact,
                              Double_t step, Double_t *safe) const
{
// Compute distance from a given point outside to the shape.
   if (iact<3 && safe) {
      // compute safe distance
      *safe = Safety(point,kFALSE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   Double_t local[3], master[3], ldir[3], rdir[3];
   memcpy(&master[0], point, 3*sizeof(Double_t));
   Int_t i;
   Double_t d1, d2, snxt=0.;
   fRightMat->MasterToLocal(point, &local[0]);
   fLeftMat->MasterToLocalVect(dir, &ldir[0]);
   fRightMat->MasterToLocalVect(dir, &rdir[0]);
   // check if inside '-'
   Bool_t inside = fRight->Contains(&local[0]);
   Double_t epsil = 0.;
   while (1) {
      if (inside) {
         // propagate to outside of '-'
         node->SetSelected(2);
         d1 = fRight->DistFromInside(&local[0], &rdir[0], iact, step, safe);
         snxt += d1+epsil;
         for (i=0; i<3; i++) master[i] += (d1+1E-8)*dir[i];
         epsil = 1.E-8;
         // now master outside '-'; check if inside '+'
         fLeftMat->MasterToLocal(&master[0], &local[0]);
         if (fLeft->Contains(&local[0])) return snxt;
      } 
      // master outside '-' and outside '+' ;  find distances to both
      node->SetSelected(1);
      fLeftMat->MasterToLocal(&master[0], &local[0]);
      d2 = fLeft->DistFromOutside(&local[0], &ldir[0], iact, step, safe);
      if (d2>1E20) return TGeoShape::Big();
      
      fRightMat->MasterToLocal(&master[0], &local[0]);
      d1 = fRight->DistFromOutside(&local[0], &rdir[0], iact, step, safe);
      if (d2<d1-TGeoShape::Tolerance()) {
         snxt += d2+epsil;
         return snxt;
      }   
      // propagate to '-'
      snxt += d1+epsil;
      for (i=0; i<3; i++) master[i] += (d1+1E-8)*dir[i];
      epsil = 1.E-8;
      // now inside '-' and not inside '+'
      fRightMat->MasterToLocal(&master[0], &local[0]);
      inside = kTRUE;
   }
}
//-----------------------------------------------------------------------------
Int_t TGeoSubtraction::GetNpoints() const
{
// Returns number of vertices for the composite shape described by this subtraction.
   return 0;
}
//-----------------------------------------------------------------------------
Double_t TGeoSubtraction::Safety(Double_t *point, Bool_t in) const
{
// Compute safety distance for a union node;
   Double_t local1[3], local2[3];
   fLeftMat->MasterToLocal(point,local1);
   Bool_t in1 = fLeft->Contains(local1);
   fRightMat->MasterToLocal(point,local2);
   Bool_t in2 = fRight->Contains(local2);
   Bool_t intrue = in1 && (!in2);
   if (in^intrue) return 0.0;
   Double_t saf1 = fLeft->Safety(local1, in1);
   Double_t saf2 = fRight->Safety(local2, in2);
   if (in1 && in2) return saf2;
   if (in1)        return TMath::Min(saf1,saf2);
   if (in2)        return TMath::Max(saf1,saf2);
   return saf1;
}   
//_____________________________________________________________________________
void TGeoSubtraction::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
// Save a primitive as a C++ statement(s) on output stream "out".
   TGeoBoolNode::SavePrimitive(out,option);
   out << "   pBoolNode = new TGeoSubtraction(";
   out << fLeft->GetPointerName() << ",";
   out << fRight->GetPointerName() << ",";
   if (!fLeftMat->IsIdentity()) out << fLeftMat->GetPointerName() << ",";
   else                         out << "0,";
   if (!fRightMat->IsIdentity()) out << fRightMat->GetPointerName() << ");" << endl;
   else                         out << "0);" << endl;
}   
//-----------------------------------------------------------------------------
void TGeoSubtraction::SetPoints(Double_t * /*points*/) const
{
// Fill buffer with shape vertices.
}
//-----------------------------------------------------------------------------
void TGeoSubtraction::SetPoints(Float_t * /*points*/) const
{
// Fill buffer with shape vertices.
}
//-----------------------------------------------------------------------------
void TGeoSubtraction::Sizeof3D() const
{
// Register 3D size of this shape.
   TGeoBoolNode::Sizeof3D();
}
   


ClassImp(TGeoIntersection)

//-----------------------------------------------------------------------------
void TGeoIntersection::Paint(Option_t *option)
{
// Paint method.
   TVirtualViewer3D *viewer = gPad->GetViewer3D();

   if (!viewer) {
      Error("Paint", "gPad->GetViewer3D() returned 0, cannot work with composite!\n");
      return;
   }

   viewer->AddCompositeOp(TBuffer3D::kCSIntersection);

   TGeoBoolNode::Paint(option);
}

//-----------------------------------------------------------------------------
TGeoIntersection::TGeoIntersection()
{
// Default constructor
}

//-----------------------------------------------------------------------------
TGeoIntersection::TGeoIntersection(const char *expr1, const char *expr2)
          :TGeoBoolNode(expr1, expr2)
{
// Constructor
}

//-----------------------------------------------------------------------------
TGeoIntersection::TGeoIntersection(TGeoShape *left, TGeoShape *right, TGeoMatrix *lmat, TGeoMatrix *rmat)
                 :TGeoBoolNode(left,right,lmat,rmat)
{
// Constructor providing pointers to components
   Bool_t hs1 = (fLeft->TestShapeBit(TGeoShape::kGeoHalfSpace))?kTRUE:kFALSE;
   Bool_t hs2 = (fRight->TestShapeBit(TGeoShape::kGeoHalfSpace))?kTRUE:kFALSE;
   if (hs1 && hs2) Fatal("ctor", "cannot intersect two half-spaces: %s * %s", left->GetName(), right->GetName());
}

//-----------------------------------------------------------------------------
TGeoIntersection::~TGeoIntersection()
{
// Destructor
// --- deletion of components handled by TGeoManager class.
}

//-----------------------------------------------------------------------------
void TGeoIntersection::ComputeBBox(Double_t &dx, Double_t &dy, Double_t &dz, Double_t *origin)
{
// Compute bounding box corresponding to a intersection of two shapes.
   Bool_t hs1 = (fLeft->TestShapeBit(TGeoShape::kGeoHalfSpace))?kTRUE:kFALSE;
   Bool_t hs2 = (fRight->TestShapeBit(TGeoShape::kGeoHalfSpace))?kTRUE:kFALSE;
   Double_t vert[48];
   Double_t pt[3];
   Int_t i;
   Double_t xmin1, xmax1, ymin1, ymax1, zmin1, zmax1;
   Double_t xmin2, xmax2, ymin2, ymax2, zmin2, zmax2;
   xmin1 = ymin1 = zmin1 = xmin2 = ymin2 = zmin2 = TGeoShape::Big();
   xmax1 = ymax1 = zmax1 = xmax2 = ymax2 = zmax2 =  -TGeoShape::Big();
   if (!hs1) {
      if (((TGeoBBox*)fLeft)->IsNullBox()) fLeft->ComputeBBox();
      ((TGeoBBox*)fLeft)->SetBoxPoints(&vert[0]);
      for (i=0; i<8; i++) {
         fLeftMat->LocalToMaster(&vert[3*i], &pt[0]);
         if (pt[0]<xmin1) xmin1=pt[0];
         if (pt[0]>xmax1) xmax1=pt[0];
         if (pt[1]<ymin1) ymin1=pt[1];
         if (pt[1]>ymax1) ymax1=pt[1];
         if (pt[2]<zmin1) zmin1=pt[2];
         if (pt[2]>zmax1) zmax1=pt[2];
      }   
   }   
   if (!hs2) {
      if (((TGeoBBox*)fRight)->IsNullBox()) fRight->ComputeBBox();
      ((TGeoBBox*)fRight)->SetBoxPoints(&vert[24]);
      for (i=8; i<16; i++) {
         fRightMat->LocalToMaster(&vert[3*i], &pt[0]);
         if (pt[0]<xmin2) xmin2=pt[0];
         if (pt[0]>xmax2) xmax2=pt[0];
         if (pt[1]<ymin2) ymin2=pt[1];
         if (pt[1]>ymax2) ymax2=pt[1];
         if (pt[2]<zmin2) zmin2=pt[2];
         if (pt[2]>zmax2) zmax2=pt[2];
      }
   }      
   if (hs1) {
      dx = 0.5*(xmax2-xmin2);
      origin[0] = 0.5*(xmax2+xmin2);   
      dy = 0.5*(ymax2-ymin2);
      origin[1] = 0.5*(ymax2+ymin2);   
      dz = 0.5*(zmax2-zmin2);
      origin[2] = 0.5*(zmax2+zmin2);   
      return;
   }            
   if (hs2) {
      dx = 0.5*(xmax1-xmin1);
      origin[0] = 0.5*(xmax1+xmin1);   
      dy = 0.5*(ymax1-ymin1);
      origin[1] = 0.5*(ymax1+ymin1);   
      dz = 0.5*(zmax1-zmin1);
      origin[2] = 0.5*(zmax1+zmin1);   
      return;
   }   
   Double_t sort[4];
   Int_t isort[4];
   sort[0] = xmin1;
   sort[1] = xmax1;
   sort[2] = xmin2;
   sort[3] = xmax2;
   TMath::Sort(4, &sort[0], &isort[0], kFALSE);
   if (isort[1]%2) {
      Warning("ComputeBBox", "shapes %s and %s do not intersect", fLeft->GetName(), fRight->GetName());
      dx = dy = dz = 0;
      memset(origin, 0, 3*sizeof(Double_t));
      return;
   }
   dx = 0.5*(sort[isort[2]]-sort[isort[1]]);
   origin[0] = 0.5*(sort[isort[1]]+sort[isort[2]]);   
   sort[0] = ymin1;
   sort[1] = ymax1;
   sort[2] = ymin2;
   sort[3] = ymax2;
   TMath::Sort(4, &sort[0], &isort[0], kFALSE);
   if (isort[1]%2) {
      Warning("ComputeBBox", "shapes %s and %s do not intersect", fLeft->GetName(), fRight->GetName());
      dx = dy = dz = 0;
      memset(origin, 0, 3*sizeof(Double_t));
      return;
   }
   dy = 0.5*(sort[isort[2]]-sort[isort[1]]);
   origin[1] = 0.5*(sort[isort[1]]+sort[isort[2]]);   
   sort[0] = zmin1;
   sort[1] = zmax1;
   sort[2] = zmin2;
   sort[3] = zmax2;
   TMath::Sort(4, &sort[0], &isort[0], kFALSE);
   if (isort[1]%2) {
      Warning("ComputeBBox", "shapes %s and %s do not intersect", fLeft->GetName(), fRight->GetName());
      dx = dy = dz = 0;
      memset(origin, 0, 3*sizeof(Double_t));
      return;
   }
   dz = 0.5*(sort[isort[2]]-sort[isort[1]]);
   origin[2] = 0.5*(sort[isort[1]]+sort[isort[2]]);   
}   

//_____________________________________________________________________________
void TGeoIntersection::ComputeNormal(Double_t *point, Double_t *dir, Double_t *norm)
{
// Normal computation in POINT. The orientation is chosen so that DIR.dot.NORM>0.
   Double_t local[3], ldir[3], lnorm[3];
   norm[0] = norm[1] = 0.;
   norm[2] = 1.;
   if (fSelected == 1) {
      fLeftMat->MasterToLocal(point, local);
      fLeftMat->MasterToLocalVect(dir, ldir);
      fLeft->ComputeNormal(local,ldir,lnorm);
      fLeftMat->LocalToMasterVect(lnorm, norm);
      return;
   }
   if (fSelected == 2) {
      fRightMat->MasterToLocal(point, local);
      fRightMat->MasterToLocalVect(dir, ldir);
      fRight->ComputeNormal(local,ldir,lnorm);
      fRightMat->LocalToMasterVect(lnorm, norm);
      return;
   }            
   fLeftMat->MasterToLocal(point,local);
   if (!fLeft->Contains(local)) {
      fLeftMat->MasterToLocalVect(dir,ldir);
      fLeft->ComputeNormal(local,ldir,lnorm);
      fLeftMat->LocalToMasterVect(lnorm,norm);
      return;
   }
   fRightMat->MasterToLocal(point,local);
   if (!fRight->Contains(local)) {
      fRightMat->MasterToLocalVect(dir,ldir);
      fRight->ComputeNormal(local,ldir,lnorm);
      fRightMat->LocalToMasterVect(lnorm,norm);
      return;
   }
   // point is inside intersection.
   local[0] = point[0] + 1E-5*dir[0];   
   local[1] = point[1] + 1E-5*dir[1];   
   local[2] = point[2] + 1E-5*dir[2];
   if (Contains(local)) {
      local[0] = point[0] - 1E-5*dir[0];   
      local[1] = point[1] - 1E-5*dir[1];   
      local[2] = point[2] - 1E-5*dir[2];
      if (Contains(local)) return;
   }
   ComputeNormal(local,dir,norm);   
}

//-----------------------------------------------------------------------------
Bool_t TGeoIntersection::Contains(Double_t *point) const
{
// Find if a intersection of two shapes contains a given point
   Double_t local[3];
   fLeftMat->MasterToLocal(point, &local[0]);
   Bool_t inside = fLeft->Contains(&local[0]);
   if (!inside) return kFALSE;
   fRightMat->MasterToLocal(point, &local[0]);
   inside = fRight->Contains(&local[0]);
   return inside;
}
//-----------------------------------------------------------------------------
Int_t TGeoIntersection::DistanceToPrimitive(Int_t /*px*/, Int_t /*py*/)
{
// Compute minimum distance to shape vertices
   return 9999;
}
//-----------------------------------------------------------------------------
Double_t TGeoIntersection::DistFromInside(Double_t *point, Double_t *dir, Int_t iact,
                              Double_t step, Double_t *safe) const
{
// Compute distance from a given point inside to the shape boundary.
   if (iact<3 && safe) {
      // compute safe distance
      *safe = Safety(point,kTRUE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   Double_t local[3], ldir[3], rdir[3];
   Double_t d1, d2, snxt=0.;
   fLeftMat->MasterToLocal(point, &local[0]);
   fLeftMat->MasterToLocalVect(dir, &ldir[0]);
   fRightMat->MasterToLocalVect(dir, &rdir[0]);
   d1 = fLeft->DistFromInside(&local[0], &ldir[0], iact, step, safe);
   fRightMat->MasterToLocal(point, &local[0]);
   d2 = fRight->DistFromInside(&local[0], &rdir[0], iact, step, safe);
   if (d1<d2) {
      snxt = d1;
      node->SetSelected(1);
   } else {
      snxt = d2;
      node->SetSelected(2);
   }      
   return snxt;
}   
//-----------------------------------------------------------------------------
Double_t TGeoIntersection::DistFromOutside(Double_t *point, Double_t *dir, Int_t iact,
                              Double_t step, Double_t *safe) const
{
// Compute distance from a given point outside to the shape.
   if (iact<3 && safe) {
      // compute safe distance
      *safe = Safety(point,kFALSE);
      if (iact==0) return TGeoShape::Big();
      if (iact==1 && step<*safe) return TGeoShape::Big();
   }
   TGeoBoolNode *node = (TGeoBoolNode*)this;
   Double_t lpt[3], rpt[3], master[3], ldir[3], rdir[3];
   memcpy(master, point, 3*sizeof(Double_t));
   Int_t i;
   Double_t d1 = 0.;
   Double_t d2 = 0.;
   fLeftMat->MasterToLocal(point, lpt);
   fRightMat->MasterToLocal(point, rpt);
   fLeftMat->MasterToLocalVect(dir, ldir);
   fRightMat->MasterToLocalVect(dir, rdir);
   Bool_t inleft = fLeft->Contains(lpt);
   Bool_t inright = fRight->Contains(rpt);
   node->SetSelected(0);
   Double_t snext = 0.0;
   if (inleft && inright) return snext;

   while (1) {
      d1 = d2 = 0.0;
      if (!inleft)  {
         d1 = fLeft->DistFromOutside(lpt,ldir,3);
         if (d1>1E20) return TGeoShape::Big();
      }
      if (!inright) {  
         d2 = fRight->DistFromOutside(rpt,rdir,3);
         if (d2>1E20) return TGeoShape::Big();
      }
   
      if (d1>d2) {
         // propagate to left shape
         snext += d1;
         node->SetSelected(1);
         inleft = kTRUE;
         for (i=0; i<3; i++) master[i] += d1*dir[i];
         fRightMat->MasterToLocal(master,rpt);
         // check if propagated point is inside right shape
         inright = fRight->Contains(rpt);
         if (inright) return snext;
         // here inleft=true, inright=false         
      } else {
         // propagate to right shape
         snext += d2;
         node->SetSelected(2);
         inright = kTRUE;
         for (i=0; i<3; i++) master[i] += d2*dir[i];
         fLeftMat->MasterToLocal(master,lpt);
         // check if propagated point is inside left shape
         inleft = fLeft->Contains(lpt);
         if (inleft) return snext;
         // here inleft=false, inright=true
      }            
   }   
   return snext;
}      

//-----------------------------------------------------------------------------
Int_t TGeoIntersection::GetNpoints() const
{
// Returns number of vertices for the composite shape described by this intersection.
   return 0;
}
//-----------------------------------------------------------------------------
Double_t TGeoIntersection::Safety(Double_t *point, Bool_t in) const
{
// Compute safety distance for a union node;
   Double_t local1[3], local2[3];
   fLeftMat->MasterToLocal(point,local1);
   Bool_t in1 = fLeft->Contains(local1);
   fRightMat->MasterToLocal(point,local2);
   Bool_t in2 = fRight->Contains(local2);
   Bool_t intrue = in1 & in2;
   if (in^intrue) return 0.0;
   Double_t saf1 = fLeft->Safety(local1, in1);
   Double_t saf2 = fRight->Safety(local2, in2);
   if (in1 && in2) return TMath::Min(saf1, saf2);
   if (in1)        return saf2;
   if (in2)        return saf1;
   return TMath::Max(saf1,saf2);
}   
//_____________________________________________________________________________
void TGeoIntersection::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
// Save a primitive as a C++ statement(s) on output stream "out".
   TGeoBoolNode::SavePrimitive(out,option);
   out << "   pBoolNode = new TGeoIntersection(";
   out << fLeft->GetPointerName() << ",";
   out << fRight->GetPointerName() << ",";
   if (!fLeftMat->IsIdentity()) out << fLeftMat->GetPointerName() << ",";
   else                         out << "0,";
   if (!fRightMat->IsIdentity()) out << fRightMat->GetPointerName() << ");" << endl;
   else                         out << "0);" << endl;
}   
//-----------------------------------------------------------------------------
void TGeoIntersection::SetPoints(Double_t * /*points*/) const
{
// Fill buffer with shape vertices.
}
//-----------------------------------------------------------------------------
void TGeoIntersection::SetPoints(Float_t * /*points*/) const
{
// Fill buffer with shape vertices.
}
//-----------------------------------------------------------------------------
void TGeoIntersection::Sizeof3D() const
{
// Register 3D size of this shape.
   TGeoBoolNode::Sizeof3D();
}
 TGeoBoolNode.cxx:1
 TGeoBoolNode.cxx:2
 TGeoBoolNode.cxx:3
 TGeoBoolNode.cxx:4
 TGeoBoolNode.cxx:5
 TGeoBoolNode.cxx:6
 TGeoBoolNode.cxx:7
 TGeoBoolNode.cxx:8
 TGeoBoolNode.cxx:9
 TGeoBoolNode.cxx:10
 TGeoBoolNode.cxx:11
 TGeoBoolNode.cxx:12
 TGeoBoolNode.cxx:13
 TGeoBoolNode.cxx:14
 TGeoBoolNode.cxx:15
 TGeoBoolNode.cxx:16
 TGeoBoolNode.cxx:17
 TGeoBoolNode.cxx:18
 TGeoBoolNode.cxx:19
 TGeoBoolNode.cxx:20
 TGeoBoolNode.cxx:21
 TGeoBoolNode.cxx:22
 TGeoBoolNode.cxx:23
 TGeoBoolNode.cxx:24
 TGeoBoolNode.cxx:25
 TGeoBoolNode.cxx:26
 TGeoBoolNode.cxx:27
 TGeoBoolNode.cxx:28
 TGeoBoolNode.cxx:29
 TGeoBoolNode.cxx:30
 TGeoBoolNode.cxx:31
 TGeoBoolNode.cxx:32
 TGeoBoolNode.cxx:33
 TGeoBoolNode.cxx:34
 TGeoBoolNode.cxx:35
 TGeoBoolNode.cxx:36
 TGeoBoolNode.cxx:37
 TGeoBoolNode.cxx:38
 TGeoBoolNode.cxx:39
 TGeoBoolNode.cxx:40
 TGeoBoolNode.cxx:41
 TGeoBoolNode.cxx:42
 TGeoBoolNode.cxx:43
 TGeoBoolNode.cxx:44
 TGeoBoolNode.cxx:45
 TGeoBoolNode.cxx:46
 TGeoBoolNode.cxx:47
 TGeoBoolNode.cxx:48
 TGeoBoolNode.cxx:49
 TGeoBoolNode.cxx:50
 TGeoBoolNode.cxx:51
 TGeoBoolNode.cxx:52
 TGeoBoolNode.cxx:53
 TGeoBoolNode.cxx:54
 TGeoBoolNode.cxx:55
 TGeoBoolNode.cxx:56
 TGeoBoolNode.cxx:57
 TGeoBoolNode.cxx:58
 TGeoBoolNode.cxx:59
 TGeoBoolNode.cxx:60
 TGeoBoolNode.cxx:61
 TGeoBoolNode.cxx:62
 TGeoBoolNode.cxx:63
 TGeoBoolNode.cxx:64
 TGeoBoolNode.cxx:65
 TGeoBoolNode.cxx:66
 TGeoBoolNode.cxx:67
 TGeoBoolNode.cxx:68
 TGeoBoolNode.cxx:69
 TGeoBoolNode.cxx:70
 TGeoBoolNode.cxx:71
 TGeoBoolNode.cxx:72
 TGeoBoolNode.cxx:73
 TGeoBoolNode.cxx:74
 TGeoBoolNode.cxx:75
 TGeoBoolNode.cxx:76
 TGeoBoolNode.cxx:77
 TGeoBoolNode.cxx:78
 TGeoBoolNode.cxx:79
 TGeoBoolNode.cxx:80
 TGeoBoolNode.cxx:81
 TGeoBoolNode.cxx:82
 TGeoBoolNode.cxx:83
 TGeoBoolNode.cxx:84
 TGeoBoolNode.cxx:85
 TGeoBoolNode.cxx:86
 TGeoBoolNode.cxx:87
 TGeoBoolNode.cxx:88
 TGeoBoolNode.cxx:89
 TGeoBoolNode.cxx:90
 TGeoBoolNode.cxx:91
 TGeoBoolNode.cxx:92
 TGeoBoolNode.cxx:93
 TGeoBoolNode.cxx:94
 TGeoBoolNode.cxx:95
 TGeoBoolNode.cxx:96
 TGeoBoolNode.cxx:97
 TGeoBoolNode.cxx:98
 TGeoBoolNode.cxx:99
 TGeoBoolNode.cxx:100
 TGeoBoolNode.cxx:101
 TGeoBoolNode.cxx:102
 TGeoBoolNode.cxx:103
 TGeoBoolNode.cxx:104
 TGeoBoolNode.cxx:105
 TGeoBoolNode.cxx:106
 TGeoBoolNode.cxx:107
 TGeoBoolNode.cxx:108
 TGeoBoolNode.cxx:109
 TGeoBoolNode.cxx:110
 TGeoBoolNode.cxx:111
 TGeoBoolNode.cxx:112
 TGeoBoolNode.cxx:113
 TGeoBoolNode.cxx:114
 TGeoBoolNode.cxx:115
 TGeoBoolNode.cxx:116
 TGeoBoolNode.cxx:117
 TGeoBoolNode.cxx:118
 TGeoBoolNode.cxx:119
 TGeoBoolNode.cxx:120
 TGeoBoolNode.cxx:121
 TGeoBoolNode.cxx:122
 TGeoBoolNode.cxx:123
 TGeoBoolNode.cxx:124
 TGeoBoolNode.cxx:125
 TGeoBoolNode.cxx:126
 TGeoBoolNode.cxx:127
 TGeoBoolNode.cxx:128
 TGeoBoolNode.cxx:129
 TGeoBoolNode.cxx:130
 TGeoBoolNode.cxx:131
 TGeoBoolNode.cxx:132
 TGeoBoolNode.cxx:133
 TGeoBoolNode.cxx:134
 TGeoBoolNode.cxx:135
 TGeoBoolNode.cxx:136
 TGeoBoolNode.cxx:137
 TGeoBoolNode.cxx:138
 TGeoBoolNode.cxx:139
 TGeoBoolNode.cxx:140
 TGeoBoolNode.cxx:141
 TGeoBoolNode.cxx:142
 TGeoBoolNode.cxx:143
 TGeoBoolNode.cxx:144
 TGeoBoolNode.cxx:145
 TGeoBoolNode.cxx:146
 TGeoBoolNode.cxx:147
 TGeoBoolNode.cxx:148
 TGeoBoolNode.cxx:149
 TGeoBoolNode.cxx:150
 TGeoBoolNode.cxx:151
 TGeoBoolNode.cxx:152
 TGeoBoolNode.cxx:153
 TGeoBoolNode.cxx:154
 TGeoBoolNode.cxx:155
 TGeoBoolNode.cxx:156
 TGeoBoolNode.cxx:157
 TGeoBoolNode.cxx:158
 TGeoBoolNode.cxx:159
 TGeoBoolNode.cxx:160
 TGeoBoolNode.cxx:161
 TGeoBoolNode.cxx:162
 TGeoBoolNode.cxx:163
 TGeoBoolNode.cxx:164
 TGeoBoolNode.cxx:165
 TGeoBoolNode.cxx:166
 TGeoBoolNode.cxx:167
 TGeoBoolNode.cxx:168
 TGeoBoolNode.cxx:169
 TGeoBoolNode.cxx:170
 TGeoBoolNode.cxx:171
 TGeoBoolNode.cxx:172
 TGeoBoolNode.cxx:173
 TGeoBoolNode.cxx:174
 TGeoBoolNode.cxx:175
 TGeoBoolNode.cxx:176
 TGeoBoolNode.cxx:177
 TGeoBoolNode.cxx:178
 TGeoBoolNode.cxx:179
 TGeoBoolNode.cxx:180
 TGeoBoolNode.cxx:181
 TGeoBoolNode.cxx:182
 TGeoBoolNode.cxx:183
 TGeoBoolNode.cxx:184
 TGeoBoolNode.cxx:185
 TGeoBoolNode.cxx:186
 TGeoBoolNode.cxx:187
 TGeoBoolNode.cxx:188
 TGeoBoolNode.cxx:189
 TGeoBoolNode.cxx:190
 TGeoBoolNode.cxx:191
 TGeoBoolNode.cxx:192
 TGeoBoolNode.cxx:193
 TGeoBoolNode.cxx:194
 TGeoBoolNode.cxx:195
 TGeoBoolNode.cxx:196
 TGeoBoolNode.cxx:197
 TGeoBoolNode.cxx:198
 TGeoBoolNode.cxx:199
 TGeoBoolNode.cxx:200
 TGeoBoolNode.cxx:201
 TGeoBoolNode.cxx:202
 TGeoBoolNode.cxx:203
 TGeoBoolNode.cxx:204
 TGeoBoolNode.cxx:205
 TGeoBoolNode.cxx:206
 TGeoBoolNode.cxx:207
 TGeoBoolNode.cxx:208
 TGeoBoolNode.cxx:209
 TGeoBoolNode.cxx:210
 TGeoBoolNode.cxx:211
 TGeoBoolNode.cxx:212
 TGeoBoolNode.cxx:213
 TGeoBoolNode.cxx:214
 TGeoBoolNode.cxx:215
 TGeoBoolNode.cxx:216
 TGeoBoolNode.cxx:217
 TGeoBoolNode.cxx:218
 TGeoBoolNode.cxx:219
 TGeoBoolNode.cxx:220
 TGeoBoolNode.cxx:221
 TGeoBoolNode.cxx:222
 TGeoBoolNode.cxx:223
 TGeoBoolNode.cxx:224
 TGeoBoolNode.cxx:225
 TGeoBoolNode.cxx:226
 TGeoBoolNode.cxx:227
 TGeoBoolNode.cxx:228
 TGeoBoolNode.cxx:229
 TGeoBoolNode.cxx:230
 TGeoBoolNode.cxx:231
 TGeoBoolNode.cxx:232
 TGeoBoolNode.cxx:233
 TGeoBoolNode.cxx:234
 TGeoBoolNode.cxx:235
 TGeoBoolNode.cxx:236
 TGeoBoolNode.cxx:237
 TGeoBoolNode.cxx:238
 TGeoBoolNode.cxx:239
 TGeoBoolNode.cxx:240
 TGeoBoolNode.cxx:241
 TGeoBoolNode.cxx:242
 TGeoBoolNode.cxx:243
 TGeoBoolNode.cxx:244
 TGeoBoolNode.cxx:245
 TGeoBoolNode.cxx:246
 TGeoBoolNode.cxx:247
 TGeoBoolNode.cxx:248
 TGeoBoolNode.cxx:249
 TGeoBoolNode.cxx:250
 TGeoBoolNode.cxx:251
 TGeoBoolNode.cxx:252
 TGeoBoolNode.cxx:253
 TGeoBoolNode.cxx:254
 TGeoBoolNode.cxx:255
 TGeoBoolNode.cxx:256
 TGeoBoolNode.cxx:257
 TGeoBoolNode.cxx:258
 TGeoBoolNode.cxx:259
 TGeoBoolNode.cxx:260
 TGeoBoolNode.cxx:261
 TGeoBoolNode.cxx:262
 TGeoBoolNode.cxx:263
 TGeoBoolNode.cxx:264
 TGeoBoolNode.cxx:265
 TGeoBoolNode.cxx:266
 TGeoBoolNode.cxx:267
 TGeoBoolNode.cxx:268
 TGeoBoolNode.cxx:269
 TGeoBoolNode.cxx:270
 TGeoBoolNode.cxx:271
 TGeoBoolNode.cxx:272
 TGeoBoolNode.cxx:273
 TGeoBoolNode.cxx:274
 TGeoBoolNode.cxx:275
 TGeoBoolNode.cxx:276
 TGeoBoolNode.cxx:277
 TGeoBoolNode.cxx:278
 TGeoBoolNode.cxx:279
 TGeoBoolNode.cxx:280
 TGeoBoolNode.cxx:281
 TGeoBoolNode.cxx:282
 TGeoBoolNode.cxx:283
 TGeoBoolNode.cxx:284
 TGeoBoolNode.cxx:285
 TGeoBoolNode.cxx:286
 TGeoBoolNode.cxx:287
 TGeoBoolNode.cxx:288
 TGeoBoolNode.cxx:289
 TGeoBoolNode.cxx:290
 TGeoBoolNode.cxx:291
 TGeoBoolNode.cxx:292
 TGeoBoolNode.cxx:293
 TGeoBoolNode.cxx:294
 TGeoBoolNode.cxx:295
 TGeoBoolNode.cxx:296
 TGeoBoolNode.cxx:297
 TGeoBoolNode.cxx:298
 TGeoBoolNode.cxx:299
 TGeoBoolNode.cxx:300
 TGeoBoolNode.cxx:301
 TGeoBoolNode.cxx:302
 TGeoBoolNode.cxx:303
 TGeoBoolNode.cxx:304
 TGeoBoolNode.cxx:305
 TGeoBoolNode.cxx:306
 TGeoBoolNode.cxx:307
 TGeoBoolNode.cxx:308
 TGeoBoolNode.cxx:309
 TGeoBoolNode.cxx:310
 TGeoBoolNode.cxx:311
 TGeoBoolNode.cxx:312
 TGeoBoolNode.cxx:313
 TGeoBoolNode.cxx:314
 TGeoBoolNode.cxx:315
 TGeoBoolNode.cxx:316
 TGeoBoolNode.cxx:317
 TGeoBoolNode.cxx:318
 TGeoBoolNode.cxx:319
 TGeoBoolNode.cxx:320
 TGeoBoolNode.cxx:321
 TGeoBoolNode.cxx:322
 TGeoBoolNode.cxx:323
 TGeoBoolNode.cxx:324
 TGeoBoolNode.cxx:325
 TGeoBoolNode.cxx:326
 TGeoBoolNode.cxx:327
 TGeoBoolNode.cxx:328
 TGeoBoolNode.cxx:329
 TGeoBoolNode.cxx:330
 TGeoBoolNode.cxx:331
 TGeoBoolNode.cxx:332
 TGeoBoolNode.cxx:333
 TGeoBoolNode.cxx:334
 TGeoBoolNode.cxx:335
 TGeoBoolNode.cxx:336
 TGeoBoolNode.cxx:337
 TGeoBoolNode.cxx:338
 TGeoBoolNode.cxx:339
 TGeoBoolNode.cxx:340
 TGeoBoolNode.cxx:341
 TGeoBoolNode.cxx:342
 TGeoBoolNode.cxx:343
 TGeoBoolNode.cxx:344
 TGeoBoolNode.cxx:345
 TGeoBoolNode.cxx:346
 TGeoBoolNode.cxx:347
 TGeoBoolNode.cxx:348
 TGeoBoolNode.cxx:349
 TGeoBoolNode.cxx:350
 TGeoBoolNode.cxx:351
 TGeoBoolNode.cxx:352
 TGeoBoolNode.cxx:353
 TGeoBoolNode.cxx:354
 TGeoBoolNode.cxx:355
 TGeoBoolNode.cxx:356
 TGeoBoolNode.cxx:357
 TGeoBoolNode.cxx:358
 TGeoBoolNode.cxx:359
 TGeoBoolNode.cxx:360
 TGeoBoolNode.cxx:361
 TGeoBoolNode.cxx:362
 TGeoBoolNode.cxx:363
 TGeoBoolNode.cxx:364
 TGeoBoolNode.cxx:365
 TGeoBoolNode.cxx:366
 TGeoBoolNode.cxx:367
 TGeoBoolNode.cxx:368
 TGeoBoolNode.cxx:369
 TGeoBoolNode.cxx:370
 TGeoBoolNode.cxx:371
 TGeoBoolNode.cxx:372
 TGeoBoolNode.cxx:373
 TGeoBoolNode.cxx:374
 TGeoBoolNode.cxx:375
 TGeoBoolNode.cxx:376
 TGeoBoolNode.cxx:377
 TGeoBoolNode.cxx:378
 TGeoBoolNode.cxx:379
 TGeoBoolNode.cxx:380
 TGeoBoolNode.cxx:381
 TGeoBoolNode.cxx:382
 TGeoBoolNode.cxx:383
 TGeoBoolNode.cxx:384
 TGeoBoolNode.cxx:385
 TGeoBoolNode.cxx:386
 TGeoBoolNode.cxx:387
 TGeoBoolNode.cxx:388
 TGeoBoolNode.cxx:389
 TGeoBoolNode.cxx:390
 TGeoBoolNode.cxx:391
 TGeoBoolNode.cxx:392
 TGeoBoolNode.cxx:393
 TGeoBoolNode.cxx:394
 TGeoBoolNode.cxx:395
 TGeoBoolNode.cxx:396
 TGeoBoolNode.cxx:397
 TGeoBoolNode.cxx:398
 TGeoBoolNode.cxx:399
 TGeoBoolNode.cxx:400
 TGeoBoolNode.cxx:401
 TGeoBoolNode.cxx:402
 TGeoBoolNode.cxx:403
 TGeoBoolNode.cxx:404
 TGeoBoolNode.cxx:405
 TGeoBoolNode.cxx:406
 TGeoBoolNode.cxx:407
 TGeoBoolNode.cxx:408
 TGeoBoolNode.cxx:409
 TGeoBoolNode.cxx:410
 TGeoBoolNode.cxx:411
 TGeoBoolNode.cxx:412
 TGeoBoolNode.cxx:413
 TGeoBoolNode.cxx:414
 TGeoBoolNode.cxx:415
 TGeoBoolNode.cxx:416
 TGeoBoolNode.cxx:417
 TGeoBoolNode.cxx:418
 TGeoBoolNode.cxx:419
 TGeoBoolNode.cxx:420
 TGeoBoolNode.cxx:421
 TGeoBoolNode.cxx:422
 TGeoBoolNode.cxx:423
 TGeoBoolNode.cxx:424
 TGeoBoolNode.cxx:425
 TGeoBoolNode.cxx:426
 TGeoBoolNode.cxx:427
 TGeoBoolNode.cxx:428
 TGeoBoolNode.cxx:429
 TGeoBoolNode.cxx:430
 TGeoBoolNode.cxx:431
 TGeoBoolNode.cxx:432
 TGeoBoolNode.cxx:433
 TGeoBoolNode.cxx:434
 TGeoBoolNode.cxx:435
 TGeoBoolNode.cxx:436
 TGeoBoolNode.cxx:437
 TGeoBoolNode.cxx:438
 TGeoBoolNode.cxx:439
 TGeoBoolNode.cxx:440
 TGeoBoolNode.cxx:441
 TGeoBoolNode.cxx:442
 TGeoBoolNode.cxx:443
 TGeoBoolNode.cxx:444
 TGeoBoolNode.cxx:445
 TGeoBoolNode.cxx:446
 TGeoBoolNode.cxx:447
 TGeoBoolNode.cxx:448
 TGeoBoolNode.cxx:449
 TGeoBoolNode.cxx:450
 TGeoBoolNode.cxx:451
 TGeoBoolNode.cxx:452
 TGeoBoolNode.cxx:453
 TGeoBoolNode.cxx:454
 TGeoBoolNode.cxx:455
 TGeoBoolNode.cxx:456
 TGeoBoolNode.cxx:457
 TGeoBoolNode.cxx:458
 TGeoBoolNode.cxx:459
 TGeoBoolNode.cxx:460
 TGeoBoolNode.cxx:461
 TGeoBoolNode.cxx:462
 TGeoBoolNode.cxx:463
 TGeoBoolNode.cxx:464
 TGeoBoolNode.cxx:465
 TGeoBoolNode.cxx:466
 TGeoBoolNode.cxx:467
 TGeoBoolNode.cxx:468
 TGeoBoolNode.cxx:469
 TGeoBoolNode.cxx:470
 TGeoBoolNode.cxx:471
 TGeoBoolNode.cxx:472
 TGeoBoolNode.cxx:473
 TGeoBoolNode.cxx:474
 TGeoBoolNode.cxx:475
 TGeoBoolNode.cxx:476
 TGeoBoolNode.cxx:477
 TGeoBoolNode.cxx:478
 TGeoBoolNode.cxx:479
 TGeoBoolNode.cxx:480
 TGeoBoolNode.cxx:481
 TGeoBoolNode.cxx:482
 TGeoBoolNode.cxx:483
 TGeoBoolNode.cxx:484
 TGeoBoolNode.cxx:485
 TGeoBoolNode.cxx:486
 TGeoBoolNode.cxx:487
 TGeoBoolNode.cxx:488
 TGeoBoolNode.cxx:489
 TGeoBoolNode.cxx:490
 TGeoBoolNode.cxx:491
 TGeoBoolNode.cxx:492
 TGeoBoolNode.cxx:493
 TGeoBoolNode.cxx:494
 TGeoBoolNode.cxx:495
 TGeoBoolNode.cxx:496
 TGeoBoolNode.cxx:497
 TGeoBoolNode.cxx:498
 TGeoBoolNode.cxx:499
 TGeoBoolNode.cxx:500
 TGeoBoolNode.cxx:501
 TGeoBoolNode.cxx:502
 TGeoBoolNode.cxx:503
 TGeoBoolNode.cxx:504
 TGeoBoolNode.cxx:505
 TGeoBoolNode.cxx:506
 TGeoBoolNode.cxx:507
 TGeoBoolNode.cxx:508
 TGeoBoolNode.cxx:509
 TGeoBoolNode.cxx:510
 TGeoBoolNode.cxx:511
 TGeoBoolNode.cxx:512
 TGeoBoolNode.cxx:513
 TGeoBoolNode.cxx:514
 TGeoBoolNode.cxx:515
 TGeoBoolNode.cxx:516
 TGeoBoolNode.cxx:517
 TGeoBoolNode.cxx:518
 TGeoBoolNode.cxx:519
 TGeoBoolNode.cxx:520
 TGeoBoolNode.cxx:521
 TGeoBoolNode.cxx:522
 TGeoBoolNode.cxx:523
 TGeoBoolNode.cxx:524
 TGeoBoolNode.cxx:525
 TGeoBoolNode.cxx:526
 TGeoBoolNode.cxx:527
 TGeoBoolNode.cxx:528
 TGeoBoolNode.cxx:529
 TGeoBoolNode.cxx:530
 TGeoBoolNode.cxx:531
 TGeoBoolNode.cxx:532
 TGeoBoolNode.cxx:533
 TGeoBoolNode.cxx:534
 TGeoBoolNode.cxx:535
 TGeoBoolNode.cxx:536
 TGeoBoolNode.cxx:537
 TGeoBoolNode.cxx:538
 TGeoBoolNode.cxx:539
 TGeoBoolNode.cxx:540
 TGeoBoolNode.cxx:541
 TGeoBoolNode.cxx:542
 TGeoBoolNode.cxx:543
 TGeoBoolNode.cxx:544
 TGeoBoolNode.cxx:545
 TGeoBoolNode.cxx:546
 TGeoBoolNode.cxx:547
 TGeoBoolNode.cxx:548
 TGeoBoolNode.cxx:549
 TGeoBoolNode.cxx:550
 TGeoBoolNode.cxx:551
 TGeoBoolNode.cxx:552
 TGeoBoolNode.cxx:553
 TGeoBoolNode.cxx:554
 TGeoBoolNode.cxx:555
 TGeoBoolNode.cxx:556
 TGeoBoolNode.cxx:557
 TGeoBoolNode.cxx:558
 TGeoBoolNode.cxx:559
 TGeoBoolNode.cxx:560
 TGeoBoolNode.cxx:561
 TGeoBoolNode.cxx:562
 TGeoBoolNode.cxx:563
 TGeoBoolNode.cxx:564
 TGeoBoolNode.cxx:565
 TGeoBoolNode.cxx:566
 TGeoBoolNode.cxx:567
 TGeoBoolNode.cxx:568
 TGeoBoolNode.cxx:569
 TGeoBoolNode.cxx:570
 TGeoBoolNode.cxx:571
 TGeoBoolNode.cxx:572
 TGeoBoolNode.cxx:573
 TGeoBoolNode.cxx:574
 TGeoBoolNode.cxx:575
 TGeoBoolNode.cxx:576
 TGeoBoolNode.cxx:577
 TGeoBoolNode.cxx:578
 TGeoBoolNode.cxx:579
 TGeoBoolNode.cxx:580
 TGeoBoolNode.cxx:581
 TGeoBoolNode.cxx:582
 TGeoBoolNode.cxx:583
 TGeoBoolNode.cxx:584
 TGeoBoolNode.cxx:585
 TGeoBoolNode.cxx:586
 TGeoBoolNode.cxx:587
 TGeoBoolNode.cxx:588
 TGeoBoolNode.cxx:589
 TGeoBoolNode.cxx:590
 TGeoBoolNode.cxx:591
 TGeoBoolNode.cxx:592
 TGeoBoolNode.cxx:593
 TGeoBoolNode.cxx:594
 TGeoBoolNode.cxx:595
 TGeoBoolNode.cxx:596
 TGeoBoolNode.cxx:597
 TGeoBoolNode.cxx:598
 TGeoBoolNode.cxx:599
 TGeoBoolNode.cxx:600
 TGeoBoolNode.cxx:601
 TGeoBoolNode.cxx:602
 TGeoBoolNode.cxx:603
 TGeoBoolNode.cxx:604
 TGeoBoolNode.cxx:605
 TGeoBoolNode.cxx:606
 TGeoBoolNode.cxx:607
 TGeoBoolNode.cxx:608
 TGeoBoolNode.cxx:609
 TGeoBoolNode.cxx:610
 TGeoBoolNode.cxx:611
 TGeoBoolNode.cxx:612
 TGeoBoolNode.cxx:613
 TGeoBoolNode.cxx:614
 TGeoBoolNode.cxx:615
 TGeoBoolNode.cxx:616
 TGeoBoolNode.cxx:617
 TGeoBoolNode.cxx:618
 TGeoBoolNode.cxx:619
 TGeoBoolNode.cxx:620
 TGeoBoolNode.cxx:621
 TGeoBoolNode.cxx:622
 TGeoBoolNode.cxx:623
 TGeoBoolNode.cxx:624
 TGeoBoolNode.cxx:625
 TGeoBoolNode.cxx:626
 TGeoBoolNode.cxx:627
 TGeoBoolNode.cxx:628
 TGeoBoolNode.cxx:629
 TGeoBoolNode.cxx:630
 TGeoBoolNode.cxx:631
 TGeoBoolNode.cxx:632
 TGeoBoolNode.cxx:633
 TGeoBoolNode.cxx:634
 TGeoBoolNode.cxx:635
 TGeoBoolNode.cxx:636
 TGeoBoolNode.cxx:637
 TGeoBoolNode.cxx:638
 TGeoBoolNode.cxx:639
 TGeoBoolNode.cxx:640
 TGeoBoolNode.cxx:641
 TGeoBoolNode.cxx:642
 TGeoBoolNode.cxx:643
 TGeoBoolNode.cxx:644
 TGeoBoolNode.cxx:645
 TGeoBoolNode.cxx:646
 TGeoBoolNode.cxx:647
 TGeoBoolNode.cxx:648
 TGeoBoolNode.cxx:649
 TGeoBoolNode.cxx:650
 TGeoBoolNode.cxx:651
 TGeoBoolNode.cxx:652
 TGeoBoolNode.cxx:653
 TGeoBoolNode.cxx:654
 TGeoBoolNode.cxx:655
 TGeoBoolNode.cxx:656
 TGeoBoolNode.cxx:657
 TGeoBoolNode.cxx:658
 TGeoBoolNode.cxx:659
 TGeoBoolNode.cxx:660
 TGeoBoolNode.cxx:661
 TGeoBoolNode.cxx:662
 TGeoBoolNode.cxx:663
 TGeoBoolNode.cxx:664
 TGeoBoolNode.cxx:665
 TGeoBoolNode.cxx:666
 TGeoBoolNode.cxx:667
 TGeoBoolNode.cxx:668
 TGeoBoolNode.cxx:669
 TGeoBoolNode.cxx:670
 TGeoBoolNode.cxx:671
 TGeoBoolNode.cxx:672
 TGeoBoolNode.cxx:673
 TGeoBoolNode.cxx:674
 TGeoBoolNode.cxx:675
 TGeoBoolNode.cxx:676
 TGeoBoolNode.cxx:677
 TGeoBoolNode.cxx:678
 TGeoBoolNode.cxx:679
 TGeoBoolNode.cxx:680
 TGeoBoolNode.cxx:681
 TGeoBoolNode.cxx:682
 TGeoBoolNode.cxx:683
 TGeoBoolNode.cxx:684
 TGeoBoolNode.cxx:685
 TGeoBoolNode.cxx:686
 TGeoBoolNode.cxx:687
 TGeoBoolNode.cxx:688
 TGeoBoolNode.cxx:689
 TGeoBoolNode.cxx:690
 TGeoBoolNode.cxx:691
 TGeoBoolNode.cxx:692
 TGeoBoolNode.cxx:693
 TGeoBoolNode.cxx:694
 TGeoBoolNode.cxx:695
 TGeoBoolNode.cxx:696
 TGeoBoolNode.cxx:697
 TGeoBoolNode.cxx:698
 TGeoBoolNode.cxx:699
 TGeoBoolNode.cxx:700
 TGeoBoolNode.cxx:701
 TGeoBoolNode.cxx:702
 TGeoBoolNode.cxx:703
 TGeoBoolNode.cxx:704
 TGeoBoolNode.cxx:705
 TGeoBoolNode.cxx:706
 TGeoBoolNode.cxx:707
 TGeoBoolNode.cxx:708
 TGeoBoolNode.cxx:709
 TGeoBoolNode.cxx:710
 TGeoBoolNode.cxx:711
 TGeoBoolNode.cxx:712
 TGeoBoolNode.cxx:713
 TGeoBoolNode.cxx:714
 TGeoBoolNode.cxx:715
 TGeoBoolNode.cxx:716
 TGeoBoolNode.cxx:717
 TGeoBoolNode.cxx:718
 TGeoBoolNode.cxx:719
 TGeoBoolNode.cxx:720
 TGeoBoolNode.cxx:721
 TGeoBoolNode.cxx:722
 TGeoBoolNode.cxx:723
 TGeoBoolNode.cxx:724
 TGeoBoolNode.cxx:725
 TGeoBoolNode.cxx:726
 TGeoBoolNode.cxx:727
 TGeoBoolNode.cxx:728
 TGeoBoolNode.cxx:729
 TGeoBoolNode.cxx:730
 TGeoBoolNode.cxx:731
 TGeoBoolNode.cxx:732
 TGeoBoolNode.cxx:733
 TGeoBoolNode.cxx:734
 TGeoBoolNode.cxx:735
 TGeoBoolNode.cxx:736
 TGeoBoolNode.cxx:737
 TGeoBoolNode.cxx:738
 TGeoBoolNode.cxx:739
 TGeoBoolNode.cxx:740
 TGeoBoolNode.cxx:741
 TGeoBoolNode.cxx:742
 TGeoBoolNode.cxx:743
 TGeoBoolNode.cxx:744
 TGeoBoolNode.cxx:745
 TGeoBoolNode.cxx:746
 TGeoBoolNode.cxx:747
 TGeoBoolNode.cxx:748
 TGeoBoolNode.cxx:749
 TGeoBoolNode.cxx:750
 TGeoBoolNode.cxx:751
 TGeoBoolNode.cxx:752
 TGeoBoolNode.cxx:753
 TGeoBoolNode.cxx:754
 TGeoBoolNode.cxx:755
 TGeoBoolNode.cxx:756
 TGeoBoolNode.cxx:757
 TGeoBoolNode.cxx:758
 TGeoBoolNode.cxx:759
 TGeoBoolNode.cxx:760
 TGeoBoolNode.cxx:761
 TGeoBoolNode.cxx:762
 TGeoBoolNode.cxx:763
 TGeoBoolNode.cxx:764
 TGeoBoolNode.cxx:765
 TGeoBoolNode.cxx:766
 TGeoBoolNode.cxx:767
 TGeoBoolNode.cxx:768
 TGeoBoolNode.cxx:769
 TGeoBoolNode.cxx:770
 TGeoBoolNode.cxx:771
 TGeoBoolNode.cxx:772
 TGeoBoolNode.cxx:773
 TGeoBoolNode.cxx:774
 TGeoBoolNode.cxx:775
 TGeoBoolNode.cxx:776
 TGeoBoolNode.cxx:777
 TGeoBoolNode.cxx:778
 TGeoBoolNode.cxx:779
 TGeoBoolNode.cxx:780
 TGeoBoolNode.cxx:781
 TGeoBoolNode.cxx:782
 TGeoBoolNode.cxx:783
 TGeoBoolNode.cxx:784
 TGeoBoolNode.cxx:785
 TGeoBoolNode.cxx:786
 TGeoBoolNode.cxx:787
 TGeoBoolNode.cxx:788
 TGeoBoolNode.cxx:789
 TGeoBoolNode.cxx:790
 TGeoBoolNode.cxx:791
 TGeoBoolNode.cxx:792
 TGeoBoolNode.cxx:793
 TGeoBoolNode.cxx:794
 TGeoBoolNode.cxx:795
 TGeoBoolNode.cxx:796
 TGeoBoolNode.cxx:797
 TGeoBoolNode.cxx:798
 TGeoBoolNode.cxx:799
 TGeoBoolNode.cxx:800
 TGeoBoolNode.cxx:801
 TGeoBoolNode.cxx:802
 TGeoBoolNode.cxx:803
 TGeoBoolNode.cxx:804
 TGeoBoolNode.cxx:805
 TGeoBoolNode.cxx:806
 TGeoBoolNode.cxx:807
 TGeoBoolNode.cxx:808
 TGeoBoolNode.cxx:809
 TGeoBoolNode.cxx:810
 TGeoBoolNode.cxx:811
 TGeoBoolNode.cxx:812
 TGeoBoolNode.cxx:813
 TGeoBoolNode.cxx:814
 TGeoBoolNode.cxx:815
 TGeoBoolNode.cxx:816
 TGeoBoolNode.cxx:817
 TGeoBoolNode.cxx:818
 TGeoBoolNode.cxx:819
 TGeoBoolNode.cxx:820
 TGeoBoolNode.cxx:821
 TGeoBoolNode.cxx:822
 TGeoBoolNode.cxx:823
 TGeoBoolNode.cxx:824
 TGeoBoolNode.cxx:825
 TGeoBoolNode.cxx:826
 TGeoBoolNode.cxx:827
 TGeoBoolNode.cxx:828
 TGeoBoolNode.cxx:829
 TGeoBoolNode.cxx:830
 TGeoBoolNode.cxx:831
 TGeoBoolNode.cxx:832
 TGeoBoolNode.cxx:833
 TGeoBoolNode.cxx:834
 TGeoBoolNode.cxx:835
 TGeoBoolNode.cxx:836
 TGeoBoolNode.cxx:837
 TGeoBoolNode.cxx:838
 TGeoBoolNode.cxx:839
 TGeoBoolNode.cxx:840
 TGeoBoolNode.cxx:841
 TGeoBoolNode.cxx:842
 TGeoBoolNode.cxx:843
 TGeoBoolNode.cxx:844
 TGeoBoolNode.cxx:845
 TGeoBoolNode.cxx:846
 TGeoBoolNode.cxx:847
 TGeoBoolNode.cxx:848
 TGeoBoolNode.cxx:849
 TGeoBoolNode.cxx:850
 TGeoBoolNode.cxx:851
 TGeoBoolNode.cxx:852
 TGeoBoolNode.cxx:853
 TGeoBoolNode.cxx:854
 TGeoBoolNode.cxx:855
 TGeoBoolNode.cxx:856
 TGeoBoolNode.cxx:857
 TGeoBoolNode.cxx:858
 TGeoBoolNode.cxx:859
 TGeoBoolNode.cxx:860
 TGeoBoolNode.cxx:861
 TGeoBoolNode.cxx:862
 TGeoBoolNode.cxx:863
 TGeoBoolNode.cxx:864
 TGeoBoolNode.cxx:865
 TGeoBoolNode.cxx:866
 TGeoBoolNode.cxx:867
 TGeoBoolNode.cxx:868
 TGeoBoolNode.cxx:869
 TGeoBoolNode.cxx:870
 TGeoBoolNode.cxx:871
 TGeoBoolNode.cxx:872
 TGeoBoolNode.cxx:873
 TGeoBoolNode.cxx:874
 TGeoBoolNode.cxx:875
 TGeoBoolNode.cxx:876
 TGeoBoolNode.cxx:877
 TGeoBoolNode.cxx:878
 TGeoBoolNode.cxx:879
 TGeoBoolNode.cxx:880
 TGeoBoolNode.cxx:881
 TGeoBoolNode.cxx:882
 TGeoBoolNode.cxx:883
 TGeoBoolNode.cxx:884
 TGeoBoolNode.cxx:885
 TGeoBoolNode.cxx:886
 TGeoBoolNode.cxx:887
 TGeoBoolNode.cxx:888
 TGeoBoolNode.cxx:889
 TGeoBoolNode.cxx:890
 TGeoBoolNode.cxx:891
 TGeoBoolNode.cxx:892
 TGeoBoolNode.cxx:893
 TGeoBoolNode.cxx:894
 TGeoBoolNode.cxx:895
 TGeoBoolNode.cxx:896
 TGeoBoolNode.cxx:897
 TGeoBoolNode.cxx:898
 TGeoBoolNode.cxx:899
 TGeoBoolNode.cxx:900
 TGeoBoolNode.cxx:901
 TGeoBoolNode.cxx:902
 TGeoBoolNode.cxx:903
 TGeoBoolNode.cxx:904
 TGeoBoolNode.cxx:905
 TGeoBoolNode.cxx:906
 TGeoBoolNode.cxx:907
 TGeoBoolNode.cxx:908
 TGeoBoolNode.cxx:909
 TGeoBoolNode.cxx:910
 TGeoBoolNode.cxx:911
 TGeoBoolNode.cxx:912
 TGeoBoolNode.cxx:913
 TGeoBoolNode.cxx:914
 TGeoBoolNode.cxx:915
 TGeoBoolNode.cxx:916
 TGeoBoolNode.cxx:917
 TGeoBoolNode.cxx:918
 TGeoBoolNode.cxx:919
 TGeoBoolNode.cxx:920
 TGeoBoolNode.cxx:921
 TGeoBoolNode.cxx:922
 TGeoBoolNode.cxx:923
 TGeoBoolNode.cxx:924
 TGeoBoolNode.cxx:925
 TGeoBoolNode.cxx:926
 TGeoBoolNode.cxx:927
 TGeoBoolNode.cxx:928
 TGeoBoolNode.cxx:929
 TGeoBoolNode.cxx:930
 TGeoBoolNode.cxx:931
 TGeoBoolNode.cxx:932
 TGeoBoolNode.cxx:933
 TGeoBoolNode.cxx:934
 TGeoBoolNode.cxx:935
 TGeoBoolNode.cxx:936
 TGeoBoolNode.cxx:937
 TGeoBoolNode.cxx:938
 TGeoBoolNode.cxx:939
 TGeoBoolNode.cxx:940
 TGeoBoolNode.cxx:941
 TGeoBoolNode.cxx:942
 TGeoBoolNode.cxx:943
 TGeoBoolNode.cxx:944
 TGeoBoolNode.cxx:945
 TGeoBoolNode.cxx:946
 TGeoBoolNode.cxx:947
 TGeoBoolNode.cxx:948
 TGeoBoolNode.cxx:949
 TGeoBoolNode.cxx:950
 TGeoBoolNode.cxx:951
 TGeoBoolNode.cxx:952
 TGeoBoolNode.cxx:953
 TGeoBoolNode.cxx:954
 TGeoBoolNode.cxx:955
 TGeoBoolNode.cxx:956
 TGeoBoolNode.cxx:957
 TGeoBoolNode.cxx:958
 TGeoBoolNode.cxx:959
 TGeoBoolNode.cxx:960
 TGeoBoolNode.cxx:961
 TGeoBoolNode.cxx:962
 TGeoBoolNode.cxx:963
 TGeoBoolNode.cxx:964
 TGeoBoolNode.cxx:965
 TGeoBoolNode.cxx:966
 TGeoBoolNode.cxx:967
 TGeoBoolNode.cxx:968
 TGeoBoolNode.cxx:969
 TGeoBoolNode.cxx:970
 TGeoBoolNode.cxx:971
 TGeoBoolNode.cxx:972
 TGeoBoolNode.cxx:973
 TGeoBoolNode.cxx:974
 TGeoBoolNode.cxx:975
 TGeoBoolNode.cxx:976
 TGeoBoolNode.cxx:977
 TGeoBoolNode.cxx:978
 TGeoBoolNode.cxx:979
 TGeoBoolNode.cxx:980
 TGeoBoolNode.cxx:981
 TGeoBoolNode.cxx:982
 TGeoBoolNode.cxx:983
 TGeoBoolNode.cxx:984
 TGeoBoolNode.cxx:985
 TGeoBoolNode.cxx:986
 TGeoBoolNode.cxx:987
 TGeoBoolNode.cxx:988
 TGeoBoolNode.cxx:989
 TGeoBoolNode.cxx:990
 TGeoBoolNode.cxx:991
 TGeoBoolNode.cxx:992
 TGeoBoolNode.cxx:993
 TGeoBoolNode.cxx:994
 TGeoBoolNode.cxx:995
 TGeoBoolNode.cxx:996
 TGeoBoolNode.cxx:997
 TGeoBoolNode.cxx:998
 TGeoBoolNode.cxx:999
 TGeoBoolNode.cxx:1000
 TGeoBoolNode.cxx:1001
 TGeoBoolNode.cxx:1002
 TGeoBoolNode.cxx:1003
 TGeoBoolNode.cxx:1004
 TGeoBoolNode.cxx:1005
 TGeoBoolNode.cxx:1006
 TGeoBoolNode.cxx:1007
 TGeoBoolNode.cxx:1008
 TGeoBoolNode.cxx:1009
 TGeoBoolNode.cxx:1010
 TGeoBoolNode.cxx:1011
 TGeoBoolNode.cxx:1012
 TGeoBoolNode.cxx:1013
 TGeoBoolNode.cxx:1014
 TGeoBoolNode.cxx:1015
 TGeoBoolNode.cxx:1016
 TGeoBoolNode.cxx:1017
 TGeoBoolNode.cxx:1018
 TGeoBoolNode.cxx:1019
 TGeoBoolNode.cxx:1020
 TGeoBoolNode.cxx:1021
 TGeoBoolNode.cxx:1022
 TGeoBoolNode.cxx:1023
 TGeoBoolNode.cxx:1024
 TGeoBoolNode.cxx:1025
 TGeoBoolNode.cxx:1026
 TGeoBoolNode.cxx:1027
 TGeoBoolNode.cxx:1028
 TGeoBoolNode.cxx:1029
 TGeoBoolNode.cxx:1030
 TGeoBoolNode.cxx:1031
 TGeoBoolNode.cxx:1032
 TGeoBoolNode.cxx:1033
 TGeoBoolNode.cxx:1034
 TGeoBoolNode.cxx:1035
 TGeoBoolNode.cxx:1036
 TGeoBoolNode.cxx:1037
 TGeoBoolNode.cxx:1038
 TGeoBoolNode.cxx:1039
 TGeoBoolNode.cxx:1040
 TGeoBoolNode.cxx:1041
 TGeoBoolNode.cxx:1042
 TGeoBoolNode.cxx:1043
 TGeoBoolNode.cxx:1044
 TGeoBoolNode.cxx:1045
 TGeoBoolNode.cxx:1046
 TGeoBoolNode.cxx:1047
 TGeoBoolNode.cxx:1048
 TGeoBoolNode.cxx:1049
 TGeoBoolNode.cxx:1050
 TGeoBoolNode.cxx:1051
 TGeoBoolNode.cxx:1052
 TGeoBoolNode.cxx:1053
 TGeoBoolNode.cxx:1054
 TGeoBoolNode.cxx:1055
 TGeoBoolNode.cxx:1056
 TGeoBoolNode.cxx:1057
 TGeoBoolNode.cxx:1058
 TGeoBoolNode.cxx:1059
 TGeoBoolNode.cxx:1060
 TGeoBoolNode.cxx:1061
 TGeoBoolNode.cxx:1062
 TGeoBoolNode.cxx:1063
 TGeoBoolNode.cxx:1064
 TGeoBoolNode.cxx:1065
 TGeoBoolNode.cxx:1066
 TGeoBoolNode.cxx:1067
 TGeoBoolNode.cxx:1068
 TGeoBoolNode.cxx:1069
 TGeoBoolNode.cxx:1070
 TGeoBoolNode.cxx:1071
 TGeoBoolNode.cxx:1072
 TGeoBoolNode.cxx:1073
 TGeoBoolNode.cxx:1074
 TGeoBoolNode.cxx:1075
 TGeoBoolNode.cxx:1076
 TGeoBoolNode.cxx:1077
 TGeoBoolNode.cxx:1078
 TGeoBoolNode.cxx:1079
 TGeoBoolNode.cxx:1080
 TGeoBoolNode.cxx:1081
 TGeoBoolNode.cxx:1082
 TGeoBoolNode.cxx:1083
 TGeoBoolNode.cxx:1084
 TGeoBoolNode.cxx:1085
 TGeoBoolNode.cxx:1086
 TGeoBoolNode.cxx:1087
 TGeoBoolNode.cxx:1088
 TGeoBoolNode.cxx:1089
 TGeoBoolNode.cxx:1090
 TGeoBoolNode.cxx:1091
 TGeoBoolNode.cxx:1092
 TGeoBoolNode.cxx:1093
 TGeoBoolNode.cxx:1094
 TGeoBoolNode.cxx:1095
 TGeoBoolNode.cxx:1096
 TGeoBoolNode.cxx:1097
 TGeoBoolNode.cxx:1098
 TGeoBoolNode.cxx:1099
 TGeoBoolNode.cxx:1100
 TGeoBoolNode.cxx:1101
 TGeoBoolNode.cxx:1102
 TGeoBoolNode.cxx:1103
 TGeoBoolNode.cxx:1104
 TGeoBoolNode.cxx:1105
 TGeoBoolNode.cxx:1106
 TGeoBoolNode.cxx:1107
 TGeoBoolNode.cxx:1108
 TGeoBoolNode.cxx:1109
 TGeoBoolNode.cxx:1110
 TGeoBoolNode.cxx:1111
 TGeoBoolNode.cxx:1112
 TGeoBoolNode.cxx:1113
 TGeoBoolNode.cxx:1114
 TGeoBoolNode.cxx:1115
 TGeoBoolNode.cxx:1116
 TGeoBoolNode.cxx:1117
 TGeoBoolNode.cxx:1118
 TGeoBoolNode.cxx:1119
 TGeoBoolNode.cxx:1120
 TGeoBoolNode.cxx:1121
 TGeoBoolNode.cxx:1122
 TGeoBoolNode.cxx:1123
 TGeoBoolNode.cxx:1124
 TGeoBoolNode.cxx:1125
 TGeoBoolNode.cxx:1126
 TGeoBoolNode.cxx:1127
 TGeoBoolNode.cxx:1128
 TGeoBoolNode.cxx:1129
 TGeoBoolNode.cxx:1130
 TGeoBoolNode.cxx:1131
 TGeoBoolNode.cxx:1132
 TGeoBoolNode.cxx:1133
 TGeoBoolNode.cxx:1134
 TGeoBoolNode.cxx:1135
 TGeoBoolNode.cxx:1136
 TGeoBoolNode.cxx:1137
 TGeoBoolNode.cxx:1138
 TGeoBoolNode.cxx:1139
 TGeoBoolNode.cxx:1140
 TGeoBoolNode.cxx:1141
 TGeoBoolNode.cxx:1142
 TGeoBoolNode.cxx:1143
 TGeoBoolNode.cxx:1144
 TGeoBoolNode.cxx:1145
 TGeoBoolNode.cxx:1146
 TGeoBoolNode.cxx:1147
 TGeoBoolNode.cxx:1148
 TGeoBoolNode.cxx:1149
 TGeoBoolNode.cxx:1150
 TGeoBoolNode.cxx:1151
 TGeoBoolNode.cxx:1152
 TGeoBoolNode.cxx:1153
 TGeoBoolNode.cxx:1154
 TGeoBoolNode.cxx:1155
 TGeoBoolNode.cxx:1156
 TGeoBoolNode.cxx:1157
 TGeoBoolNode.cxx:1158
 TGeoBoolNode.cxx:1159
 TGeoBoolNode.cxx:1160
 TGeoBoolNode.cxx:1161
 TGeoBoolNode.cxx:1162
 TGeoBoolNode.cxx:1163
 TGeoBoolNode.cxx:1164
 TGeoBoolNode.cxx:1165
 TGeoBoolNode.cxx:1166
 TGeoBoolNode.cxx:1167
 TGeoBoolNode.cxx:1168
 TGeoBoolNode.cxx:1169
 TGeoBoolNode.cxx:1170
 TGeoBoolNode.cxx:1171
 TGeoBoolNode.cxx:1172
 TGeoBoolNode.cxx:1173
 TGeoBoolNode.cxx:1174
 TGeoBoolNode.cxx:1175
 TGeoBoolNode.cxx:1176
 TGeoBoolNode.cxx:1177
 TGeoBoolNode.cxx:1178
 TGeoBoolNode.cxx:1179
 TGeoBoolNode.cxx:1180
 TGeoBoolNode.cxx:1181
 TGeoBoolNode.cxx:1182
 TGeoBoolNode.cxx:1183
 TGeoBoolNode.cxx:1184
 TGeoBoolNode.cxx:1185
 TGeoBoolNode.cxx:1186
 TGeoBoolNode.cxx:1187
 TGeoBoolNode.cxx:1188
 TGeoBoolNode.cxx:1189
 TGeoBoolNode.cxx:1190
 TGeoBoolNode.cxx:1191
 TGeoBoolNode.cxx:1192
 TGeoBoolNode.cxx:1193
 TGeoBoolNode.cxx:1194
 TGeoBoolNode.cxx:1195
 TGeoBoolNode.cxx:1196
 TGeoBoolNode.cxx:1197
 TGeoBoolNode.cxx:1198
 TGeoBoolNode.cxx:1199
 TGeoBoolNode.cxx:1200
 TGeoBoolNode.cxx:1201
 TGeoBoolNode.cxx:1202
 TGeoBoolNode.cxx:1203
 TGeoBoolNode.cxx:1204
 TGeoBoolNode.cxx:1205
 TGeoBoolNode.cxx:1206
 TGeoBoolNode.cxx:1207
 TGeoBoolNode.cxx:1208
 TGeoBoolNode.cxx:1209
 TGeoBoolNode.cxx:1210
 TGeoBoolNode.cxx:1211
 TGeoBoolNode.cxx:1212
 TGeoBoolNode.cxx:1213
 TGeoBoolNode.cxx:1214
 TGeoBoolNode.cxx:1215
 TGeoBoolNode.cxx:1216
 TGeoBoolNode.cxx:1217
 TGeoBoolNode.cxx:1218
 TGeoBoolNode.cxx:1219
 TGeoBoolNode.cxx:1220
 TGeoBoolNode.cxx:1221
 TGeoBoolNode.cxx:1222
 TGeoBoolNode.cxx:1223
 TGeoBoolNode.cxx:1224
 TGeoBoolNode.cxx:1225
 TGeoBoolNode.cxx:1226
 TGeoBoolNode.cxx:1227
 TGeoBoolNode.cxx:1228
 TGeoBoolNode.cxx:1229
 TGeoBoolNode.cxx:1230
 TGeoBoolNode.cxx:1231
 TGeoBoolNode.cxx:1232
 TGeoBoolNode.cxx:1233
 TGeoBoolNode.cxx:1234