G4VIntersectionLocator.hh

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 G4VIntersectionLocator 
00031 //
00032 // class description:
00033 // 
00034 // Base class for the calculation of the intersection point with a boundary 
00035 // when PropagationInField is used.
00036 // Gives possibility to choose the method of intersection; concrete locators
00037 // implemented are: G4SimpleLocator, G4MultiLevelLocator, G4BrentLocator.
00038 //
00039 // Key Method: EstimateIntersectionPoint()
00040 
00041 // History:
00042 // -------
00043 // 27.10.08 - John Apostolakis, Tatiana Nikitina: Design and implementation 
00044 // ---------------------------------------------------------------------------
00045 
00046 #ifndef G4VINTERSECTIONLOCATOR_HH
00047 #define G4VINTERSECTIONLOCATOR_HH
00048 
00049 #include "G4Types.hh" 
00050 #include "G4ThreeVector.hh"
00051 #include "G4FieldTrack.hh"
00052 
00053 #include "G4Navigator.hh"
00054 #include "G4ChordFinder.hh"
00055 
00056 class G4VIntersectionLocator
00057  {
00058    public:  // with description 
00059  
00060      G4VIntersectionLocator(G4Navigator *theNavigator);
00061        // Constructor
00062      virtual ~G4VIntersectionLocator();
00063        // Default destructor
00064      
00065      virtual G4bool EstimateIntersectionPoint( 
00066          const  G4FieldTrack&       curveStartPointTangent,  // A
00067          const  G4FieldTrack&       curveEndPointTangent,    // B
00068          const  G4ThreeVector&      trialPoint,              // E
00069                 G4FieldTrack&       intersectPointTangent,   // Output
00070                 G4bool&             recalculatedEndPoint,    // Out
00071                 G4double&           fPreviousSafety,         // In/Out
00072                 G4ThreeVector&      fPreviousSftOrigin) = 0; // In/Out   
00073        // If such an intersection exists, this function calculates the
00074        // intersection point of the true path of the particle with the surface
00075        // of the current volume (or of one of its daughters). 
00076        // Should use lateral displacement as measure of convergence
00077        // NOTE: changes the safety!
00078 
00079      void printStatus( const G4FieldTrack& startFT,
00080                        const G4FieldTrack& currentFT, 
00081                              G4double      requestStep, 
00082                              G4double      safety,
00083                              G4int         step);
00084        // Print Method, useful mostly for debugging
00085 
00086      inline G4bool IntersectChord( const G4ThreeVector&  StartPointA,
00087                                    const G4ThreeVector&  EndPointB,
00088                                    G4double      &NewSafety,
00089                                    G4double      &PreviousSafety,    // In/Out
00090                                    G4ThreeVector &PreviousSftOrigin, // In/Out
00091                                    G4double      &LinearStepLength,
00092                                    G4ThreeVector &IntersectionPoint,
00093                                    G4bool        *calledNavigator=0 );
00094        // Intersect the chord from StartPointA to EndPointB and return
00095        // whether an intersection occurred. NOTE: changes the Safety!
00096 
00097      inline void    SetEpsilonStepFor( G4double EpsilonStep );
00098      inline void    SetDeltaIntersectionFor( G4double deltaIntersection );
00099      inline void    SetNavigatorFor( G4Navigator *fNavigator );
00100      inline void    SetChordFinderFor(G4ChordFinder *fCFinder );
00101        // These parameters must be set at each step, in case they were changed
00102 
00103        // Note: This simple approach ensures that all scenarios are considered. 
00104        //   [ Future refinement may identify which are invariant during a 
00105        //      track, run or event ]
00106 
00107     inline void     SetVerboseFor(G4int fVerbose);
00108     inline G4int    GetVerboseFor();
00109        // Controling verbosity enables checking of the locating of intersections
00110 
00111   public:  // without description
00112     // Additional inline Set/Get methods for parameters, dependent objects
00113 
00114     inline G4double       GetDeltaIntersectionFor();
00115     inline G4double       GetEpsilonStepFor();
00116     inline G4Navigator*   GetNavigatorFor();
00117     inline G4ChordFinder* GetChordFinderFor();
00118 
00119     inline void   SetSafetyParametersFor(G4bool UseSafety );
00120 
00121     inline void   AddAdjustementOfFoundIntersection(G4bool UseCorrection);
00122     inline G4bool GetAdjustementOfFoundIntersection();
00123       // Methods to be made Obsolete - replaced by methods below
00124     inline void   AdjustIntersections(G4bool UseCorrection); 
00125     inline G4bool AreIntersectionsAdjusted(){ return fUseNormalCorrection; }  
00126       // Change adjustment flag  ( New Interface ) 
00127 
00128   protected:  // with description
00129 
00130     G4FieldTrack ReEstimateEndpoint( const G4FieldTrack &CurrentStateA,  
00131                                      const G4FieldTrack &EstimtdEndStateB,
00132                                            G4double      linearDistSq,
00133                                            G4double      curveDist);
00134       // Return new estimate for state after curveDist starting from
00135       // CurrentStateA, to replace EstimtdEndStateB, and report displacement
00136       // (if field is compiled verbose)
00137 
00138     G4ThreeVector GetSurfaceNormal(const G4ThreeVector &CurrentInt_Point,
00139                                          G4bool        &validNormal); // const
00140       // Position *must* be the intersection point from last call
00141       // to G4Navigator's ComputeStep  (via IntersectChord )
00142       // Will try to use cached (last) value in Navigator for speed, 
00143       // if it was kept and valid.
00144       // Value returned is in global coordinates.
00145       // It does NOT guarantee to obtain Normal. This can happen eg if:
00146       //  - the "Intersection" Point is not on a surface, potentially due to
00147       //  - inaccuracies in the transformations used, or
00148       //  - issues with the Solid.
00149 
00150     G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector &CurrentE_Point,
00151                                                G4bool        &validNormal);
00152       // Return the SurfaceNormal of Intersecting Solid in global coordinates
00153       // Costlier then GetSurfaceNormal
00154 
00155     G4bool AdjustmentOfFoundIntersection(const G4ThreeVector &A,
00156                                          const G4ThreeVector &CurrentE_Point, 
00157                                          const G4ThreeVector &CurrentF_Point,
00158                                          const G4ThreeVector &MomentumDir,
00159                                          const G4bool         IntersectAF, 
00160                                                G4ThreeVector &IntersectionPoint,
00161                                                G4double      &NewSafety,
00162                                                G4double      &fPrevSafety,
00163                                                G4ThreeVector &fPrevSftOrigin );
00164       // Optional method for adjustment of located intersection point
00165       // using the surface-normal
00166   
00167     void ReportTrialStep( G4int step_no, 
00168                           const G4ThreeVector& ChordAB_v,
00169                           const G4ThreeVector& ChordEF_v,
00170                           const G4ThreeVector& NewMomentumDir,
00171                           const G4ThreeVector& NormalAtEntry,
00172                           G4bool validNormal   );
00173       // Print a three-line report on the current "sub-step", ie trial intersection
00174 
00175   private:  // no description
00176 
00177     G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector &CurrentE_Point,
00178                                               G4bool &validNormal);
00179       // Return the SurfaceNormal of Intersecting Solid  in local coordinates
00180 
00181     G4ThreeVector GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
00182                                         G4bool        &validNormal) const; 
00183       // Position *must* be the intersection point from last call
00184       // to G4Navigator's ComputeStep  (via IntersectChord )
00185       // Temporary - will use the same method in the Navigator
00186 
00187   protected:
00188 
00189     G4double kCarTolerance;         // Constant
00190 
00191     G4int    fVerboseLevel;          // For debugging
00192     G4bool   fUseNormalCorrection;   // Configuration parameter
00193 
00194     G4Navigator   *fiNavigator;
00195 
00196     G4ChordFinder *fiChordFinder;
00197     G4double       fiEpsilonStep;
00198     G4double       fiDeltaIntersection;
00199     G4bool         fiUseSafety;
00200       // Parameters set at each physical step by calling method 
00201       // by G4PropagatorInField
00202 
00203     G4Navigator *fHelpingNavigator;
00204       // Helper for location
00205 
00206     G4TouchableHistory *fpTouchable;
00207       // Touchable history hook
00208 };
00209 
00210 #include "G4VIntersectionLocator.icc"
00211 
00212 #endif

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