Back to home page

EIC code displayed by LXR

 
 

    


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

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