G4VoxelNavigation.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 //
00030 // class G4VoxelNavigation Implementation
00031 //
00032 // Author: P.Kent, 1996
00033 //
00034 // --------------------------------------------------------------------
00035 
00036 #include "G4VoxelNavigation.hh"
00037 #include "G4GeometryTolerance.hh"
00038 #include "G4VoxelSafety.hh"
00039 
00040 // ********************************************************************
00041 // Constructor
00042 // ********************************************************************
00043 //
00044 G4VoxelNavigation::G4VoxelNavigation()
00045   : fBList(), fVoxelDepth(-1),
00046     fVoxelAxisStack(kNavigatorVoxelStackMax,kXAxis),
00047     fVoxelNoSlicesStack(kNavigatorVoxelStackMax,0),
00048     fVoxelSliceWidthStack(kNavigatorVoxelStackMax,0.),
00049     fVoxelNodeNoStack(kNavigatorVoxelStackMax,0),
00050     fVoxelHeaderStack(kNavigatorVoxelStackMax,(G4SmartVoxelHeader*)0),
00051     fVoxelNode(0), fpVoxelSafety(0), fCheck(false), fBestSafety(false)
00052 {
00053   fLogger = new G4NavigationLogger("G4VoxelNavigation");
00054   fpVoxelSafety = new G4VoxelSafety (); 
00055 }
00056 
00057 // ********************************************************************
00058 // Destructor
00059 // ********************************************************************
00060 //
00061 G4VoxelNavigation::~G4VoxelNavigation()
00062 {
00063   delete fpVoxelSafety;
00064   delete fLogger;
00065 }
00066 
00067 // ********************************************************************
00068 // ComputeStep
00069 // ********************************************************************
00070 //
00071 G4double
00072 G4VoxelNavigation::ComputeStep( const G4ThreeVector& localPoint,
00073                                 const G4ThreeVector& localDirection,
00074                                 const G4double currentProposedStepLength,
00075                                       G4double& newSafety,
00076                                       G4NavigationHistory& history,
00077                                       G4bool& validExitNormal,
00078                                       G4ThreeVector& exitNormal,
00079                                       G4bool& exiting,
00080                                       G4bool& entering,
00081                                       G4VPhysicalVolume *(*pBlockedPhysical),
00082                                       G4int& blockedReplicaNo )
00083 {
00084   G4VPhysicalVolume *motherPhysical, *samplePhysical, *blockedExitedVol=0;
00085   G4LogicalVolume *motherLogical;
00086   G4VSolid *motherSolid;
00087   G4ThreeVector sampleDirection;
00088   G4double ourStep=currentProposedStepLength, motherSafety, ourSafety;
00089   G4int localNoDaughters, sampleNo;
00090 
00091   G4bool initialNode, noStep;
00092   G4SmartVoxelNode *curVoxelNode;
00093   G4int curNoVolumes, contentNo;
00094   G4double voxelSafety;
00095 
00096   motherPhysical = history.GetTopVolume();
00097   motherLogical = motherPhysical->GetLogicalVolume();
00098   motherSolid = motherLogical->GetSolid();
00099 
00100   //
00101   // Compute mother safety
00102   //
00103 
00104   motherSafety = motherSolid->DistanceToOut(localPoint);
00105   ourSafety = motherSafety;                 // Working isotropic safety
00106   
00107 #ifdef G4VERBOSE
00108   if ( fCheck )
00109   {
00110     fLogger->PreComputeStepLog (motherPhysical, motherSafety, localPoint);
00111   }
00112 #endif
00113 
00114   //
00115   // Compute daughter safeties & intersections
00116   //
00117 
00118   // Exiting normal optimisation
00119   //
00120   if ( exiting && validExitNormal )
00121   {
00122     if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
00123     {
00124       // Block exited daughter volume
00125       //
00126       blockedExitedVol = *pBlockedPhysical;
00127       ourSafety = 0;
00128     }
00129   }
00130   exiting = false;
00131   entering = false;
00132 
00133   localNoDaughters = motherLogical->GetNoDaughters();
00134 
00135   fBList.Enlarge(localNoDaughters);
00136   fBList.Reset();
00137 
00138   initialNode = true;
00139   noStep = true;
00140 
00141   while (noStep)
00142   {
00143     curVoxelNode = fVoxelNode;
00144     curNoVolumes = curVoxelNode->GetNoContained();
00145     for (contentNo=curNoVolumes-1; contentNo>=0; contentNo--)
00146     {
00147       sampleNo = curVoxelNode->GetVolume(contentNo);
00148       if ( !fBList.IsBlocked(sampleNo) )
00149       {
00150         fBList.BlockVolume(sampleNo);
00151         samplePhysical = motherLogical->GetDaughter(sampleNo);
00152         if ( samplePhysical!=blockedExitedVol )
00153         {
00154           G4AffineTransform sampleTf(samplePhysical->GetRotation(),
00155                                      samplePhysical->GetTranslation());
00156           sampleTf.Invert();
00157           const G4ThreeVector samplePoint =
00158                      sampleTf.TransformPoint(localPoint);
00159           const G4VSolid *sampleSolid     =
00160                      samplePhysical->GetLogicalVolume()->GetSolid();
00161           const G4double sampleSafety     =
00162                      sampleSolid->DistanceToIn(samplePoint);
00163 #ifdef G4VERBOSE
00164           if( fCheck )
00165           {
00166             fLogger->PrintDaughterLog(sampleSolid,samplePoint,sampleSafety,0);
00167           }
00168 #endif
00169           if ( sampleSafety<ourSafety )
00170           {
00171             ourSafety = sampleSafety;
00172           }
00173           if ( sampleSafety<=ourStep )
00174           {
00175             sampleDirection = sampleTf.TransformAxis(localDirection);
00176             G4double sampleStep =
00177                      sampleSolid->DistanceToIn(samplePoint, sampleDirection);
00178 #ifdef G4VERBOSE
00179             if( fCheck )
00180             {
00181               fLogger->PrintDaughterLog(sampleSolid, samplePoint,
00182                                         sampleSafety, sampleStep);
00183             }
00184 #endif
00185             if ( sampleStep<=ourStep )
00186             {
00187               ourStep = sampleStep;
00188               entering = true;
00189               exiting = false;
00190               *pBlockedPhysical = samplePhysical;
00191               blockedReplicaNo = -1;
00192 #ifdef G4VERBOSE
00193               // Check to see that the resulting point is indeed in/on volume.
00194               // This check could eventually be made only for successful
00195               // candidate.
00196 
00197               if ( fCheck )
00198               {
00199                 fLogger->AlongComputeStepLog (sampleSolid, samplePoint,
00200                   sampleDirection, localDirection, sampleSafety, sampleStep);
00201               }
00202 #endif
00203             }
00204           }
00205         }
00206       }
00207     }
00208     if (initialNode)
00209     {
00210       initialNode = false;
00211       voxelSafety = ComputeVoxelSafety(localPoint);
00212       if ( voxelSafety<ourSafety )
00213       {
00214         ourSafety = voxelSafety;
00215       }
00216       if ( currentProposedStepLength<ourSafety )
00217       {
00218         // Guaranteed physics limited
00219         //      
00220         noStep = false;
00221         entering = false;
00222         exiting = false;
00223         *pBlockedPhysical = 0;
00224         ourStep = kInfinity;
00225       }
00226       else
00227       {
00228         //
00229         // Compute mother intersection if required
00230         //
00231         if ( motherSafety<=ourStep )
00232         {
00233           G4double motherStep =
00234               motherSolid->DistanceToOut(localPoint,
00235                                          localDirection,
00236                                          true, &validExitNormal, &exitNormal);
00237 #ifdef G4VERBOSE
00238           if ( fCheck )
00239           {
00240             fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
00241                                         motherStep, motherSafety);
00242           }
00243 #endif
00244           if ( motherStep<=ourStep )
00245           {
00246             ourStep = motherStep;
00247             exiting = true;
00248             entering = false;
00249             if ( validExitNormal )
00250             {
00251               const G4RotationMatrix *rot = motherPhysical->GetRotation();
00252               if (rot)
00253               {
00254                 exitNormal *= rot->inverse();
00255               }
00256             }  
00257           }
00258           else
00259           {
00260             validExitNormal = false;
00261           }
00262         }
00263       }
00264       newSafety = ourSafety;
00265     }
00266     if (noStep)
00267     {
00268       noStep = LocateNextVoxel(localPoint, localDirection, ourStep);
00269     }
00270   }  // end -while (noStep)- loop
00271 
00272   return ourStep;
00273 }
00274 
00275 // ********************************************************************
00276 // ComputeVoxelSafety
00277 //
00278 // Computes safety from specified point to voxel boundaries
00279 // using already located point
00280 // o collected boundaries for most derived level
00281 // o adjacent boundaries for previous levels
00282 // ********************************************************************
00283 //
00284 G4double
00285 G4VoxelNavigation::ComputeVoxelSafety(const G4ThreeVector& localPoint) const
00286 {
00287   G4SmartVoxelHeader *curHeader;
00288   G4double voxelSafety, curNodeWidth;
00289   G4double curNodeOffset, minCurCommonDelta, maxCurCommonDelta;
00290   G4int minCurNodeNoDelta, maxCurNodeNoDelta;
00291   G4int localVoxelDepth, curNodeNo;
00292   EAxis curHeaderAxis;
00293 
00294   localVoxelDepth = fVoxelDepth;
00295 
00296   curHeader = fVoxelHeaderStack[localVoxelDepth];
00297   curHeaderAxis = fVoxelAxisStack[localVoxelDepth];
00298   curNodeNo = fVoxelNodeNoStack[localVoxelDepth];
00299   curNodeWidth = fVoxelSliceWidthStack[localVoxelDepth];
00300   
00301   // Compute linear intersection distance to boundaries of max/min
00302   // to collected nodes at current level
00303   //
00304   curNodeOffset = curNodeNo*curNodeWidth;
00305   maxCurNodeNoDelta = fVoxelNode->GetMaxEquivalentSliceNo()-curNodeNo;
00306   minCurNodeNoDelta = curNodeNo-fVoxelNode->GetMinEquivalentSliceNo();
00307   minCurCommonDelta = localPoint(curHeaderAxis)
00308                       - curHeader->GetMinExtent() - curNodeOffset;
00309   maxCurCommonDelta = curNodeWidth-minCurCommonDelta;
00310 
00311   if ( minCurNodeNoDelta<maxCurNodeNoDelta )
00312   {
00313     voxelSafety = minCurNodeNoDelta*curNodeWidth;
00314     voxelSafety += minCurCommonDelta;
00315   }
00316   else if (maxCurNodeNoDelta < minCurNodeNoDelta)
00317   {
00318     voxelSafety = maxCurNodeNoDelta*curNodeWidth;
00319     voxelSafety += maxCurCommonDelta;
00320   }
00321   else    // (maxCurNodeNoDelta == minCurNodeNoDelta)
00322   {
00323     voxelSafety = minCurNodeNoDelta*curNodeWidth;
00324     voxelSafety += std::min(minCurCommonDelta,maxCurCommonDelta);
00325   }
00326 
00327   // Compute isotropic safety to boundaries of previous levels
00328   // [NOT to collected boundaries]
00329   //
00330   while ( (localVoxelDepth>0) && (voxelSafety>0) )
00331   {
00332     localVoxelDepth--;
00333     curHeader = fVoxelHeaderStack[localVoxelDepth];
00334     curHeaderAxis = fVoxelAxisStack[localVoxelDepth];
00335     curNodeNo = fVoxelNodeNoStack[localVoxelDepth];
00336     curNodeWidth = fVoxelSliceWidthStack[localVoxelDepth];
00337     curNodeOffset = curNodeNo*curNodeWidth;
00338     minCurCommonDelta = localPoint(curHeaderAxis)
00339                         - curHeader->GetMinExtent() - curNodeOffset;
00340     maxCurCommonDelta = curNodeWidth-minCurCommonDelta;
00341     
00342     if ( minCurCommonDelta<voxelSafety )
00343     {
00344       voxelSafety = minCurCommonDelta;
00345     }
00346     if ( maxCurCommonDelta<voxelSafety )
00347     {
00348       voxelSafety = maxCurCommonDelta;
00349     }
00350   }
00351   if ( voxelSafety<0 )
00352   {
00353     voxelSafety = 0;
00354   }
00355 
00356   return voxelSafety;
00357 }
00358 
00359 // ********************************************************************
00360 // LocateNextVoxel
00361 //
00362 // Finds the next voxel from the current voxel and point
00363 // in the specified direction
00364 //
00365 // Returns false if all voxels considered
00366 //              [current Step ends inside same voxel or leaves all voxels]
00367 //         true  otherwise
00368 //              [the information on the next voxel is put into the set of
00369 //               fVoxel* variables & "stacks"] 
00370 // ********************************************************************
00371 // 
00372 G4bool
00373 G4VoxelNavigation::LocateNextVoxel(const G4ThreeVector& localPoint,
00374                                    const G4ThreeVector& localDirection,
00375                                    const G4double currentStep)
00376 {
00377   G4SmartVoxelHeader *workHeader=0, *newHeader=0;
00378   G4SmartVoxelProxy *newProxy=0;
00379   G4SmartVoxelNode *newVoxelNode=0;
00380   G4ThreeVector targetPoint, voxelPoint;
00381   G4double workNodeWidth, workMinExtent, workCoord;
00382   G4double minVal, maxVal, newDistance=0.;
00383   G4double newHeaderMin, newHeaderNodeWidth;
00384   G4int depth=0, newDepth=0, workNodeNo=0, newNodeNo=0, newHeaderNoSlices=0;
00385   EAxis workHeaderAxis, newHeaderAxis;
00386   G4bool isNewVoxel=false;
00387   
00388   G4double currentDistance = currentStep;
00389   static const G4double sigma = 0.5*G4GeometryTolerance::GetInstance()
00390                                     ->GetSurfaceTolerance();
00391 
00392   // Determine if end of Step within current voxel
00393   //
00394   for (depth=0; depth<fVoxelDepth; depth++)
00395   {
00396     targetPoint = localPoint+localDirection*currentDistance;
00397     newDistance = currentDistance;
00398     workHeader = fVoxelHeaderStack[depth];
00399     workHeaderAxis = fVoxelAxisStack[depth];
00400     workNodeNo = fVoxelNodeNoStack[depth];
00401     workNodeWidth = fVoxelSliceWidthStack[depth];
00402     workMinExtent = workHeader->GetMinExtent();
00403     workCoord = targetPoint(workHeaderAxis);
00404     minVal = workMinExtent+workNodeNo*workNodeWidth;
00405 
00406     if ( minVal<=workCoord+sigma )
00407     {
00408       maxVal = minVal+workNodeWidth;
00409       if ( maxVal<=workCoord-sigma )
00410       {
00411         // Must consider next voxel
00412         //
00413         newNodeNo = workNodeNo+1;
00414         newHeader = workHeader;
00415         newDistance = (maxVal-localPoint(workHeaderAxis))
00416                     / localDirection(workHeaderAxis);
00417         isNewVoxel = true;
00418         newDepth = depth;
00419       }
00420     }
00421     else
00422     {
00423       newNodeNo = workNodeNo-1;
00424       newHeader = workHeader;
00425       newDistance = (minVal-localPoint(workHeaderAxis))
00426                   / localDirection(workHeaderAxis);
00427       isNewVoxel = true;
00428       newDepth = depth;
00429     }
00430     currentDistance = newDistance;
00431   }
00432   targetPoint = localPoint+localDirection*currentDistance;
00433 
00434   // Check if end of Step within collected boundaries of current voxel
00435   //
00436   depth = fVoxelDepth;
00437   {
00438     workHeader = fVoxelHeaderStack[depth];
00439     workHeaderAxis = fVoxelAxisStack[depth];
00440     workNodeNo = fVoxelNodeNoStack[depth];
00441     workNodeWidth = fVoxelSliceWidthStack[depth];
00442     workMinExtent = workHeader->GetMinExtent();
00443     workCoord = targetPoint(workHeaderAxis);
00444     minVal = workMinExtent+fVoxelNode->GetMinEquivalentSliceNo()*workNodeWidth;
00445 
00446     if ( minVal<=workCoord+sigma )
00447     {
00448       maxVal = workMinExtent+(fVoxelNode->GetMaxEquivalentSliceNo()+1)
00449                             *workNodeWidth;
00450       if ( maxVal<=workCoord-sigma )
00451       {
00452         newNodeNo = fVoxelNode->GetMaxEquivalentSliceNo()+1;
00453         newHeader = workHeader;
00454         newDistance = (maxVal-localPoint(workHeaderAxis))
00455                     / localDirection(workHeaderAxis);
00456         isNewVoxel = true;
00457         newDepth = depth;
00458       }
00459     }
00460     else
00461     {
00462       newNodeNo = fVoxelNode->GetMinEquivalentSliceNo()-1;
00463       newHeader = workHeader;
00464       newDistance = (minVal-localPoint(workHeaderAxis))
00465                   / localDirection(workHeaderAxis);
00466       isNewVoxel = true;
00467       newDepth = depth;
00468     }
00469     currentDistance = newDistance;
00470   }
00471   if (isNewVoxel)
00472   {
00473     // Compute new voxel & adjust voxel stack
00474     //
00475     // newNodeNo=Candidate node no at 
00476     // newDepth =refinement depth of crossed voxel boundary
00477     // newHeader=Header for crossed voxel
00478     // newDistance=distance to crossed voxel boundary (along the track)
00479     //
00480     if ( (newNodeNo<0) || (newNodeNo>=newHeader->GetNoSlices()))
00481     {
00482       // Leaving mother volume
00483       //
00484       isNewVoxel = false;
00485     }
00486     else
00487     {
00488       // Compute intersection point on the least refined
00489       // voxel boundary that is hit
00490       //
00491       voxelPoint = localPoint+localDirection*newDistance;
00492       fVoxelNodeNoStack[newDepth] = newNodeNo;
00493       fVoxelDepth = newDepth;
00494       newVoxelNode = 0;
00495       while ( !newVoxelNode )
00496       {
00497         newProxy = newHeader->GetSlice(newNodeNo);
00498         if (newProxy->IsNode())
00499         {
00500           newVoxelNode = newProxy->GetNode();
00501         }
00502         else
00503         {
00504           fVoxelDepth++;
00505           newHeader = newProxy->GetHeader();
00506           newHeaderAxis = newHeader->GetAxis();
00507           newHeaderNoSlices = newHeader->GetNoSlices();
00508           newHeaderMin = newHeader->GetMinExtent();
00509           newHeaderNodeWidth = (newHeader->GetMaxExtent()-newHeaderMin)
00510                              / newHeaderNoSlices;
00511           newNodeNo = G4int( (voxelPoint(newHeaderAxis)-newHeaderMin)
00512                              / newHeaderNodeWidth );
00513           // Rounding protection
00514           //
00515           if ( newNodeNo<0 )
00516           {
00517             newNodeNo=0;
00518           }
00519           else if ( newNodeNo>=newHeaderNoSlices )
00520                {
00521                  newNodeNo = newHeaderNoSlices-1;
00522                }
00523           // Stack info for stepping
00524           //
00525           fVoxelAxisStack[fVoxelDepth] = newHeaderAxis;
00526           fVoxelNoSlicesStack[fVoxelDepth] = newHeaderNoSlices;
00527           fVoxelSliceWidthStack[fVoxelDepth] = newHeaderNodeWidth;
00528           fVoxelNodeNoStack[fVoxelDepth] = newNodeNo;
00529           fVoxelHeaderStack[fVoxelDepth] = newHeader;
00530         }
00531       }
00532       fVoxelNode = newVoxelNode;
00533     }
00534   }
00535   return isNewVoxel;        
00536 }
00537 
00538 // ********************************************************************
00539 // ComputeSafety
00540 //
00541 // Calculates the isotropic distance to the nearest boundary from the
00542 // specified point in the local coordinate system. 
00543 // The localpoint utilised must be within the current volume.
00544 // ********************************************************************
00545 //
00546 G4double
00547 G4VoxelNavigation::ComputeSafety(const G4ThreeVector& localPoint,
00548                                  const G4NavigationHistory& history,
00549                                  const G4double       maxLength)
00550 {
00551   G4VPhysicalVolume *motherPhysical, *samplePhysical;
00552   G4LogicalVolume *motherLogical;
00553   G4VSolid *motherSolid;
00554   G4double motherSafety, ourSafety;
00555   G4int sampleNo;
00556   G4SmartVoxelNode *curVoxelNode;
00557   G4int curNoVolumes, contentNo;
00558   G4double voxelSafety;
00559 
00560   motherPhysical = history.GetTopVolume();
00561   motherLogical = motherPhysical->GetLogicalVolume();
00562   motherSolid = motherLogical->GetSolid();
00563 
00564   if( fBestSafety )
00565   { 
00566     return fpVoxelSafety->ComputeSafety( localPoint,*motherPhysical,maxLength );
00567   }
00568 
00569   //
00570   // Compute mother safety
00571   //
00572 
00573   motherSafety = motherSolid->DistanceToOut(localPoint);
00574   ourSafety = motherSafety;                 // Working isotropic safety
00575 
00576   if( motherSafety == 0.0 )
00577   {
00578 #ifdef G4DEBUG_NAVIGATION
00579     // Check that point is inside mother volume
00580     EInside  insideMother= motherSolid->Inside(localPoint);
00581 
00582     if( insideMother == kOutside )
00583     {
00584       G4ExceptionDescription message;
00585       message << "Safety method called for location outside current Volume." << G4endl
00586          << "Location for safety is Outside this volume. " << G4endl
00587          << "The approximate distance to the solid "
00588          << "(safety from outside) is: "
00589          << motherSolid->DistanceToIn( localPoint ) << G4endl;
00590       message << "  Problem occurred with physical volume: "
00591          << " Name: " << motherPhysical->GetName()
00592          << " Copy No: " << motherPhysical->GetCopyNo() << G4endl
00593          << "    Local Point = " << localPoint << G4endl;
00594       message << "  Description of solid: " << G4endl
00595             << *motherSolid << G4endl;
00596       G4Exception("G4VoxelNavigation::ComputeSafety()", "GeomNav0003",
00597                   JustWarning,  // FatalException,
00598                   message);
00599     }
00600 
00601     // Following check is NOT for an issue - it is only for information
00602     //  It is allowed that a solid gives approximate safety - even zero.
00603     //
00604     if( insideMother == kInside ) // && fVerbose )
00605     {
00606       G4ExceptionDescription messageIn;
00607       
00608       messageIn << " Point is Inside, but safety is Zero ."  << G4endl;
00609       messageIn << " Inexact safety for volume " << motherPhysical->GetName() << G4endl
00610              << "  Solid: Name= " << motherSolid->GetName()
00611              << "   Type= " << motherSolid->GetEntityType() << G4endl;
00612       messageIn << "  Local point= " << localPoint << G4endl;
00613       messageIn << "  Solid parameters: " << G4endl << *motherSolid << G4endl;
00614       G4Exception("G4VoxelNavigation::ComputeSafety()", "GeomNav0003",
00615                   JustWarning, messageIn);
00616     }
00617 #endif
00618     // if( insideMother != kInside )
00619     return 0.0;
00620   }
00621    
00622 #ifdef G4VERBOSE
00623   if( fCheck )
00624   {
00625     fLogger->ComputeSafetyLog (motherSolid, localPoint, motherSafety, true);
00626   }
00627 #endif
00628   //
00629   // Compute daughter safeties
00630   //
00631   // Look only inside the current Voxel only (in the first version).
00632   //
00633   curVoxelNode = fVoxelNode;
00634   curNoVolumes = curVoxelNode->GetNoContained();
00635 
00636   for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
00637   {
00638     sampleNo = curVoxelNode->GetVolume(contentNo);
00639     samplePhysical = motherLogical->GetDaughter(sampleNo);
00640 
00641     G4AffineTransform sampleTf(samplePhysical->GetRotation(),
00642                                samplePhysical->GetTranslation());
00643     sampleTf.Invert();
00644     const G4ThreeVector samplePoint =
00645                           sampleTf.TransformPoint(localPoint);
00646     const G4VSolid *sampleSolid     =
00647                           samplePhysical->GetLogicalVolume()->GetSolid();
00648     G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
00649     if ( sampleSafety<ourSafety )
00650     {
00651       ourSafety = sampleSafety;
00652     }
00653 #ifdef G4VERBOSE
00654     if( fCheck )
00655     {
00656       fLogger->ComputeSafetyLog (sampleSolid,samplePoint,sampleSafety,false);
00657     }
00658 #endif
00659   }
00660   voxelSafety = ComputeVoxelSafety(localPoint);
00661   if ( voxelSafety<ourSafety )
00662   {
00663     ourSafety = voxelSafety;
00664   }
00665   return ourSafety;
00666 }
00667 
00668 // ********************************************************************
00669 // SetVerboseLevel
00670 // ********************************************************************
00671 //
00672 void  G4VoxelNavigation::SetVerboseLevel(G4int level)
00673 {
00674   if( fLogger )       fLogger->SetVerboseLevel(level);
00675   if( fpVoxelSafety)  fpVoxelSafety->SetVerboseLevel( level ); 
00676 }

Generated on Mon May 27 17:50:18 2013 for Geant4 by  doxygen 1.4.7