Back to home page

EIC code displayed by LXR

 
 

    


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

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