Back to home page

EIC code displayed by LXR

 
 

    


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

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 G4PropagatorInField 
0027 //
0028 // class description:
0029 // 
0030 // This class performs the navigation/propagation of a particle/track 
0031 // in a magnetic field. The field is in general non-uniform.
0032 // For the calculation of the path, it relies on the class G4ChordFinder.
0033 
0034 // History:
0035 // -------
0036 // 25.10.96 John Apostolakis,  design and implementation 
0037 // 25.03.97 John Apostolakis,  adaptation for G4Transportation and cleanup
0038 //  8.11.02 John Apostolakis,  changes to enable use of safety in intersecting
0039 // ---------------------------------------------------------------------------
0040 #ifndef G4PropagatorInField_hh 
0041 #define G4PropagatorInField_hh  1
0042 
0043 #include "G4Types.hh"
0044 
0045 #include <vector>
0046 
0047 #include "G4FieldTrack.hh"
0048 #include "G4FieldManager.hh"
0049 #include "G4VIntersectionLocator.hh"
0050 
0051 class G4ChordFinder; 
0052 
0053 class G4Navigator;
0054 class G4VPhysicalVolume;
0055 class G4VCurvedTrajectoryFilter;
0056 
0057 class G4PropagatorInField
0058 {
0059 
0060  public:  // with description
0061 
0062    G4PropagatorInField( G4Navigator* theNavigator, 
0063                         G4FieldManager* detectorFieldMgr,
0064                         G4VIntersectionLocator* vLocator = nullptr );
0065   ~G4PropagatorInField();
0066 
0067    G4double ComputeStep( G4FieldTrack& pFieldTrack,
0068                          G4double pCurrentProposedStepLength,
0069                          G4double& pNewSafety, 
0070                          G4VPhysicalVolume* pPhysVol = nullptr,
0071                          G4bool canRelaxDeltaChord = false);
0072      // Compute the next geometric Step
0073 
0074    inline G4ThreeVector EndPosition() const;       
0075    inline G4ThreeVector EndMomentumDir() const;
0076    inline G4bool        IsParticleLooping() const;
0077      // Return the state after the Step
0078 
0079    inline G4double GetEpsilonStep() const;
0080      // Relative accuracy for current Step (Calc.)
0081    inline void     SetEpsilonStep(G4double newEps);
0082      // The ratio DeltaOneStep()/h_current_step
0083 
0084    G4FieldManager* FindAndSetFieldManager(G4VPhysicalVolume* pCurrentPhysVol);
0085      // Set (and return) the correct field manager (global or local), 
0086      // if it exists.
0087      // Should be called before ComputeStep is called;
0088      // Currently, ComputeStep will call it, if it has not been called.
0089  
0090    inline G4ChordFinder* GetChordFinder();
0091 
0092           G4int SetVerboseLevel( G4int verbose );
0093    inline G4int GetVerboseLevel() const;
0094    inline G4int Verbose() const;
0095    inline void CheckMode(G4bool mode);
0096 
0097    inline void   SetVerboseTrace( G4bool enable );
0098    inline G4bool GetVerboseTrace();
0099    // Tracing key parts of Compute Step
0100    
0101    inline G4int GetMaxLoopCount() const;
0102    inline void  SetMaxLoopCount( G4int new_max );
0103      // A maximum for the number of substeps that a particle can take.
0104      //   Above this number it is signaled as 'looping'.
0105 
0106    void printStatus( const G4FieldTrack&      startFT,
0107                      const G4FieldTrack&      currentFT, 
0108                            G4double           requestStep, 
0109                            G4double           safety,
0110                            G4int              step, 
0111                            G4VPhysicalVolume* startVolume);
0112      // Print Method - useful mostly for debugging.
0113 
0114    inline G4FieldTrack GetEndState() const;
0115 
0116    inline G4double GetMinimumEpsilonStep() const; // Min for relative accuracy
0117    inline void     SetMinimumEpsilonStep( G4double newEpsMin ); // of any step
0118    inline G4double GetMaximumEpsilonStep() const;
0119    inline void     SetMaximumEpsilonStep( G4double newEpsMax );
0120      // The 4 above methods are now obsolescent but *for now* will work 
0121      // They are being replaced by same-name methods in G4FieldManager,
0122      // allowing the specialisation in different volumes. 
0123      // Their new behaviour is to change the values for the global field
0124      // manager
0125 
0126    void     SetLargestAcceptableStep( G4double newBigDist );
0127    G4double GetLargestAcceptableStep();
0128    void     ResetLargestAcceptableStep();
0129      // Obtain / change the size of the largest step the method will undertake
0130      // Reset method uses the world volume's 
0131 
0132    G4double GetMaxStepSizeMultiplier();
0133    void     SetMaxStepSizeMultiplier(G4double vm);
0134      // Control extra Multiplier parameter for limiting long steps.
0135    G4double GetMinBigDistance();
0136    void     SetMinBigDistance(G4double val);
0137      // Control minimum 'directional' distance in case of too-large step
0138 
0139    void SetTrajectoryFilter(G4VCurvedTrajectoryFilter* filter);
0140      // Set the filter that examines & stores 'intermediate' 
0141      // curved trajectory points.  Currently only position is stored.
0142 
0143    std::vector<G4ThreeVector>* GimmeTrajectoryVectorAndForgetIt() const;
0144      // Access the points which have passed by the filter.
0145      // Responsibility for deleting the points lies with the client.
0146      // This method MUST BE called exactly ONCE per step. 
0147 
0148    void ClearPropagatorState();
0149      // Clear all the State of this class and its current associates
0150      // --> the current field manager & chord finder will also be called
0151 
0152    inline void SetDetectorFieldManager( G4FieldManager* newGlobalFieldManager );
0153      // Update this (dangerous) state -- for the time being
0154   
0155    inline void   SetUseSafetyForOptimization( G4bool );
0156    inline G4bool GetUseSafetyForOptimization();
0157      // Toggle & view parameter for using safety to discard 
0158      // unneccesary calls to navigator (thus 'optimising' performance)
0159    inline G4bool IntersectChord( const G4ThreeVector& StartPointA,
0160                                  const G4ThreeVector& EndPointB,
0161                                        G4double&      NewSafety,
0162                                        G4double&      LinearStepLength,
0163                                        G4ThreeVector& IntersectionPoint);
0164      // Intersect the chord from StartPointA to EndPointB
0165      // and return whether an intersection occurred
0166      // NOTE: Safety is changed!
0167 
0168    inline G4bool IsFirstStepInVolume();
0169    inline G4bool IsLastStepInVolume();
0170    inline void PrepareNewTrack();
0171       
0172    inline G4VIntersectionLocator* GetIntersectionLocator();
0173    inline void SetIntersectionLocator(G4VIntersectionLocator* pLocator );
0174      // Change or get the object which calculates the exact 
0175      // intersection point with the next boundary
0176 
0177    inline G4int GetIterationsToIncreaseChordDistance() const;
0178    inline void  SetIterationsToIncreaseChordDistance(G4int numIters);
0179      // Control the parameter which enables the temporary 'relaxation'
0180      //   which ensures that chord segments are short enough so that
0181      //   their sagitta is small than delta-chord parameter.
0182      // The Set method increases the value of delta-chord temporarily,
0183      //   doubling it once the number of iterations substeps reach
0184      //   value of 'IncreaseChordDistanceThreshold'.  It is also doubled
0185      //   again every time the iteration count reaches a multiple of this
0186      //   value.
0187      // Note: delta-chord is reset to its original value at the end of
0188      //   each call to ComputeStep.
0189 
0190  public:  // without description
0191 
0192    inline G4double GetDeltaIntersection() const;
0193    inline G4double GetDeltaOneStep() const;
0194 
0195    inline G4FieldManager* GetCurrentFieldManager();
0196    inline G4EquationOfMotion* GetCurrentEquationOfMotion();
0197       // Auxiliary methods - their results can/will change during propagation
0198 
0199    inline void SetNavigatorForPropagating(G4Navigator* SimpleOrMultiNavigator); 
0200    inline G4Navigator* GetNavigatorForPropagating();
0201 
0202    inline void SetThresholdNoZeroStep( G4int noAct,
0203                                        G4int noHarsh,
0204                                        G4int noAbandon );
0205    inline G4int GetThresholdNoZeroSteps( G4int i ); 
0206 
0207    inline G4double GetZeroStepThreshold(); 
0208    inline void     SetZeroStepThreshold( G4double newLength ); 
0209    
0210    void RefreshIntersectionLocator(); 
0211      // Update the Locator with parameters from this class
0212      // and from current field manager
0213 
0214  protected:  // without description
0215 
0216    void PrintStepLengthDiagnostic( G4double      currentProposedStepLength,
0217                                    G4double      decreaseFactor,
0218                                    G4double      stepTrial,
0219                              const G4FieldTrack& aFieldTrack);
0220 
0221    void ReportLoopingParticle( G4int count,  G4double StepTaken,
0222                                G4double stepRequest, const char* methodName,
0223                                const G4ThreeVector&      momentumVec,
0224                                G4VPhysicalVolume* physVol);
0225    void ReportStuckParticle(G4int noZeroSteps, G4double proposedStep,
0226                             G4double lastTriedStep, G4VPhysicalVolume* physVol);
0227 
0228  private:
0229 
0230    // ----------------------------------------------------------------------
0231    //  DATA Members
0232    // ----------------------------------------------------------------------
0233 
0234    //  ==================================================================
0235    //  INVARIANTS - Must not change during tracking
0236 
0237    //  ** PARAMETERS -----------
0238    G4int fMax_loop_count = 1000;
0239      // Limit for the number of sub-steps taken in one call to ComputeStep
0240    G4int fIncreaseChordDistanceThreshold = 100;
0241    G4bool fUseSafetyForOptimisation = true;
0242      // (false) is less sensitive to incorrect safety
0243 
0244    //  Thresholds for identifying "abnormal" cases - which cause looping
0245    //
0246    G4int fActionThreshold_NoZeroSteps = 2;        // Threshold # - above it act
0247    G4int fSevereActionThreshold_NoZeroSteps = 10; // Threshold # to act harshly
0248    G4int fAbandonThreshold_NoZeroSteps = 50;      // Threshold # to abandon
0249    G4double fZeroStepThreshold = 0.0; 
0250      // Threshold *length* for counting of tiny or 'zero' steps 
0251 
0252    // Parameters related to handling of very large steps which
0253    //   occur typically in large volumes with vacuum or very thin gas
0254    G4double fLargestAcceptableStep;
0255      // Maximum size of a step - for optimization (and to avoid problems)
0256    G4double fMaxStepSizeMultiplier = 3;
0257      // Multiplier for directional exit distance used as extra long-step limit 
0258    G4double fMinBigDistance= 100. ; // * CLHEP::mm
0259      // Minimum distance added to directional exit distance
0260    //  ** End of PARAMETERS -----
0261 
0262    G4double kCarTolerance;
0263        // Geometrical tolerance defining surface thickness
0264 
0265    G4bool fAllocatedLocator;                    //  Book-keeping
0266 
0267    //  --------------------------------------------------------
0268    //  ** Dependent Objects - to which work is delegated 
0269 
0270    G4FieldManager* fDetectorFieldMgr; 
0271        // The  Field Manager of the whole Detector.  (default)
0272 
0273    G4VIntersectionLocator* fIntersectionLocator;
0274      // Refines candidate intersection
0275 
0276    G4VCurvedTrajectoryFilter* fpTrajectoryFilter = nullptr;
0277      // The filter encapsulates the algorithm which selects which
0278      // intermediate points should be stored in a trajectory. 
0279      // When it is NULL, no intermediate points will be stored.
0280      // Else PIF::ComputeStep must submit (all) intermediate
0281      // points it calculates, to this filter.  (jacek 04/11/2002)
0282 
0283    G4Navigator* fNavigator;
0284      // Set externally - only by tracking / run manager
0285    //
0286    //  ** End of Dependent Objects ----------------------------
0287 
0288    //  End of INVARIANTS 
0289    //  ==================================================================
0290 
0291    //  STATE information
0292    //  -----------------
0293    G4FieldManager* fCurrentFieldMgr;
0294      // The  Field Manager of the current volume (may be the global)
0295    G4bool fSetFieldMgr = false;  // Has it been set for the current step?
0296 
0297    // Parameters of current step
0298    G4double fEpsilonStep;            // Relative accuracy of current Step
0299    G4FieldTrack End_PointAndTangent; // End point storage
0300    G4bool fParticleIsLooping = false;
0301    G4int fNoZeroStep = 0;            // Count of zero Steps
0302 
0303    // State used for Optimisation
0304    G4double fFull_CurveLen_of_LastAttempt = -1; 
0305    G4double fLast_ProposedStepLength = -1; 
0306      // Previous step information -- for use in adjust step size
0307    G4ThreeVector fPreviousSftOrigin;
0308    G4double fPreviousSafety = 0.0; 
0309      // Last safety origin & value: for optimisation
0310 
0311    G4int fVerboseLevel = 0;
0312    G4bool fVerbTracePiF = false;
0313    G4bool fCheck = false;
0314      // For debugging purposes
0315 
0316    G4bool fFirstStepInVolume = true; 
0317    G4bool fLastStepInVolume = true; 
0318    G4bool fNewTrack = true;
0319 };
0320 
0321 // Inline methods
0322 //
0323 #include "G4PropagatorInField.icc"
0324 
0325 #endif