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 // G4BogackiShampine45
0027 //
0028 // Class description:
0029 //
0030 //  An implementation of the embedded RK method from the following paper 
0031 //  by P. Bogacki and L. F. Shampine:
0032 //     "An efficient Runge-Kutta (4,5) pair"
0033 //     Comput. Math. with Appl., vol. 32, no. 6, pp. 15-28, Sep. 1996.
0034 //
0035 //  An interpolation method provides the value of an intermediate
0036 //  point in a step -- if a step was sucessful.
0037 //
0038 //  This version can provide the FSAL property of the method,
0039 //  which allows the reuse of the last derivative in the next step,
0040 //  but only by using the additional method GetLastDyDx() (an alternative
0041 //  interface for simpler use of FSAL is under development).
0042 
0043 // Created: Somnath Banerjee, Google Summer of Code 2015, 25 May 2015
0044 // Supervision: John Apostolakis, CERN
0045 // --------------------------------------------------------------------
0046 #ifndef BOGACKI_SHAMPINE_45_HH
0047 #define BOGACKI_SHAMPINE_45_HH
0048 
0049 #include "G4MagIntegratorStepper.hh"
0050 
0051 class G4BogackiShampine45 : public G4MagIntegratorStepper
0052 {
0053   public:
0054 
0055     G4BogackiShampine45(G4EquationOfMotion* EqRhs,
0056                         G4int numberOfVariables = 6,
0057                         G4bool primary =  true);
0058     ~G4BogackiShampine45() override;
0059     
0060     G4BogackiShampine45(const G4BogackiShampine45&) = delete;
0061     G4BogackiShampine45& operator=(const G4BogackiShampine45&) = delete;
0062 
0063     void Stepper( const G4double y[],
0064                   const G4double dydx[],
0065                         G4double h,
0066                         G4double yout[],
0067                         G4double yerr[] ) override ;
0068 
0069     // This Stepper provides 'dense output'. After a successful
0070     // step, it is possible to obtain an estimate of the value
0071     // of the function at an intermediate point of the interval.
0072     // This requires only two additional evaluations of the
0073     // derivative (and thus the field).
0074 
0075     inline void SetupInterpolation()
0076     {
0077       SetupInterpolationHigh(); // ( yInput, dydx, Step);
0078     }
0079     
0080     // For calculating the output at the tau fraction of Step
0081     //
0082     inline void Interpolate( G4double tau, 
0083                              G4double yOut[] ) // Output value
0084     {
0085       InterpolateHigh( tau, yOut);       
0086       // InterpolateHigh( yInput, dydx, Step, yOut, tau);
0087     }
0088     
0089     void SetupInterpolationHigh();
0090     
0091     // For calculating the output at the tau fraction of Step
0092     //
0093     void InterpolateHigh( G4double tau,
0094                           G4double yOut[] ) const;
0095     
0096     G4double  DistChord()   const override;
0097     G4int IntegratorOrder() const override { return 4; }
0098 
0099     void GetLastDydx( G4double dyDxLast[] );
0100 
0101     void PrepareConstants();  // Initialise the values of the bi[][] array
0102    
0103   private:
0104     
0105     G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8,
0106              *ak9, *ak10, *ak11, *yTemp, *yIn;
0107 
0108     G4double *p[6];
0109 
0110     G4double fLastStepLength = -1.0;
0111     G4double *fLastInitialVector, *fLastFinalVector, *fLastDyDx,
0112              *fMidVector, *fMidError;
0113       // For DistChord calculations
0114     
0115     G4BogackiShampine45* fAuxStepper = nullptr;
0116       // For chord - until interpolation is proven
0117     G4bool fPreparedInterpolation = false;
0118 
0119     // Class constants
0120     static G4bool fPreparedConstants;   
0121     static G4double bi[12][7];
0122 };
0123 
0124 #endif