|
||||
File indexing completed on 2025-01-18 09:58:10
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 // G4DormandPrinceRK78 0027 // 0028 // Class description: 0029 // 0030 // Dormand-Prince 8(7)13M non-FSAL RK method, a 13 stage embedded 0031 // explicit Runge-Kutta method, using a pair of 7th and 8th order formulae. 0032 // 0033 // Paper proposing this RK scheme: 0034 // P.J. Prince, J.R. Dormand, "High order embedded Runge-Kutta formulae", 0035 // Journal of Computational and Applied Mathematics, Volume 7, Issue 1, 1981, 0036 // Pages 67-75, ISSN 0377-0427, DOI: 10.1016/0771-050X(81)90010-3 0037 0038 // Created: Somnath Banerjee, Google Summer of Code 2015, 28 June 2015 0039 // Supervision: John Apostolakis, CERN 0040 // -------------------------------------------------------------------- 0041 #ifndef G4DORMAND_PRINCE_RK78_HH 0042 #define G4DORMAND_PRINCE_RK78_HH 0043 0044 #include "G4MagIntegratorStepper.hh" 0045 0046 class G4DormandPrinceRK78 : public G4MagIntegratorStepper 0047 { 0048 public: 0049 0050 G4DormandPrinceRK78(G4EquationOfMotion* EqRhs, 0051 G4int numberOfVariables = 6, 0052 G4bool primary = true); 0053 ~G4DormandPrinceRK78() override; 0054 0055 G4DormandPrinceRK78(const G4DormandPrinceRK78&) = delete; 0056 G4DormandPrinceRK78& operator=(const G4DormandPrinceRK78&) = delete; 0057 0058 void Stepper( const G4double y[], 0059 const G4double dydx[], 0060 G4double h, 0061 G4double yout[], 0062 G4double yerr[]) override ; 0063 0064 G4double DistChord() const override; 0065 inline G4int IntegratorOrder() const override { return 7; } 0066 0067 private : 0068 0069 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8, 0070 *ak9, *ak10, *ak11, *ak12, *ak13, 0071 *yTemp, *yIn; 0072 0073 G4double fLastStepLength = -1.0; 0074 G4double *fLastInitialVector, *fLastFinalVector, 0075 *fLastDyDx, *fMidVector, *fMidError; 0076 // For DistChord calculations 0077 0078 G4DormandPrinceRK78* fAuxStepper = nullptr; 0079 }; 0080 0081 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |