Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/magneticfield/include/G4BogackiShampine23.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 // G4BogackiShampine23
 27 //
 28 // Class description:
 29 //
 30 //  Bogacki-Shampine - 4 - 3(2) non-FSAL implementation 
 31 //
 32 //  An implementation of the embedded RK method from the paper 
 33 //  [1] P. Bogacki and L. F. Shampine,
 34 //     "A 3(2) pair of Runge - Kutta formulas"
 35 //     Appl. Math. Lett., vol. 2, no. 4, pp. 321-325, Jan. 1989.
 36 //
 37 //  This version does not utilise the FSAL property of the method,
 38 //  which would allow the reuse of the last derivative in the next step.
 39 //  (Alternative FSAL implementation created with revised interface)
 40 
 41 // Created: Somnath Banerjee, Google Summer of Code 2015, 20 May 2015
 42 // Supervision: John Apostolakis, CERN
 43 // --------------------------------------------------------------------
 44 #ifndef G4BOGACKI_SHAMPINE23_HH
 45 #define G4BOGACKI_SHAMPINE23_HH
 46 
 47 #include "G4MagIntegratorStepper.hh"
 48 #include "G4FieldTrack.hh"
 49 
 50 class G4BogackiShampine23 : public G4MagIntegratorStepper
 51 {
 52   public:
 53 
 54     G4BogackiShampine23(G4EquationOfMotion* EqRhs,
 55                         G4int numberOfVariables = 6);
 56 
 57     void Stepper(const G4double yInput[],
 58                  const G4double dydx[],
 59                        G4double hstep,
 60                        G4double yOutput[],
 61                        G4double yError[]) override;
 62 
 63     void Stepper(const G4double yInput[],
 64                  const G4double dydx[],
 65                        G4double hstep,
 66                        G4double yOutput[],
 67                        G4double yError[],
 68                        G4double dydxOutput[]);
 69 
 70     G4BogackiShampine23(const G4BogackiShampine23&) = delete;
 71     G4BogackiShampine23& operator = (const G4BogackiShampine23&) = delete;
 72 
 73     G4double DistChord() const override;
 74     G4int IntegratorOrder() const  override { return 3; }
 75 
 76   private:
 77 
 78     void makeStep(const G4double yInput[],
 79                   const G4double dydx[],
 80                   const G4double hstep,
 81                         G4double yOutput[],
 82                         G4double* dydxOutput = nullptr,
 83                         G4double* yError = nullptr) const;
 84 
 85     G4double fyIn[G4FieldTrack::ncompSVEC],
 86              fdydx[G4FieldTrack::ncompSVEC],
 87              fyOut[G4FieldTrack::ncompSVEC],
 88              fdydxOut[G4FieldTrack::ncompSVEC];
 89     G4double fhstep = -1.0;
 90 };
 91 
 92 #endif
 93