Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/magneticfield/include/G4DormandPrince745.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 // G4DormandPrince745
 27 //
 28 // Class desription:
 29 //
 30 //  An implementation of the 5th order embedded RK method from the paper:
 31 //  J. R. Dormand and P. J. Prince, "A family of embedded Runge-Kutta formulae"
 32 //  Journal of computational and applied Math., vol.6, no.1, pp.19-26, 1980.
 33 //
 34 //  DormandPrince7 - 5(4) embedded RK method
 35 //
 36 
 37 // Created: Somnath Banerjee, Google Summer of Code 2015, 25 May 2015
 38 // Supervision: John Apostolakis, CERN
 39 // --------------------------------------------------------------------
 40 #ifndef G4DORMAND_PRINCE_745_HH
 41 #define G4DORMAND_PRINCE_745_HH
 42 
 43 #include "G4MagIntegratorStepper.hh"
 44 #include "G4FieldUtils.hh"
 45 
 46 class G4DormandPrince745 : public G4MagIntegratorStepper
 47 {
 48   public:
 49 
 50     G4DormandPrince745(G4EquationOfMotion* equation,
 51                        G4int numberOfVariables = 6);
 52 
 53     void Stepper(const G4double yInput[],
 54                  const G4double dydx[],
 55                        G4double hstep,
 56                        G4double yOutput[],
 57                        G4double yError[]) override;
 58 
 59     void Stepper(const G4double yInput[],
 60                  const G4double dydx[],
 61                        G4double hstep,
 62                        G4double yOutput[],
 63                        G4double yError[],
 64                        G4double dydxOutput[]);
 65 
 66     inline void SetupInterpolation() {}
 67 
 68     inline void Interpolate(G4double tau, G4double yOut[]) const
 69     {
 70       Interpolate4thOrder(yOut, tau);
 71     }
 72       // For calculating the output at the tau fraction of Step
 73 
 74     G4double DistChord() const override;
 75 
 76     G4int IntegratorOrder() const override { return 4; }
 77 
 78     const G4String& StepperType() const;
 79     const G4String& StepperDescription() const;
 80    
 81     const field_utils::State& GetYOut() const { return fyOut; }
 82 
 83     void Interpolate4thOrder(G4double yOut[], G4double tau) const;
 84 
 85     void SetupInterpolation5thOrder();
 86     void Interpolate5thOrder(G4double yOut[], G4double tau) const;
 87 
 88     G4EquationOfMotion* GetSpecificEquation() { return GetEquationOfMotion(); }
 89 
 90   private:
 91 
 92     field_utils::State ak2, ak3, ak4, ak5, ak6, ak7, ak8, ak9;
 93     field_utils::State fyIn, fyOut, fdydxIn;
 94 
 95     G4double fLastStepLength = -1.0;
 96 };
 97 
 98 #endif
 99