Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:54

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 G4PathFinder 
0027 //
0028 // Class description:
0029 // 
0030 // This class directs the lock-stepped propagation of a track in the 
0031 // 'mass' and other parallel geometries.  It ensures that tracking 
0032 // in a magnetic field sees these parallel geometries at each trial step, 
0033 // and that the earliest boundary limits the step.
0034 // 
0035 // For the movement in field, it relies on the class G4PropagatorInField
0036 
0037 // History:
0038 // -------
0039 //  7.10.05 John Apostolakis,  Draft design 
0040 // 26.04.06 John Apostolakis,  Revised design and first implementation 
0041 // ---------------------------------------------------------------------------
0042 #ifndef G4PATHFINDER_HH 
0043 #define G4PATHFINDER_HH  1
0044 
0045 #include <vector>
0046 #include "G4Types.hh"
0047 
0048 #include "G4FieldTrack.hh"
0049 
0050 class G4TransportationManager; 
0051 class G4Navigator;
0052 
0053 #include "G4MultiNavigator.hh"
0054 #include "G4TouchableHandle.hh"
0055 
0056 class G4PropagatorInField;
0057 
0058 class G4PathFinder
0059 {
0060 
0061  public:  // with description
0062 
0063    static G4PathFinder* GetInstance();
0064      // Retrieve singleton instance and create it if not existing.
0065 
0066    static G4PathFinder* GetInstanceIfExist();
0067      // Retrieve singleton instance pointer.
0068 
0069    G4double ComputeStep( const G4FieldTrack& pFieldTrack,
0070                          G4double  pCurrentProposedStepLength,
0071                          G4int     navigatorId, // Identifies the geometry
0072                          G4int     stepNo,      // See next step/check 
0073                          G4double& pNewSafety,  // Only for this geometry
0074                          ELimited& limitedStep,      
0075                          G4FieldTrack& EndState, 
0076                          G4VPhysicalVolume* currentVolume );
0077      // Compute the next geometric Step -- curved or linear
0078      // If it is called with a larger 'stepNo' it will execute a new step;
0079      // if 'stepNo' is same as last call, then the results for 
0080      // the geometry with Id. number 'navigatorId' will be returned. 
0081 
0082    void Locate( const G4ThreeVector& position, 
0083                 const G4ThreeVector& direction,
0084                       G4bool relativeSearch = true); 
0085      // Make primary relocation of global point in all navigators,
0086      // and update them.
0087 
0088    void ReLocate( const G4ThreeVector& position ); 
0089      // Make secondary relocation of global point (within safety only) 
0090      // in all navigators, and update them.
0091 
0092    void PrepareNewTrack( const G4ThreeVector& position,
0093                          const G4ThreeVector& direction,
0094                                G4VPhysicalVolume* massStartVol = nullptr); 
0095      // Check and cache set of active navigators.
0096 
0097    void EndTrack();
0098      // Signal end of tracking of current track.  
0099      //   Reset internal state
0100      //   Inform TransportationManager to use 'ordinary' Navigator
0101 
0102     G4TouchableHandle CreateTouchableHandle( G4int navId ) const;
0103     inline G4VPhysicalVolume* GetLocatedVolume( G4int navId ) const;
0104 
0105    // -----------------------------------------------------------------
0106   
0107    inline G4bool IsParticleLooping() const;
0108 
0109    inline G4double GetCurrentSafety() const;
0110      // Minimum value of safety after last ComputeStep
0111    inline G4double GetMinimumStep() const;      
0112      // Get the minimum step size from the last ComputeStep call
0113      //   - in case full step is taken, this is kInfinity
0114    inline unsigned int  GetNumberGeometriesLimitingStep() const; 
0115 
0116    G4double ComputeSafety( const G4ThreeVector& globalPoint); 
0117      // Recompute safety for the relevant point the endpoint of the last step!!
0118      // Maintain vector of individual safety values (for next method)
0119 
0120    G4double ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint );
0121      // Obtain safety for navigator/geometry navId for last point 'computed'
0122      //   --> last point for which ComputeSafety was called
0123      //   Returns the point (center) for which this safety is valid
0124 
0125    void EnableParallelNavigation( G4bool enableChoice = true ); 
0126      // Must call it to ensure that PathFinder is prepared,  
0127      // especially for curved tracks. If true it switches PropagatorInField
0128      // to use MultiNavigator. Must call it with false to undo (=PiF use
0129      // Navigator for tracking!)
0130 
0131    inline G4int  SetVerboseLevel(G4int lev = -1);
0132 
0133  public:  // with description
0134 
0135    inline G4int GetMaxLoopCount() const;
0136    inline void  SetMaxLoopCount( G4int new_max );
0137      // A maximum for the number of steps that a (looping) particle can take
0138 
0139  public:  // without description
0140 
0141    inline void MovePoint();
0142      // Signal that location will be moved -- internal use primarily
0143 
0144    // To provide best compatibility between Coupled and Old Transportation
0145    // the next two methods are provided:
0146 
0147    G4double LastPreSafety( G4int navId, G4ThreeVector& globalCenterPoint,
0148                            G4double& minSafety ); 
0149      // Obtain last safety needed in ComputeStep (for geometry navId)
0150      // --> last point at which ComputeStep recalculated safety
0151      //     Returns the point (center) for which this safety is valid
0152      //     and also the minimum safety over all navigators (i.e. full)
0153 
0154    void PushPostSafetyToPreSafety(); 
0155      // Tell PathFinder to copy PostStep Safety to PreSafety
0156      // (for use at next step)
0157 
0158    G4String& LimitedString( ELimited lim );
0159      // Convert ELimited to string
0160 
0161    ~G4PathFinder();
0162      // Destructor
0163 
0164  protected:  // without description
0165 
0166    G4double DoNextLinearStep( const G4FieldTrack& FieldTrack,
0167                                     G4double      proposedStepLength); 
0168 
0169    G4double DoNextCurvedStep( const G4FieldTrack& FieldTrack,
0170                                     G4double      proposedStepLength,
0171                                     G4VPhysicalVolume* pCurrentPhysVolume); 
0172 
0173   void WhichLimited();
0174   void PrintLimited();
0175     // Print key details out for debugging
0176 
0177   inline G4bool UseSafetyForOptimization( G4bool );
0178     // Whether use safety to discard unneccesary calls to navigator
0179 
0180   void ReportMove( const G4ThreeVector& OldV,
0181                    const G4ThreeVector& NewV,
0182                    const G4String& Quantity ) const;
0183     // Helper method to report movement (likely of initial point)
0184 
0185  protected:
0186 
0187    G4PathFinder();  //  Singleton 
0188 
0189    inline G4Navigator* GetNavigator(G4int n) const;
0190 
0191  private:
0192 
0193    // ----------------------------------------------------------------------
0194    //  DATA Members
0195    // ----------------------------------------------------------------------
0196 
0197    G4MultiNavigator* fpMultiNavigator; 
0198      // Object that enables G4PropagatorInField to see many geometries
0199 
0200    G4int fNoActiveNavigators = 0; 
0201    G4bool fNewTrack = false; // Flag a new track (ensure first step)
0202 
0203    static const G4int fMaxNav = 16;
0204 
0205    // Global state (retained during stepping for one track)
0206 
0207    G4Navigator*  fpNavigator[fMaxNav];
0208 
0209    // State changed in a step computation
0210 
0211    ELimited      fLimitedStep[fMaxNav];
0212    G4bool        fLimitTruth[fMaxNav];
0213    G4double      fCurrentStepSize[fMaxNav]; 
0214    G4int fNoGeometriesLimiting = 0;  // How many processes contribute to limit
0215 
0216    G4ThreeVector fPreSafetyLocation;
0217      // last initial position for which safety evaluated
0218    G4double      fPreSafetyMinValue = -1.0;
0219      // - corresponding value of full safety
0220    G4double      fPreSafetyValues[ fMaxNav ];
0221      // Safeties for the above point
0222 
0223    // This part of the state can be retained for several calls --> CARE
0224 
0225    G4ThreeVector fPreStepLocation;
0226      // point where last ComputeStep called
0227    G4double fMinSafety_PreStepPt = -1.0;
0228      // - corresponding value of full safety
0229    G4double fCurrentPreStepSafety[ fMaxNav ];
0230      // Safeties for the above point.
0231      // This changes at each step, so it can differ when steps
0232      // inside min-safety are made
0233 
0234    G4bool fPreStepCenterRenewed = false;
0235      // Whether PreSafety coincides with PreStep point 
0236 
0237    G4double fMinStep = -1.0;  // As reported by Navigators -- can be kInfinity
0238    G4double fTrueMinStep = -1.0;  // Corrected in case >= proposed
0239 
0240    // State after calling 'locate'
0241    //
0242    G4VPhysicalVolume* fLocatedVolume[fMaxNav];
0243    G4ThreeVector      fLastLocatedPosition; 
0244 
0245    // State after calling 'ComputeStep'
0246    // (others member variables will be affected)
0247    //
0248    G4FieldTrack fEndState;  // Point, velocity, ... at proposed step end
0249    G4bool fFieldExertedForce = false; // In current proposed step
0250 
0251    G4bool fRelocatedPoint = false; // Signals that point was or is being moved 
0252                                    // from the position of the last location or
0253                                    // the endpoint resulting from ComputeStep()
0254                                    // -- invalidates fEndState
0255 
0256    // State for 'ComputeSafety' and related methods
0257    //
0258    G4ThreeVector fSafetyLocation;
0259      // point where ComputeSafety is called
0260    G4double      fMinSafety_atSafLocation = -1.0;
0261      // - corresponding value of safety
0262    G4double      fNewSafetyComputed[ fMaxNav ];
0263      // Safeties for last ComputeSafety
0264 
0265    // State for Step numbers
0266    //
0267    G4int fLastStepNo = -1, fCurrentStepNo = -1; 
0268 
0269    G4int fVerboseLevel = 0;  // For debugging purposes
0270 
0271    G4TransportationManager* fpTransportManager; // Cache for frequent use
0272    G4PropagatorInField* fpFieldPropagator;
0273 
0274    G4double kCarTolerance;
0275 
0276    static G4ThreadLocal G4PathFinder* fpPathFinder;
0277 };
0278 
0279 // ********************************************************************
0280 // Inline methods.
0281 // ********************************************************************
0282 
0283 inline G4VPhysicalVolume* G4PathFinder::GetLocatedVolume( G4int navId ) const
0284 {
0285   G4VPhysicalVolume* vol = nullptr;  
0286   if( (navId < fMaxNav) && (navId >= 0) ) { vol= fLocatedVolume[navId]; }
0287   return vol; 
0288 }
0289 
0290 inline G4int G4PathFinder::SetVerboseLevel(G4int newLevel)
0291 {
0292   G4int old = fVerboseLevel;
0293   fVerboseLevel = newLevel;
0294   return old;
0295 }
0296 
0297 inline G4double G4PathFinder::GetMinimumStep() const
0298 { 
0299   return fMinStep; 
0300 } 
0301 
0302 inline unsigned int G4PathFinder::GetNumberGeometriesLimitingStep() const
0303 {
0304   unsigned int noGeometries = fNoGeometriesLimiting;
0305   return noGeometries;
0306 }
0307 
0308 inline G4double G4PathFinder::GetCurrentSafety() const
0309 {
0310   return fMinSafety_PreStepPt;
0311 }
0312 
0313 inline void G4PathFinder::MovePoint()
0314 {
0315   fRelocatedPoint = true;
0316 }
0317 
0318 inline G4Navigator* G4PathFinder::GetNavigator(G4int n) const
0319 { 
0320   if( (n>fNoActiveNavigators) || (n<0) ) { n=0; }
0321   return fpNavigator[n];
0322 }
0323 
0324 inline G4double
0325 G4PathFinder::ObtainSafety( G4int navId, G4ThreeVector& globalCenterPoint )
0326 {
0327   globalCenterPoint = fSafetyLocation; 
0328   return fNewSafetyComputed[ navId ];
0329 }
0330 
0331 inline G4double
0332 G4PathFinder::LastPreSafety( G4int navId, G4ThreeVector& globalCenterPoint, 
0333                              G4double& minSafety )
0334 {
0335   globalCenterPoint = fPreSafetyLocation;
0336   minSafety =         fPreSafetyMinValue;
0337   return  fPreSafetyValues[ navId ];
0338 }
0339 
0340 #endif