Back to home page

EIC code displayed by LXR

 
 

    


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

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 // G4BogackiShampine23
0027 //
0028 // Class description:
0029 //
0030 //  Bogacki-Shampine - 4 - 3(2) non-FSAL implementation 
0031 //
0032 //  An implementation of the embedded RK method from the paper 
0033 //  [1] P. Bogacki and L. F. Shampine,
0034 //     "A 3(2) pair of Runge - Kutta formulas"
0035 //     Appl. Math. Lett., vol. 2, no. 4, pp. 321-325, Jan. 1989.
0036 //
0037 //  This version does not utilise the FSAL property of the method,
0038 //  which would allow the reuse of the last derivative in the next step.
0039 //  (Alternative FSAL implementation created with revised interface)
0040 
0041 // Created: Somnath Banerjee, Google Summer of Code 2015, 20 May 2015
0042 // Supervision: John Apostolakis, CERN
0043 // --------------------------------------------------------------------
0044 #ifndef G4BOGACKI_SHAMPINE23_HH
0045 #define G4BOGACKI_SHAMPINE23_HH
0046 
0047 #include "G4MagIntegratorStepper.hh"
0048 #include "G4FieldTrack.hh"
0049 
0050 class G4BogackiShampine23 : public G4MagIntegratorStepper
0051 {
0052   public:
0053 
0054     G4BogackiShampine23(G4EquationOfMotion* EqRhs,
0055                         G4int numberOfVariables = 6);
0056 
0057     void Stepper(const G4double yInput[],
0058                  const G4double dydx[],
0059                        G4double hstep,
0060                        G4double yOutput[],
0061                        G4double yError[]) override;
0062 
0063     void Stepper(const G4double yInput[],
0064                  const G4double dydx[],
0065                        G4double hstep,
0066                        G4double yOutput[],
0067                        G4double yError[],
0068                        G4double dydxOutput[]);
0069 
0070     G4BogackiShampine23(const G4BogackiShampine23&) = delete;
0071     G4BogackiShampine23& operator = (const G4BogackiShampine23&) = delete;
0072 
0073     G4double DistChord() const override;
0074     G4int IntegratorOrder() const  override { return 3; }
0075 
0076   private:
0077 
0078     void makeStep(const G4double yInput[],
0079                   const G4double dydx[],
0080                   const G4double hstep,
0081                         G4double yOutput[],
0082                         G4double* dydxOutput = nullptr,
0083                         G4double* yError = nullptr) const;
0084 
0085     G4double fyIn[G4FieldTrack::ncompSVEC],
0086              fdydx[G4FieldTrack::ncompSVEC],
0087              fyOut[G4FieldTrack::ncompSVEC],
0088              fdydxOut[G4FieldTrack::ncompSVEC];
0089     G4double fhstep = -1.0;
0090 };
0091 
0092 #endif