G4VIntersectionLocator.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 // $Id$
00027 //
00028 // Class G4VIntersectionLocator implementation
00029 //
00030 // 27.10.08 - John Apostolakis, Tatiana Nikitina.
00031 // ---------------------------------------------------------------------------
00032  
00033 #include <iomanip>
00034 #include <sstream>
00035 
00036 #include "globals.hh"
00037 #include "G4ios.hh"
00038 #include "G4SystemOfUnits.hh"
00039 #include "G4VIntersectionLocator.hh"
00040 #include "G4GeometryTolerance.hh"
00041 
00043 //
00044 // Constructor
00045 //
00046 G4VIntersectionLocator:: G4VIntersectionLocator(G4Navigator *theNavigator): 
00047   fUseNormalCorrection(false), 
00048   fiNavigator( theNavigator ),
00049   fiChordFinder( 0 ),            // Not set - overridden at each step
00050   fiEpsilonStep( -1.0 ),         // Out of range - overridden at each step
00051   fiDeltaIntersection( -1.0 ),   // Out of range - overridden at each step
00052   fiUseSafety(false),            // Default - overridden at each step
00053   fpTouchable(0)           
00054 {
00055   kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00056   fVerboseLevel = 0;
00057   fHelpingNavigator = new G4Navigator();
00058 }      
00059 
00061 //
00062 // Destructor.
00063 //
00064 G4VIntersectionLocator::~G4VIntersectionLocator()
00065 {
00066   delete fHelpingNavigator;
00067   delete fpTouchable;
00068 }
00069 
00071 //
00072 // Dumps status of propagator.
00073 //
00074 void
00075 G4VIntersectionLocator::printStatus( const G4FieldTrack&        StartFT,
00076                                      const G4FieldTrack&        CurrentFT, 
00077                                            G4double             requestStep, 
00078                                            G4double             safety,
00079                                            G4int                stepNo)
00080 {
00081   const G4int verboseLevel= fVerboseLevel;
00082   const G4ThreeVector StartPosition       = StartFT.GetPosition();
00083   const G4ThreeVector StartUnitVelocity   = StartFT.GetMomentumDir();
00084   const G4ThreeVector CurrentPosition     = CurrentFT.GetPosition();
00085   const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir();
00086 
00087   G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength();
00088   G4int oldprc;  // cout/cerr precision settings
00089 
00090   if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) )
00091   {
00092     oldprc = G4cout.precision(4);
00093     G4cout << std::setw( 6)  << " " 
00094            << std::setw( 25) << " Current Position  and  Direction" << " "
00095            << G4endl; 
00096     G4cout << std::setw( 5) << "Step#" 
00097            << std::setw(10) << "  s  " << " "
00098            << std::setw(10) << "X(mm)" << " "
00099            << std::setw(10) << "Y(mm)" << " "  
00100            << std::setw(10) << "Z(mm)" << " "
00101            << std::setw( 7) << " N_x " << " "
00102            << std::setw( 7) << " N_y " << " "
00103            << std::setw( 7) << " N_z " << " " ;
00104     G4cout << std::setw( 7) << " Delta|N|" << " "
00105            << std::setw( 9) << "StepLen" << " "  
00106            << std::setw(12) << "StartSafety" << " "  
00107            << std::setw( 9) << "PhsStep" << " ";  
00108     G4cout << G4endl;
00109     G4cout.precision(oldprc);
00110   }
00111   if((stepNo == 0) && (verboseLevel <=3))
00112   {
00113     // Recurse to print the start values
00114     //
00115     printStatus( StartFT, StartFT, -1.0, safety, -1);
00116   }
00117   if( verboseLevel <= 3 )
00118   {
00119     if( stepNo >= 0)
00120     {
00121        G4cout << std::setw( 4) << stepNo << " ";
00122     }
00123     else
00124     {
00125        G4cout << std::setw( 5) << "Start" ;
00126     }
00127     oldprc = G4cout.precision(8);
00128     G4cout << std::setw(10) << CurrentFT.GetCurveLength() << " "; 
00129     G4cout << std::setw(10) << CurrentPosition.x() << " "
00130            << std::setw(10) << CurrentPosition.y() << " "
00131            << std::setw(10) << CurrentPosition.z() << " ";
00132     G4cout.precision(4);
00133     G4cout << std::setw( 7) << CurrentUnitVelocity.x() << " "
00134            << std::setw( 7) << CurrentUnitVelocity.y() << " "
00135            << std::setw( 7) << CurrentUnitVelocity.z() << " ";
00136     G4cout.precision(3); 
00137     G4cout << std::setw( 7)
00138            << CurrentFT.GetMomentum().mag()- StartFT.GetMomentum().mag()
00139            << " "; 
00140     G4cout << std::setw( 9) << step_len << " "; 
00141     G4cout << std::setw(12) << safety << " ";
00142     if( requestStep != -1.0 )
00143     {
00144       G4cout << std::setw( 9) << requestStep << " ";
00145     }
00146     else
00147     {
00148       G4cout << std::setw( 9) << "Init/NotKnown" << " "; 
00149     }
00150     G4cout << G4endl;
00151     G4cout.precision(oldprc);
00152   }
00153   else // if( verboseLevel > 3 )
00154   {
00155     //  Multi-line output
00156        
00157     G4cout << "Step taken was " << step_len  
00158            << " out of PhysicalStep= " <<  requestStep << G4endl;
00159     G4cout << "Final safety is: " << safety << G4endl;
00160     G4cout << "Chord length = " << (CurrentPosition-StartPosition).mag()
00161            << G4endl;
00162     G4cout << G4endl; 
00163   }
00164 }
00165 
00167 //
00168 // ReEstimateEndPoint.
00169 //
00170 G4FieldTrack G4VIntersectionLocator::
00171 ReEstimateEndpoint( const G4FieldTrack& CurrentStateA,  
00172                     const G4FieldTrack& EstimatedEndStateB,
00173                           G4double      linearDistSq,
00174                           G4double      curveDist )
00175 {  
00176   G4FieldTrack newEndPoint( CurrentStateA );
00177   G4MagInt_Driver* integrDriver= GetChordFinderFor()->GetIntegrationDriver(); 
00178 
00179   G4FieldTrack retEndPoint( CurrentStateA );
00180   G4bool goodAdvance;
00181   G4int  itrial=0;
00182   const G4int no_trials= 20;
00183 
00184   G4double endCurveLen= EstimatedEndStateB.GetCurveLength();
00185   do
00186   {
00187      G4double currentCurveLen= newEndPoint.GetCurveLength();
00188      G4double advanceLength= endCurveLen - currentCurveLen ; 
00189      if (std::abs(advanceLength)<kCarTolerance)
00190      {
00191        goodAdvance=true;
00192      }
00193      else{
00194      goodAdvance= 
00195        integrDriver->AccurateAdvance(newEndPoint, advanceLength,
00196                                      GetEpsilonStepFor());
00197      }
00198   }
00199   while( !goodAdvance && (++itrial < no_trials) );
00200 
00201   if( goodAdvance )
00202   {
00203     retEndPoint= newEndPoint; 
00204   }
00205   else
00206   {
00207     retEndPoint= EstimatedEndStateB; // Could not improve without major work !!
00208   }
00209 
00210   //  All the work is done
00211   //  below are some diagnostics only -- before the return!
00212   // 
00213   static const G4String MethodName("G4VIntersectionLocator::ReEstimateEndpoint");
00214 
00215 #ifdef G4VERBOSE
00216   G4int  latest_good_trials=0;
00217   if( itrial > 1)
00218   {
00219     if( fVerboseLevel > 0 )
00220     {
00221       G4cout << MethodName << " called - goodAdv= " << goodAdvance
00222              << " trials = " << itrial
00223              << " previous good= " << latest_good_trials
00224              << G4endl;
00225     }
00226     latest_good_trials=0; 
00227   }
00228   else
00229   {   
00230     latest_good_trials++; 
00231   }
00232 #endif
00233 
00234 #ifdef G4DEBUG_FIELD
00235   G4double lengthDone = newEndPoint.GetCurveLength() 
00236                       - CurrentStateA.GetCurveLength(); 
00237   if( !goodAdvance )
00238   {
00239     if( fVerboseLevel >= 3 )
00240     {
00241       G4cout << MethodName << "> AccurateAdvance failed " ;
00242       G4cout << " in " << itrial << " integration trials/steps. " << G4endl;
00243       G4cout << " It went only " << lengthDone << " instead of " << curveDist
00244              << " -- a difference of " << curveDist - lengthDone  << G4endl;
00245       G4cout << " ReEstimateEndpoint> Reset endPoint to original value!"
00246              << G4endl;
00247     }
00248   }
00249 
00250   static G4int noInaccuracyWarnings = 0; 
00251   G4int maxNoWarnings = 10;
00252   if (  (noInaccuracyWarnings < maxNoWarnings ) 
00253        || (fVerboseLevel > 1) )
00254     {
00255       G4cerr << "G4PropagatorInField::LocateIntersectionPoint():"
00256              << G4endl
00257              << " Warning: Integration inaccuracy requires" 
00258              <<   " an adjustment in the step's endpoint."  << G4endl
00259              << "   Two mid-points are further apart than their"
00260              <<   " curve length difference"                << G4endl 
00261              << "   Dist = "       << std::sqrt(linearDistSq)
00262              << " curve length = " << curveDist             << G4endl; 
00263       G4cerr << " Correction applied is " 
00264              << (newEndPoint.GetPosition()-EstimatedEndStateB.GetPosition()).mag()
00265              << G4endl;
00266     }
00267 #else
00268   // Statistics on the RMS value of the corrections
00269 
00270   static G4int    noCorrections=0;
00271   static G4double sumCorrectionsSq = 0;
00272   noCorrections++; 
00273   if( goodAdvance )
00274   {
00275     sumCorrectionsSq += (EstimatedEndStateB.GetPosition() - 
00276                          newEndPoint.GetPosition()).mag2();
00277   }
00278   linearDistSq -= curveDist; // To use linearDistSq ... !
00279 #endif
00280 
00281   return retEndPoint;
00282 }
00283 
00285 //
00286 // Method for finding SurfaceNormal of Intersecting Solid 
00287 //
00288 G4ThreeVector G4VIntersectionLocator::
00289 GetLocalSurfaceNormal(const G4ThreeVector& CurrentE_Point, G4bool& validNormal)
00290 {
00291   G4ThreeVector Normal(G4ThreeVector(0.0,0.0,0.0));
00292   G4VPhysicalVolume* located;
00293 
00294   validNormal = false;
00295   fHelpingNavigator->SetWorldVolume(GetNavigatorFor()->GetWorldVolume());
00296   located = fHelpingNavigator->LocateGlobalPointAndSetup( CurrentE_Point );
00297 
00298   delete fpTouchable;
00299   fpTouchable = fHelpingNavigator->CreateTouchableHistory();
00300 
00301   // To check if we can use GetGlobalExitNormal() 
00302   //
00303   G4ThreeVector localPosition = fpTouchable->GetHistory()
00304                 ->GetTopTransform().TransformPoint(CurrentE_Point);
00305 
00306   // Issue: in the case of coincident surfaces, this version does not recognise 
00307   //        which side you are located onto (can return vector with wrong sign.)
00308   // TO-DO: use direction (of chord) to identify volume we will be "entering"
00309 
00310   if( located != 0)
00311   { 
00312     G4LogicalVolume* pLogical= located->GetLogicalVolume(); 
00313     G4VSolid*        pSolid; 
00314 
00315     if( (pLogical != 0) && ( (pSolid=pLogical->GetSolid()) !=0 )  )
00316     {
00317       // G4bool     goodPoint,    nearbyPoint;   
00318       // G4int   numGoodPoints,   numNearbyPoints;  // --> use for stats
00319       if ( ( pSolid->Inside(localPosition)==kSurface )
00320            || ( pSolid->DistanceToOut(localPosition) < 1000.0 * kCarTolerance )
00321          )
00322       {
00323         Normal = pSolid->SurfaceNormal(localPosition);
00324         validNormal = true;
00325 
00326 #ifdef G4DEBUG_FIELD
00327         if( std::fabs(Normal.mag2() - 1.0 ) > perMille) 
00328         {
00329           G4cerr << "PROBLEM in G4VIntersectionLocator::GetLocalSurfaceNormal."
00330                  << G4endl;
00331           G4cerr << "  Normal is not unit - mag=" << Normal.mag() << G4endl; 
00332           G4cerr << "  at trial local point " << CurrentE_Point << G4endl;
00333           G4cerr <<  "  Solid is " << *pSolid << G4endl;
00334         }
00335 #endif
00336       }
00337     }
00338   }
00339 
00340   return Normal;
00341 }
00342 
00344 //
00345 // Adjustment of Found Intersection
00346 //
00347 G4bool G4VIntersectionLocator::
00348 AdjustmentOfFoundIntersection( const G4ThreeVector& CurrentA_Point,
00349                                const G4ThreeVector& CurrentE_Point, 
00350                                const G4ThreeVector& CurrentF_Point,
00351                                const G4ThreeVector& MomentumDir,
00352                                const G4bool         IntersectAF,
00353                                      G4ThreeVector& IntersectionPoint,  // I/O
00354                                      G4double&      NewSafety,          // I/O 
00355                                      G4double&      fPreviousSafety,    // I/O
00356                                      G4ThreeVector& fPreviousSftOrigin )// I/O
00357 {     
00358   G4double dist,lambda;
00359   G4ThreeVector Normal, NewPoint, Point_G;
00360   G4bool goodAdjust=false, Intersects_FP=false, validNormal=false;
00361 
00362   // Get SurfaceNormal of Intersecting Solid
00363   //
00364   Normal = GetGlobalSurfaceNormal(CurrentE_Point,validNormal);
00365   if(!validNormal) { return false; }
00366 
00367   // Intersection between Line and Plane
00368   //
00369   G4double n_d_m = Normal.dot(MomentumDir);
00370   if ( std::abs(n_d_m)>kCarTolerance )
00371   {
00372 #ifdef G4VERBOSE
00373     if ( fVerboseLevel>1 )
00374     {
00375       G4cerr << "WARNING - "
00376              << "G4VIntersectionLocator::AdjustementOfFoundIntersection()"
00377              << G4endl
00378              << "        No intersection. Parallels lines!" << G4endl;
00379     }
00380 #endif
00381     lambda =- Normal.dot(CurrentF_Point-CurrentE_Point)/n_d_m;
00382 
00383     // New candidate for Intersection
00384     //
00385     NewPoint = CurrentF_Point+lambda*MomentumDir;
00386 
00387     // Distance from CurrentF to Calculated Intersection
00388     //
00389     dist = std::abs(lambda);
00390 
00391     if ( dist<kCarTolerance*0.001 )  { return false; }
00392 
00393     // Calculation of new intersection point on the path.
00394     //
00395     if ( IntersectAF )  //  First part intersects
00396     {
00397       G4double stepLengthFP; 
00398       G4ThreeVector Point_P = CurrentA_Point;
00399       GetNavigatorFor()->LocateGlobalPointWithinVolume(Point_P);
00400       Intersects_FP = IntersectChord( Point_P, NewPoint, NewSafety,
00401                                       fPreviousSafety, fPreviousSftOrigin,
00402                                       stepLengthFP, Point_G );
00403 
00404     }
00405     else   // Second part intersects
00406     {      
00407       G4double stepLengthFP; 
00408       GetNavigatorFor()->LocateGlobalPointWithinVolume(CurrentF_Point );
00409       Intersects_FP = IntersectChord( CurrentF_Point, NewPoint, NewSafety,
00410                                       fPreviousSafety, fPreviousSftOrigin,
00411                                       stepLengthFP, Point_G );
00412     }
00413     if ( Intersects_FP )
00414     {
00415       goodAdjust = true;
00416       IntersectionPoint = Point_G;              
00417     }
00418   }
00419 
00420   return goodAdjust;
00421 }
00422 
00423 G4ThreeVector
00424 G4VIntersectionLocator::GetSurfaceNormal(const G4ThreeVector& CurrentInt_Point,
00425                                                G4bool& validNormal) // const
00426 {
00427   G4ThreeVector NormalAtEntry; // ( -10. , -10., -10. ); 
00428 
00429   G4ThreeVector NormalAtEntryLast, NormalAtEntryGlobal, diffNormals;
00430   G4bool validNormalLast; 
00431 
00432   // Relies on a call to Navigator::ComputeStep in IntersectChord before
00433   // this call
00434   //
00435   NormalAtEntryLast = GetLastSurfaceNormal( CurrentInt_Point, validNormalLast );
00436     // May return valid=false in cases, including
00437     //  - if the candidate volume was not found (eg exiting world), or
00438     //  - a replica was involved -- determined the step size.
00439     // (This list is not complete.) 
00440 
00441 #ifdef G4DEBUG_FIELD
00442   if  ( validNormalLast
00443    && ( std::fabs(NormalAtEntryLast.mag2() - 1.0) > perThousand ) )
00444   {
00445     std::ostringstream message; 
00446     message << "G4VIntersectionLocator::GetSurfaceNormal -- identified problem."
00447             << G4endl;
00448     message << "PROBLEM: Normal is not unit - magnitude = "
00449             << NormalAtEntryLast.mag() << G4endl; 
00450     message << "   at trial intersection point " << CurrentInt_Point << G4endl;
00451     message << "   Obtained from Get *Last* Surface Normal." << G4endl; 
00452     G4Exception("G4VIntersectionLocator::GetGlobalSurfaceNormal()",
00453                 "GeomNav1002", JustWarning, message);
00454   }
00455 #endif
00456 
00457   if( validNormalLast ) 
00458   {
00459     NormalAtEntry=NormalAtEntryLast;  
00460   }
00461   validNormal  = validNormalLast;
00462 
00463   return NormalAtEntry;
00464 }
00465 
00466 G4ThreeVector G4VIntersectionLocator::
00467 GetGlobalSurfaceNormal(const G4ThreeVector& CurrentE_Point,
00468                              G4bool& validNormal)
00469 {
00470   G4ThreeVector     localNormal=
00471       GetLocalSurfaceNormal( CurrentE_Point, validNormal );
00472   G4AffineTransform localToGlobal=           //  Must use the same Navigator !!
00473       fHelpingNavigator->GetLocalToGlobalTransform();
00474   G4ThreeVector     globalNormal =
00475     localToGlobal.TransformAxis( localNormal );
00476 
00477 #ifdef G4DEBUG_FIELD
00478   if( validNormal && ( std::fabs(globalNormal.mag2() - 1.0) > perThousand ) ) 
00479   {
00480     std::ostringstream message; 
00481     message << "**************************************************************"
00482             << G4endl;
00483     message << " Bad Normal in G4VIntersectionLocator::GetGlobalSurfaceNormal "
00484             << G4endl;
00485     message << "  * Constituents: " << G4endl;
00486     message << "    Local  Normal= " << localNormal << G4endl;
00487     message << "    Transform: " << G4endl
00488             << "      Net Translation= " << localToGlobal.NetTranslation()
00489             << G4endl
00490             << "      Net Rotation   = " << localToGlobal.NetRotation()
00491             << G4endl;
00492     message << "  * Result: " << G4endl;
00493     message << "     Global Normal= " << localNormal << G4endl;
00494     message << "**************************************************************"
00495             << G4endl;
00496     G4Exception("G4VIntersectionLocator::GetGlobalSurfaceNormal()",
00497                 "GeomNav1002", JustWarning, message);
00498   }
00499 #endif
00500 
00501   return globalNormal;
00502 }
00503 
00504 G4ThreeVector 
00505 G4VIntersectionLocator::GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
00506                                               G4bool& normalIsValid) const
00507 {
00508   G4ThreeVector normalVec;
00509   G4bool        validNorm;
00510   normalVec = fiNavigator->GetGlobalExitNormal( intersectPoint, &validNorm ); 
00511   normalIsValid= validNorm;
00512 
00513   return normalVec;
00514 }
00515 
00516 void G4VIntersectionLocator::ReportTrialStep( G4int step_no, 
00517                                         const G4ThreeVector& ChordAB_v,
00518                                         const G4ThreeVector& ChordEF_v,
00519                                         const G4ThreeVector& NewMomentumDir,
00520                                         const G4ThreeVector& NormalAtEntry,
00521                                               G4bool validNormal )
00522 {
00523   G4double       ABchord_length  = ChordAB_v.mag(); 
00524   G4double       MomDir_dot_Norm = NewMomentumDir.dot( NormalAtEntry ) ;
00525   G4double       MomDir_dot_ABchord;
00526   MomDir_dot_ABchord= (1.0 / ABchord_length) * NewMomentumDir.dot( ChordAB_v );
00527 
00528   std::ostringstream  outStream; 
00529   outStream // G4cout 
00530     << std::setw(6)  << " Step# "
00531     << std::setw(17) << " |ChordEF|(mag)" << "  "
00532     << std::setw(18) << " uMomentum.Normal" << "  "
00533     << std::setw(18) << " uMomentum.ABdir " << "  " 
00534     << std::setw(16) << " AB-dist         " << " " 
00535     << " Chord Vector (EF) " 
00536     << G4endl;
00537   outStream.precision(7); 
00538   outStream  // G4cout 
00539     << " " << std::setw(5)  << step_no           
00540     << " " << std::setw(18) << ChordEF_v.mag() 
00541     << " " << std::setw(18) << MomDir_dot_Norm    
00542     << " " << std::setw(18) << MomDir_dot_ABchord 
00543     << " " << std::setw(12) << ABchord_length     
00544     << " " << ChordEF_v
00545     << G4endl;
00546   outStream  // G4cout
00547     << " MomentumDir= " << " " << NewMomentumDir 
00548     << " Normal at Entry E= " << NormalAtEntry
00549     << " AB chord =   " << ChordAB_v
00550     << G4endl;
00551   G4cout << outStream.str();  // ostr_verbose;
00552 
00553   if( ( std::fabs(NormalAtEntry.mag2() - 1.0) > perThousand ) ) 
00554   {
00555     G4cerr << " PROBLEM in G4VIntersectionLocator::ReportTrialStep " << G4endl
00556            << "         Normal is not unit - mag=" <<  NormalAtEntry.mag() 
00557            << "         ValidNormalAtE = " << validNormal
00558            << G4endl; 
00559   }
00560   return; 
00561 }

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