Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/Geant4/G4VIntersectionLocator.hh was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // Class G4VIntersectionLocator 
0027 //
0028 // class description:
0029 //
0030 // Base class for the calculation of the intersection point with a boundary 
0031 // when PropagationInField is used.
0032 // Gives possibility to choose the method of intersection; concrete locators
0033 // implemented are: G4SimpleLocator, G4MultiLevelLocator, G4BrentLocator.
0034 //
0035 // Key Method: EstimateIntersectionPoint()
0036 
0037 // 27.10.08 - John Apostolakis, Tatiana Nikitina: Design and implementation 
0038 // ---------------------------------------------------------------------------
0039 #ifndef G4VINTERSECTIONLOCATOR_HH
0040 #define G4VINTERSECTIONLOCATOR_HH
0041 
0042 #include "G4Types.hh" 
0043 #include "G4ThreeVector.hh"
0044 #include "G4FieldTrack.hh"
0045 
0046 #include "G4Navigator.hh"
0047 #include "G4ChordFinder.hh"
0048 
0049 class G4VIntersectionLocator
0050  {
0051    public:  // with description 
0052  
0053      G4VIntersectionLocator(G4Navigator *theNavigator);
0054        // Constructor
0055      virtual ~G4VIntersectionLocator();
0056        // Default destructor
0057      
0058      virtual G4bool EstimateIntersectionPoint( 
0059          const G4FieldTrack&       curveStartPointTangent,  // A
0060          const G4FieldTrack&       curveEndPointTangent,    // B
0061          const G4ThreeVector&      trialPoint,              // E
0062                G4FieldTrack&       intersectPointTangent,   // Output
0063                G4bool&             recalculatedEndPoint,    // Out
0064                G4double&           fPreviousSafety,         // In/Out
0065                G4ThreeVector&      fPreviousSftOrigin) = 0; // In/Out   
0066        // If such an intersection exists, this function calculates the
0067        // intersection point of the true path of the particle with the surface
0068        // of the current volume (or of one of its daughters). 
0069        // Should use lateral displacement as measure of convergence
0070        // NOTE: changes the safety!
0071 
0072      void printStatus( const G4FieldTrack& startFT,
0073                        const G4FieldTrack& currentFT, 
0074                              G4double      requestStep, 
0075                              G4double      safety,
0076                              G4int         stepNum);
0077        // Print Method, useful mostly for debugging
0078 
0079      inline G4bool IntersectChord( const G4ThreeVector&  StartPointA,
0080                                    const G4ThreeVector&  EndPointB,
0081                                    G4double&      NewSafety,
0082                                    G4double&      PreviousSafety,    // In/Out
0083                                    G4ThreeVector& PreviousSftOrigin, // In/Out
0084                                    G4double&      LinearStepLength,
0085                                    G4ThreeVector& IntersectionPoint,
0086                                    G4bool*        calledNavigator = nullptr );
0087        // Intersect the chord from StartPointA to EndPointB and return
0088        // whether an intersection occurred. NOTE: changes the Safety!
0089 
0090      inline void SetEpsilonStepFor( G4double EpsilonStep );
0091      inline void SetDeltaIntersectionFor( G4double deltaIntersection );
0092      inline void SetNavigatorFor( G4Navigator* fNavigator );
0093      inline void SetChordFinderFor(G4ChordFinder* fCFinder );
0094        // These parameters must be set at each step, in case they were changed
0095        // Note: This simple approach ensures that all scenarios are considered. 
0096        //       Future refinement may identify which are invariant during a 
0097        //       track, run or event
0098 
0099      inline void  SetVerboseFor(G4int fVerbose);
0100      inline G4int GetVerboseFor();
0101        // Controling verbosity enables checking of the locating of intersections
0102 
0103   public:  // without description
0104 
0105     // Additional inline Set/Get methods for parameters, dependent objects
0106 
0107     inline G4double       GetDeltaIntersectionFor();
0108     inline G4double       GetEpsilonStepFor();
0109     inline G4Navigator*   GetNavigatorFor();
0110     inline G4ChordFinder* GetChordFinderFor();
0111 
0112     inline void SetSafetyParametersFor(G4bool UseSafety );
0113 
0114     inline void AddAdjustementOfFoundIntersection(G4bool UseCorrection);
0115     inline G4bool GetAdjustementOfFoundIntersection();
0116       // Methods to be made Obsolete - replaced by methods below
0117     inline void AdjustIntersections(G4bool UseCorrection); 
0118     inline G4bool AreIntersectionsAdjusted(){ return fUseNormalCorrection; }  
0119       // Change adjustment flag  ( New Interface ) 
0120 
0121     static void printStatus( const G4FieldTrack& startFT,
0122                              const G4FieldTrack& currentFT, 
0123                                    G4double      requestStep, 
0124                                    G4double      safety,
0125                                    G4int         stepNum,
0126                                    std::ostream& oss,
0127                                    G4int         verboseLevel );
0128       // Print Method for any ostream - e.g. cerr -- and for G4Exception
0129 
0130     inline void SetCheckMode( G4bool value ) { fCheckMode = value; }
0131     inline G4bool GetCheckMode()             { return fCheckMode; }
0132 
0133   protected:  // with description
0134 
0135     G4FieldTrack ReEstimateEndpoint( const G4FieldTrack& CurrentStateA,  
0136                                      const G4FieldTrack& EstimtdEndStateB,
0137                                            G4double linearDistSq, // not used
0138                                            G4double curveDist );  // not used
0139       // Return new estimate for state after curveDist starting from
0140       // CurrentStateA, to replace EstimtdEndStateB, and report displacement
0141       // (if field is compiled verbose)
0142 
0143     G4bool CheckAndReEstimateEndpoint( const G4FieldTrack& CurrentStartA,  
0144                                        const G4FieldTrack& EstimatedEndB,
0145                                              G4FieldTrack& RevisedEndPoint,
0146                                              G4int &       errorCode);
0147       // Check whether EndB is too far from StartA to be reached 
0148       // and if, re-estimate new value for EndB (return in RevisedEndPoint)
0149       // Report error if  EndB is before StartA (in curve length)
0150       // In that case return errorCode = 2.
0151 
0152     G4ThreeVector GetSurfaceNormal(const G4ThreeVector& CurrentInt_Point,
0153                                          G4bool& validNormal);
0154       // Position *must* be the intersection point from last call
0155       // to G4Navigator's ComputeStep  (via IntersectChord )
0156       // Will try to use cached (last) value in Navigator for speed, 
0157       // if it was kept and valid.
0158       // Value returned is in global coordinates.
0159       // It does NOT guarantee to obtain Normal. This can happen eg if:
0160       //  - the "Intersection" Point is not on a surface, potentially due to
0161       //  - inaccuracies in the transformations used, or
0162       //  - issues with the Solid.
0163 
0164     G4ThreeVector GetGlobalSurfaceNormal(const G4ThreeVector& CurrentE_Point,
0165                                                G4bool& validNormal);
0166       // Return the SurfaceNormal of Intersecting Solid in global coordinates
0167       // Costlier then GetSurfaceNormal
0168 
0169     G4bool AdjustmentOfFoundIntersection(const G4ThreeVector& A,
0170                                          const G4ThreeVector& CurrentE_Point, 
0171                                          const G4ThreeVector& CurrentF_Point,
0172                                          const G4ThreeVector& MomentumDir,
0173                                          const G4bool         IntersectAF, 
0174                                                G4ThreeVector& IntersectionPoint,
0175                                                G4double&      NewSafety,
0176                                                G4double&      fPrevSafety,
0177                                                G4ThreeVector& fPrevSftOrigin );
0178       // Optional method for adjustment of located intersection point
0179       // using the surface-normal
0180   
0181     void ReportTrialStep( G4int step_no, 
0182                           const G4ThreeVector& ChordAB_v,
0183                           const G4ThreeVector& ChordEF_v,
0184                           const G4ThreeVector& NewMomentumDir,
0185                           const G4ThreeVector& NormalAtEntry,
0186                           G4bool validNormal   );
0187       // Print a three-line report on the current "sub-step", ie trial
0188       // intersection
0189 
0190     G4bool LocateGlobalPointWithinVolumeAndCheck( const G4ThreeVector& pos );
0191       // Locate point using navigator - updates state of Navigator.
0192       // By default, it assumes that the point is inside the current volume,
0193       // and returns true.
0194       // In check mode, checks whether the point is *inside* the volume.
0195       //   If it is inside, it returns true.
0196       //   If not, issues a warning and returns false.
0197 
0198     void LocateGlobalPointWithinVolumeCheckAndReport( const G4ThreeVector& pos,
0199                                             const G4String& CodeLocationInfo,
0200                                                   G4int     CheckMode );
0201       // As above, but report information about code location.
0202       // If CheckMode > 1, report extra information.
0203 
0204   protected:  // without description
0205 
0206     // Auxiliary methods -- to report issues
0207 
0208     void ReportReversedPoints( std::ostringstream& ossMsg,
0209                                const G4FieldTrack& StartPointVel, 
0210                                const G4FieldTrack& EndPointVel,
0211                                      G4double NewSafety, G4double epsStep,
0212                                const G4FieldTrack& CurrentA_PointVelocity,
0213                                const G4FieldTrack& CurrentB_PointVelocity,
0214                                const G4FieldTrack& SubStart_PointVelocity,
0215                                const G4ThreeVector& CurrentE_Point,
0216                                const G4FieldTrack& ApproxIntersecPointV,
0217                                G4int sbstp_no, G4int sbstp_no_p, G4int depth );
0218       // Build error messsage (in ossMsg) to report that point 'B' has
0219       // gone past 'A'
0220 
0221     void ReportProgress( std::ostream& oss,
0222                          const G4FieldTrack& StartPointVel, 
0223                          const G4FieldTrack& EndPointVel,
0224                                G4int         substep_no, 
0225                          const G4FieldTrack& A_PtVel,    // G4double safetyA
0226                          const G4FieldTrack& B_PtVel,  
0227                                G4double      safetyLast,
0228                                G4int         depth= -1 );
0229       // Report the current status / progress in finding the first intersection
0230 
0231      void ReportImmediateHit( const char*          MethodName, 
0232                               const G4ThreeVector& StartPosition, 
0233                               const G4ThreeVector& TrialPoint, 
0234                                     G4double       tolerance,
0235                                  unsigned long int numCalls );
0236       // Report case: trial point is 'close' to start, within tolerance
0237     
0238   private:  // no description
0239 
0240     G4ThreeVector GetLocalSurfaceNormal(const G4ThreeVector& CurrentE_Point,
0241                                               G4bool& validNormal);
0242       // Return the SurfaceNormal of Intersecting Solid  in local coordinates
0243 
0244     G4ThreeVector GetLastSurfaceNormal( const G4ThreeVector& intersectPoint,
0245                                               G4bool& validNormal) const; 
0246       // Position *must* be the intersection point from last call
0247       // to G4Navigator's ComputeStep  (via IntersectChord )
0248       // Temporary - will use the same method in the Navigator
0249 
0250   protected:
0251 
0252     G4double kCarTolerance;                  // Constant
0253 
0254     G4int    fVerboseLevel = 0;              // For debugging
0255     G4bool   fUseNormalCorrection = false;   // Configuration parameter
0256     G4bool   fCheckMode = false;
0257     G4bool   fiUseSafety = false;    // Whether to use safety for 'fast steps'
0258    
0259     G4Navigator* fiNavigator;
0260 
0261     G4ChordFinder* fiChordFinder = nullptr;  // Overridden at each step
0262     G4double fiEpsilonStep = -1.0;           // Overridden at each step
0263     G4double fiDeltaIntersection = -1.0;     // Overridden at each step
0264       // Parameters set at each physical step by calling method 
0265       // by G4PropagatorInField
0266 
0267     G4Navigator *fHelpingNavigator;
0268       // Helper for location
0269 
0270     G4TouchableHistory *fpTouchable = nullptr;
0271       // Touchable history hook
0272 };
0273 
0274 #include "G4VIntersectionLocator.icc"
0275 
0276 #endif