|
||||
File indexing completed on 2025-01-18 09:58:23
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 // G4HelixMixedStepper 0027 // 0028 // Class description: 0029 // 0030 // G4HelixMixedStepper split the Method used for Integration in two: 0031 // 0032 // If Stepping Angle ( h / R_curve) < pi/3 : use Stepper for small step 0033 // 0034 // Else use HelixExplicitEuler Stepper 0035 // 0036 // Stepper for the small step is G4ClassicalRK4 by default, but 0037 // it possible to choose other stepper,like G4CashKarpRK45 or G4RKG3_Stepper, 0038 // by setting StepperNumber : new HelixMixedStepper(EqRhs,N) 0039 // 0040 // N=2 G4SimpleRunge; N=3 G4SimpleHeum; 0041 // N=4 G4ClassicalRK4; 0042 // N=6 G4HelixImplicitEuler; N=7 G4HelixSimpleRunge; 0043 // N=8 G4CashKarpRK45; N=9 G4ExactHelixStepper; 0044 // N=10 G4RKG3_Stepper; N=13 G4NystromRK4 0045 // N=23 BogackiShampine23 N=145 TsitourasRK45 0046 // N=45 BogackiShampine45 N=745 DormandPrince745 (ie DoPri5) 0047 // 0048 // For completeness also available are: 0049 // N=11 G4ExplicitEuler N=12 G4ImplicitEuler; -- Likely poor 0050 // N=5 G4HelixExplicitEuler (testing only) 0051 // For recommendations see comments in 'SetupStepper' method. 0052 // 0053 // Note: Like other helix steppers, only applicable in pure magnetic field 0054 0055 // Created: T.Nikitina, CERN - 18.05.2007, derived from G4ExactHelicalStepper 0056 // ------------------------------------------------------------------- 0057 #ifndef G4HELIXMIXEDSTEPPER_HH 0058 #define G4HELIXMIXEDSTEPPER_HH 0059 0060 #include "G4MagHelicalStepper.hh" 0061 0062 class G4HelixMixedStepper : public G4MagHelicalStepper 0063 { 0064 public: 0065 0066 G4HelixMixedStepper(G4Mag_EqRhs* EqRhs, 0067 G4int StepperNumber = -1, 0068 G4double Angle_threshold = -1.0); 0069 ~G4HelixMixedStepper() override; 0070 0071 void Stepper( const G4double y[], 0072 const G4double dydx[], 0073 G4double h, 0074 G4double yout[], 0075 G4double yerr[] ) override; 0076 // Step 'integration' for step size 'h' 0077 // If SteppingAngle= h/R_curve < pi/3 uses default RK stepper 0078 // else use Helix Fast Method 0079 0080 void DumbStepper( const G4double y[], 0081 G4ThreeVector Bfld, 0082 G4double h, 0083 G4double yout[]) override; 0084 0085 G4double DistChord() const override; 0086 // Estimate maximum distance of curved solution and chord ... 0087 0088 inline void SetVerbose (G4int newvalue) { fVerbose = newvalue; } 0089 0090 void PrintCalls(); 0091 G4MagIntegratorStepper* SetupStepper(G4Mag_EqRhs* EqRhs, G4int StepperName); 0092 0093 inline void SetAngleThreshold( G4double val ) { fAngle_threshold = val; } 0094 inline G4double GetAngleThreshold() { return fAngle_threshold; } 0095 inline G4int IntegratorOrder() const override { return 4; } 0096 0097 private: 0098 0099 G4MagIntegratorStepper* fRK4Stepper = nullptr; 0100 // Mixed Integration RK4 for 'small' steps 0101 G4int fStepperNumber = -1; 0102 // Int ID of RK stepper 0103 G4double fAngle_threshold = -1.0; 0104 // Threshold angle (in radians ) - above it Helical stepper is used 0105 0106 private: 0107 0108 G4int fVerbose = 0; 0109 0110 G4int fNumCallsRK4 = 0; 0111 G4int fNumCallsHelix = 0; 0112 // Used for statistic = how many calls to different steppers 0113 }; 0114 0115 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |