[985] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
[1347] | 26 | // $Id: G4VIntersectionLocator.cc,v 1.8 2010/07/13 15:59:42 gcosmo Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $ |
---|
[985] | 28 | // |
---|
| 29 | // Class G4VIntersectionLocator implementation |
---|
| 30 | // |
---|
| 31 | // 27.10.08 - John Apostolakis, Tatiana Nikitina. |
---|
| 32 | // --------------------------------------------------------------------------- |
---|
| 33 | |
---|
| 34 | #include <iomanip> |
---|
| 35 | |
---|
| 36 | #include "globals.hh" |
---|
| 37 | #include "G4ios.hh" |
---|
| 38 | #include "G4VIntersectionLocator.hh" |
---|
| 39 | #include "G4GeometryTolerance.hh" |
---|
| 40 | |
---|
| 41 | /////////////////////////////////////////////////////////////////////////// |
---|
| 42 | // |
---|
| 43 | // Constructor |
---|
| 44 | // |
---|
[1228] | 45 | G4VIntersectionLocator:: G4VIntersectionLocator(G4Navigator *theNavigator): |
---|
| 46 | fUseNormalCorrection(false), |
---|
| 47 | fiNavigator( theNavigator ), |
---|
| 48 | fiChordFinder( 0 ), // Not set - overridden at each step |
---|
| 49 | fiEpsilonStep( -1.0 ), // Out of range - overridden at each step |
---|
| 50 | fiDeltaIntersection( -1.0 ), // Out of range - overridden at each step |
---|
| 51 | fiUseSafety(false) // Default - overridden at each step |
---|
[985] | 52 | { |
---|
| 53 | kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); |
---|
| 54 | fVerboseLevel = 0; |
---|
| 55 | fHelpingNavigator = new G4Navigator(); |
---|
| 56 | } |
---|
| 57 | |
---|
| 58 | /////////////////////////////////////////////////////////////////////////// |
---|
| 59 | // |
---|
| 60 | // Destructor. |
---|
| 61 | // |
---|
| 62 | G4VIntersectionLocator::~G4VIntersectionLocator() |
---|
| 63 | { |
---|
| 64 | delete fHelpingNavigator; |
---|
| 65 | } |
---|
| 66 | |
---|
| 67 | /////////////////////////////////////////////////////////////////////////// |
---|
| 68 | // |
---|
| 69 | // Dumps status of propagator. |
---|
| 70 | // |
---|
| 71 | void |
---|
| 72 | G4VIntersectionLocator::printStatus( const G4FieldTrack& StartFT, |
---|
| 73 | const G4FieldTrack& CurrentFT, |
---|
| 74 | G4double requestStep, |
---|
| 75 | G4double safety, |
---|
| 76 | G4int stepNo) |
---|
| 77 | { |
---|
| 78 | const G4int verboseLevel= fVerboseLevel; |
---|
| 79 | const G4ThreeVector StartPosition = StartFT.GetPosition(); |
---|
| 80 | const G4ThreeVector StartUnitVelocity = StartFT.GetMomentumDir(); |
---|
| 81 | const G4ThreeVector CurrentPosition = CurrentFT.GetPosition(); |
---|
| 82 | const G4ThreeVector CurrentUnitVelocity = CurrentFT.GetMomentumDir(); |
---|
| 83 | |
---|
| 84 | G4double step_len = CurrentFT.GetCurveLength() - StartFT.GetCurveLength(); |
---|
[1347] | 85 | G4int oldprc; // cout/cerr precision settings |
---|
| 86 | |
---|
[985] | 87 | if( ((stepNo == 0) && (verboseLevel <3)) || (verboseLevel >= 3) ) |
---|
| 88 | { |
---|
[1347] | 89 | oldprc = G4cout.precision(4); |
---|
[985] | 90 | G4cout << std::setw( 6) << " " |
---|
| 91 | << std::setw( 25) << " Current Position and Direction" << " " |
---|
| 92 | << G4endl; |
---|
| 93 | G4cout << std::setw( 5) << "Step#" |
---|
| 94 | << std::setw(10) << " s " << " " |
---|
| 95 | << std::setw(10) << "X(mm)" << " " |
---|
| 96 | << std::setw(10) << "Y(mm)" << " " |
---|
| 97 | << std::setw(10) << "Z(mm)" << " " |
---|
| 98 | << std::setw( 7) << " N_x " << " " |
---|
| 99 | << std::setw( 7) << " N_y " << " " |
---|
| 100 | << std::setw( 7) << " N_z " << " " ; |
---|
[1347] | 101 | G4cout << std::setw( 7) << " Delta|N|" << " " |
---|
[985] | 102 | << std::setw( 9) << "StepLen" << " " |
---|
| 103 | << std::setw(12) << "StartSafety" << " " |
---|
| 104 | << std::setw( 9) << "PhsStep" << " "; |
---|
| 105 | G4cout << G4endl; |
---|
[1347] | 106 | G4cout.precision(oldprc); |
---|
[985] | 107 | } |
---|
| 108 | if((stepNo == 0) && (verboseLevel <=3)) |
---|
| 109 | { |
---|
| 110 | // Recurse to print the start values |
---|
| 111 | // |
---|
| 112 | printStatus( StartFT, StartFT, -1.0, safety, -1); |
---|
| 113 | } |
---|
| 114 | if( verboseLevel <= 3 ) |
---|
| 115 | { |
---|
| 116 | if( stepNo >= 0) |
---|
| 117 | { |
---|
| 118 | G4cout << std::setw( 4) << stepNo << " "; |
---|
| 119 | } |
---|
| 120 | else |
---|
| 121 | { |
---|
| 122 | G4cout << std::setw( 5) << "Start" ; |
---|
| 123 | } |
---|
[1347] | 124 | oldprc = G4cout.precision(8); |
---|
[985] | 125 | G4cout << std::setw(10) << CurrentFT.GetCurveLength() << " "; |
---|
| 126 | G4cout << std::setw(10) << CurrentPosition.x() << " " |
---|
| 127 | << std::setw(10) << CurrentPosition.y() << " " |
---|
| 128 | << std::setw(10) << CurrentPosition.z() << " "; |
---|
| 129 | G4cout.precision(4); |
---|
| 130 | G4cout << std::setw( 7) << CurrentUnitVelocity.x() << " " |
---|
| 131 | << std::setw( 7) << CurrentUnitVelocity.y() << " " |
---|
| 132 | << std::setw( 7) << CurrentUnitVelocity.z() << " "; |
---|
[1347] | 133 | G4cout.precision(3); |
---|
| 134 | G4cout << std::setw( 7) |
---|
| 135 | << CurrentFT.GetMomentum().mag()- StartFT.GetMomentum().mag() |
---|
| 136 | << " "; |
---|
| 137 | G4cout << std::setw( 9) << step_len << " "; |
---|
| 138 | G4cout << std::setw(12) << safety << " "; |
---|
| 139 | if( requestStep != -1.0 ) |
---|
| 140 | { |
---|
| 141 | G4cout << std::setw( 9) << requestStep << " "; |
---|
| 142 | } |
---|
| 143 | else |
---|
| 144 | { |
---|
| 145 | G4cout << std::setw( 9) << "Init/NotKnown" << " "; |
---|
| 146 | } |
---|
| 147 | G4cout << G4endl; |
---|
| 148 | G4cout.precision(oldprc); |
---|
| 149 | } |
---|
| 150 | else // if( verboseLevel > 3 ) |
---|
| 151 | { |
---|
| 152 | // Multi-line output |
---|
[985] | 153 | |
---|
[1347] | 154 | G4cout << "Step taken was " << step_len |
---|
| 155 | << " out of PhysicalStep= " << requestStep << G4endl; |
---|
| 156 | G4cout << "Final safety is: " << safety << G4endl; |
---|
| 157 | G4cout << "Chord length = " << (CurrentPosition-StartPosition).mag() |
---|
| 158 | << G4endl; |
---|
| 159 | G4cout << G4endl; |
---|
| 160 | } |
---|
[985] | 161 | } |
---|
| 162 | |
---|
| 163 | /////////////////////////////////////////////////////////////////////////// |
---|
| 164 | // |
---|
| 165 | // ReEstimateEndPoint. |
---|
| 166 | // |
---|
| 167 | G4FieldTrack G4VIntersectionLocator:: |
---|
| 168 | ReEstimateEndpoint( const G4FieldTrack &CurrentStateA, |
---|
| 169 | const G4FieldTrack &EstimatedEndStateB, |
---|
| 170 | G4double linearDistSq, |
---|
| 171 | G4double curveDist ) |
---|
| 172 | { |
---|
| 173 | G4FieldTrack newEndPoint( CurrentStateA ); |
---|
| 174 | G4MagInt_Driver* integrDriver= GetChordFinderFor()->GetIntegrationDriver(); |
---|
| 175 | |
---|
| 176 | G4FieldTrack retEndPoint( CurrentStateA ); |
---|
| 177 | G4bool goodAdvance; |
---|
| 178 | G4int itrial=0; |
---|
| 179 | const G4int no_trials= 20; |
---|
| 180 | |
---|
| 181 | G4double endCurveLen= EstimatedEndStateB.GetCurveLength(); |
---|
| 182 | do |
---|
| 183 | { |
---|
| 184 | G4double currentCurveLen= newEndPoint.GetCurveLength(); |
---|
| 185 | G4double advanceLength= endCurveLen - currentCurveLen ; |
---|
| 186 | if (std::abs(advanceLength)<kCarTolerance) |
---|
| 187 | { |
---|
[1228] | 188 | goodAdvance=true; |
---|
[985] | 189 | } |
---|
[1228] | 190 | else{ |
---|
[985] | 191 | goodAdvance= |
---|
| 192 | integrDriver->AccurateAdvance(newEndPoint, advanceLength, |
---|
| 193 | GetEpsilonStepFor()); |
---|
[1228] | 194 | } |
---|
[985] | 195 | } |
---|
| 196 | while( !goodAdvance && (++itrial < no_trials) ); |
---|
| 197 | |
---|
| 198 | if( goodAdvance ) |
---|
| 199 | { |
---|
| 200 | retEndPoint= newEndPoint; |
---|
| 201 | } |
---|
| 202 | else |
---|
| 203 | { |
---|
| 204 | retEndPoint= EstimatedEndStateB; // Could not improve without major work !! |
---|
| 205 | } |
---|
| 206 | |
---|
| 207 | // All the work is done |
---|
| 208 | // below are some diagnostics only -- before the return! |
---|
| 209 | // |
---|
| 210 | static const G4String MethodName("G4VIntersectionLocator::ReEstimateEndpoint"); |
---|
| 211 | |
---|
| 212 | #ifdef G4VERBOSE |
---|
| 213 | G4int latest_good_trials=0; |
---|
| 214 | if( itrial > 1) |
---|
| 215 | { |
---|
| 216 | if( fVerboseLevel > 0 ) |
---|
| 217 | { |
---|
| 218 | G4cout << MethodName << " called - goodAdv= " << goodAdvance |
---|
| 219 | << " trials = " << itrial |
---|
| 220 | << " previous good= " << latest_good_trials |
---|
| 221 | << G4endl; |
---|
| 222 | } |
---|
| 223 | latest_good_trials=0; |
---|
| 224 | } |
---|
| 225 | else |
---|
| 226 | { |
---|
| 227 | latest_good_trials++; |
---|
| 228 | } |
---|
| 229 | #endif |
---|
| 230 | |
---|
| 231 | #ifdef G4DEBUG_FIELD |
---|
| 232 | G4double lengthDone = newEndPoint.GetCurveLength() |
---|
| 233 | - CurrentStateA.GetCurveLength(); |
---|
| 234 | if( !goodAdvance ) |
---|
| 235 | { |
---|
| 236 | if( fVerboseLevel >= 3 ) |
---|
| 237 | { |
---|
| 238 | G4cout << MethodName << "> AccurateAdvance failed " ; |
---|
| 239 | G4cout << " in " << itrial << " integration trials/steps. " << G4endl; |
---|
| 240 | G4cout << " It went only " << lengthDone << " instead of " << curveDist |
---|
| 241 | << " -- a difference of " << curveDist - lengthDone << G4endl; |
---|
| 242 | G4cout << " ReEstimateEndpoint> Reset endPoint to original value!" |
---|
| 243 | << G4endl; |
---|
| 244 | } |
---|
| 245 | } |
---|
| 246 | |
---|
| 247 | static G4int noInaccuracyWarnings = 0; |
---|
| 248 | G4int maxNoWarnings = 10; |
---|
| 249 | if ( (noInaccuracyWarnings < maxNoWarnings ) |
---|
| 250 | || (fVerboseLevel > 1) ) |
---|
| 251 | { |
---|
| 252 | G4cerr << "G4PropagatorInField::LocateIntersectionPoint():" |
---|
| 253 | << G4endl |
---|
| 254 | << " Warning: Integration inaccuracy requires" |
---|
| 255 | << " an adjustment in the step's endpoint." << G4endl |
---|
| 256 | << " Two mid-points are further apart than their" |
---|
| 257 | << " curve length difference" << G4endl |
---|
| 258 | << " Dist = " << std::sqrt(linearDistSq) |
---|
| 259 | << " curve length = " << curveDist << G4endl; |
---|
| 260 | G4cerr << " Correction applied is " |
---|
| 261 | << (newEndPoint.GetPosition()-EstimatedEndStateB.GetPosition()).mag() |
---|
| 262 | << G4endl; |
---|
| 263 | } |
---|
| 264 | #else |
---|
| 265 | // Statistics on the RMS value of the corrections |
---|
| 266 | |
---|
| 267 | static G4int noCorrections=0; |
---|
| 268 | static G4double sumCorrectionsSq = 0; |
---|
| 269 | noCorrections++; |
---|
| 270 | if( goodAdvance ) |
---|
| 271 | { |
---|
| 272 | sumCorrectionsSq += (EstimatedEndStateB.GetPosition() - |
---|
| 273 | newEndPoint.GetPosition()).mag2(); |
---|
| 274 | } |
---|
| 275 | linearDistSq -= curveDist; // To use linearDistSq ... ! |
---|
| 276 | #endif |
---|
| 277 | |
---|
| 278 | return retEndPoint; |
---|
| 279 | } |
---|
| 280 | |
---|
| 281 | /////////////////////////////////////////////////////////////////////////// |
---|
| 282 | // |
---|
| 283 | // Method for finding SurfaceNormal of Intersecting Solid |
---|
| 284 | // |
---|
| 285 | G4ThreeVector G4VIntersectionLocator:: |
---|
| 286 | GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point, G4bool &validNormal) |
---|
| 287 | { |
---|
| 288 | G4ThreeVector Normal(G4ThreeVector(0,0,0)); |
---|
| 289 | G4VPhysicalVolume* located; |
---|
| 290 | |
---|
| 291 | validNormal = false; |
---|
| 292 | fHelpingNavigator->SetWorldVolume(GetNavigatorFor()->GetWorldVolume()); |
---|
| 293 | located = fHelpingNavigator->LocateGlobalPointAndSetup( CurrentE_Point ); |
---|
| 294 | G4TouchableHistoryHandle aTouchable = fHelpingNavigator |
---|
| 295 | ->CreateTouchableHistoryHandle(); |
---|
| 296 | G4ThreeVector localPosition = aTouchable->GetHistory() |
---|
| 297 | ->GetTopTransform().TransformPoint(CurrentE_Point); |
---|
| 298 | |
---|
| 299 | if( located != 0) |
---|
| 300 | { |
---|
| 301 | if (located->GetLogicalVolume() |
---|
| 302 | ->GetSolid()->Inside(localPosition)==kSurface) |
---|
| 303 | { |
---|
| 304 | Normal = located->GetLogicalVolume() |
---|
| 305 | ->GetSolid()->SurfaceNormal(localPosition); |
---|
| 306 | validNormal = true; |
---|
| 307 | } |
---|
| 308 | } |
---|
| 309 | |
---|
| 310 | return Normal; |
---|
| 311 | } |
---|
| 312 | |
---|
| 313 | /////////////////////////////////////////////////////////////////////////// |
---|
| 314 | // |
---|
| 315 | // Adjustment of Found Intersection |
---|
| 316 | // |
---|
| 317 | G4bool G4VIntersectionLocator:: |
---|
| 318 | AdjustmentOfFoundIntersection( const G4ThreeVector &CurrentA_Point, |
---|
| 319 | const G4ThreeVector &CurrentE_Point, |
---|
| 320 | const G4ThreeVector &CurrentF_Point, |
---|
| 321 | const G4ThreeVector &MomentumDir, |
---|
| 322 | const G4bool IntersectAF, |
---|
| 323 | G4ThreeVector &IntersectionPoint, // I/O |
---|
| 324 | G4double &NewSafety, // I/O |
---|
| 325 | G4double &fPreviousSafety, // I/O |
---|
| 326 | G4ThreeVector &fPreviousSftOrigin )// I/O |
---|
| 327 | { |
---|
| 328 | G4double dist,lambda; |
---|
| 329 | G4ThreeVector Normal, NewPoint, Point_G; |
---|
| 330 | G4bool goodAdjust=false, Intersects_FP=false, validNormal=false; |
---|
| 331 | |
---|
| 332 | // Get SurfaceNormal of Intersecting Solid |
---|
| 333 | // |
---|
| 334 | Normal=GetLocalSurfaceNormal(CurrentE_Point,validNormal); |
---|
| 335 | if(!validNormal) { return false; } |
---|
| 336 | |
---|
| 337 | // Intersection between Line and Plane |
---|
| 338 | // |
---|
| 339 | G4double n_d_m = Normal.dot(MomentumDir); |
---|
[1058] | 340 | if ( std::abs(n_d_m)>kCarTolerance ) |
---|
[985] | 341 | { |
---|
| 342 | if ( fVerboseLevel>1 ) |
---|
| 343 | { |
---|
| 344 | G4cerr << "WARNING - " |
---|
| 345 | << "G4VIntersectionLocator::AdjustementOfFoundIntersection()" |
---|
| 346 | << G4endl |
---|
| 347 | << " No intersection. Parallels lines!" << G4endl; |
---|
| 348 | return false; |
---|
| 349 | } |
---|
| 350 | lambda =- Normal.dot(CurrentF_Point-CurrentE_Point)/n_d_m; |
---|
| 351 | |
---|
| 352 | // New candidate for Intersection |
---|
| 353 | // |
---|
| 354 | NewPoint = CurrentF_Point+lambda*MomentumDir; |
---|
| 355 | |
---|
| 356 | // Distance from CurrentF to Calculated Intersection |
---|
| 357 | // |
---|
| 358 | dist = std::abs(lambda); |
---|
| 359 | |
---|
| 360 | if ( dist<kCarTolerance*0.001 ) { return false; } |
---|
| 361 | |
---|
| 362 | // Calculation of new intersection point on the path. |
---|
| 363 | // |
---|
| 364 | if ( IntersectAF ) // First part intersects |
---|
| 365 | { |
---|
| 366 | G4double stepLengthFP; |
---|
| 367 | G4ThreeVector Point_P = CurrentA_Point; |
---|
| 368 | GetNavigatorFor()->LocateGlobalPointWithinVolume(Point_P); |
---|
| 369 | Intersects_FP = IntersectChord( Point_P, NewPoint, NewSafety, |
---|
| 370 | fPreviousSafety, fPreviousSftOrigin, |
---|
| 371 | stepLengthFP, Point_G ); |
---|
| 372 | |
---|
| 373 | } |
---|
| 374 | else // Second part intersects |
---|
| 375 | { |
---|
| 376 | G4double stepLengthFP; |
---|
| 377 | GetNavigatorFor()->LocateGlobalPointWithinVolume(CurrentF_Point ); |
---|
| 378 | Intersects_FP = IntersectChord( CurrentF_Point, NewPoint, NewSafety, |
---|
| 379 | fPreviousSafety, fPreviousSftOrigin, |
---|
| 380 | stepLengthFP, Point_G ); |
---|
| 381 | } |
---|
| 382 | if ( Intersects_FP ) |
---|
| 383 | { |
---|
| 384 | goodAdjust = true; |
---|
| 385 | IntersectionPoint = Point_G; |
---|
| 386 | } |
---|
| 387 | } |
---|
| 388 | |
---|
| 389 | return goodAdjust; |
---|
| 390 | } |
---|