Back to home page

EIC code displayed by LXR

 
 

    


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

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 // G4FSALIntegrationDriver
0027 //
0028 // Class description:
0029 //
0030 // Driver class which controls the integration error of a Runge-Kutta stepper 
0031 
0032 // Created: D.Sorokin, 2017
0033 // --------------------------------------------------------------------
0034 #ifndef G4FSALINTEGRATIONDRIVER_HH
0035 #define G4FSALINTEGRATIONDRIVER_HH
0036 
0037 #include "G4RKIntegrationDriver.hh"
0038 #include "G4ChordFinderDelegate.hh"
0039 
0040 template <class T>
0041 class G4FSALIntegrationDriver : public G4RKIntegrationDriver<T>,
0042                        public G4ChordFinderDelegate<G4FSALIntegrationDriver<T>>
0043 {
0044   public:
0045 
0046     G4FSALIntegrationDriver(G4double hminimum,
0047                             T*       stepper,
0048                             G4int    numberOfComponents = 6,
0049                             G4int    statisticsVerbosity = 1);
0050 
0051    ~G4FSALIntegrationDriver() override;
0052 
0053     G4FSALIntegrationDriver(const G4FSALIntegrationDriver&) = delete;
0054     G4FSALIntegrationDriver& operator=(const G4FSALIntegrationDriver&) = delete;
0055 
0056     G4double AdvanceChordLimited(G4FieldTrack& track,
0057                                  G4double hstep,
0058                                  G4double eps,
0059                                  G4double chordDistance) override;
0060 
0061     void OnStartTracking() override
0062     {
0063       ChordFinderDelegate::ResetStepEstimate();
0064     }
0065 
0066     void OnComputeStep(const G4FieldTrack* /*track*/ = nullptr) override {}
0067 
0068     G4bool DoesReIntegrate() const override { return true; } 
0069    
0070     G4bool AccurateAdvance( G4FieldTrack& track,
0071                             G4double hstep,
0072                             G4double eps, // Requested y_err/hstep
0073                             G4double hinitial = 0.0) override;
0074       // Integrates ODE from current s (s=s0) to s=s0+h with accuracy eps.
0075       // On output track is replaced by value at end of interval.
0076       // The concept is similar to the odeint routine from NRC p.721-722.
0077 
0078     G4bool QuickAdvance(       G4FieldTrack& fieldTrack,
0079                          const G4double dydx[],
0080                                G4double hstep,
0081                                G4double& dchord_step,
0082                                G4double& dyerr ) override;
0083       // QuickAdvance just tries one Step - it does not ensure accuracy.
0084 
0085     void SetVerboseLevel(G4int newLevel) override;
0086     G4int GetVerboseLevel() const override;
0087 
0088     void StreamInfo( std::ostream& os ) const override;
0089      // Write out the parameters / state of the driver
0090    
0091     // Accessors
0092 
0093     G4double GetMinimumStep() const;
0094     void SetMinimumStep(G4double newval);
0095 
0096     void OneGoodStep(G4double y[],  // InOut
0097                      G4double dydx[],
0098                      G4double& curveLength,
0099                      G4double htry,
0100                      G4double eps,
0101                      G4double& hdid,
0102                      G4double& hnext);
0103       // This takes one Step that is of size htry, or as large 
0104       // as possible while satisfying the accuracy criterion of:
0105       //     yerr < eps * |y_end-y_start|
0106 
0107     G4double GetSmallestFraction() const;
0108     void SetSmallestFraction(G4double val);
0109 
0110   protected:
0111 
0112     void IncrementQuickAdvanceCalls();
0113 
0114   private:
0115 
0116     void CheckStep(const G4ThreeVector& posIn, 
0117                    const G4ThreeVector& posOut, 
0118                          G4double hdid);
0119 
0120     G4double fMinimumStep;
0121       // Minimum Step allowed in a Step (in absolute units)
0122 
0123     G4double fSmallestFraction{1e-12};
0124       // Smallest fraction of (existing) curve length - in relative units
0125       // below this fraction the current step will be the last
0126       // Expected range: smaller than 0.1 * epsilon and bigger than 5e-13
0127       //    ( Note: this range is not enforced. )
0128 
0129     G4int fVerboseLevel;
0130       // Verbosity level for printing (debug, ..)
0131       // Could be varied during tracking - to help identify issues
0132 
0133     G4int fNoQuickAvanceCalls{0};
0134     G4int fNoAccurateAdvanceCalls{0};
0135     G4int fNoAccurateAdvanceBadSteps{0};
0136     G4int fNoAccurateAdvanceGoodSteps{0};
0137 
0138     using Base = G4RKIntegrationDriver<T>;
0139     using ChordFinderDelegate = G4ChordFinderDelegate<G4FSALIntegrationDriver<T>>;
0140 };
0141 
0142 #include "G4FSALIntegrationDriver.icc"
0143 
0144 #endif