Back to home page

EIC code displayed by LXR

 
 

    


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

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 // G4FSALDormandPrince745
0027 //
0028 // Class description:
0029 //
0030 // DormandPrince7 - 5(4) FSAL stepper
0031 
0032 // Created: Somnath Banerjee, Google Summer of Code 2015, 25 May 2015
0033 // Supervision: John Apostolakis, CERN
0034 // --------------------------------------------------------------------
0035 #ifndef G4FSALDORMANDPRINCE745_HH
0036 #define G4FSALDORMANDPRINCE745_HH
0037 
0038 #include "G4VFSALIntegrationStepper.hh"
0039 
0040 class G4FSALDormandPrince745 : public G4VFSALIntegrationStepper
0041 {
0042   public:
0043 
0044     G4FSALDormandPrince745(G4EquationOfMotion* EqRhs,
0045                            G4int numberOfVariables = 6,
0046                            G4bool primary = true);
0047    ~G4FSALDormandPrince745() override;
0048 
0049     G4FSALDormandPrince745(const G4FSALDormandPrince745&) = delete;
0050     G4FSALDormandPrince745& operator=(const G4FSALDormandPrince745&) = delete;
0051 
0052     void Stepper( const G4double y[],
0053                   const G4double dydx[],
0054                         G4double h,
0055                         G4double yout[],
0056                         G4double yerr[],
0057                         G4double nextDydx[]) override ;
0058     void interpolate( const G4double yInput[],
0059                       const G4double dydx[],
0060                             G4double yOut[],
0061                             G4double Step,
0062                             G4double tau ) ;
0063 
0064     void SetupInterpolate( const G4double yInput[],
0065                            const G4double dydx[],
0066                            const G4double Step );
0067       // For higher order Interpolant
0068     
0069     void Interpolate( const G4double yInput[],
0070                       const G4double dydx[],
0071                       const G4double Step,
0072                             G4double yOut[],
0073                             G4double tau );
0074       // For calculating the output at the tau fraction of Step
0075     
0076 
0077     G4double  DistChord()   const override;
0078     inline G4int IntegratorOrder() const override { return 4; }
0079     inline G4bool isFSAL() const { return true; }
0080     
0081   private:
0082     
0083     G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7,
0084              *ak8, *ak9,         // For additional stages in the interpolant
0085              *yTemp, *yIn;
0086     
0087     G4double* pseudoDydx_for_DistChord;
0088       // Only for use with DistChord()
0089     
0090     G4double fLastStepLength = -1.0;
0091     G4double *fLastInitialVector, *fLastFinalVector,
0092              *fInitialDyDx, *fLastDyDx, *fMidVector, *fMidError;
0093       // For DistChord() calculations
0094     
0095     G4FSALDormandPrince745* fAuxStepper = nullptr;
0096 };
0097 
0098 #endif