Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:03

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 // G4RKG3_Stepper
0027 //
0028 // Class description:
0029 //
0030 // Runga-Kutta integrator stepper from Geant3.
0031 
0032 // Created: J.Apostolakis, V.Grichine - 30.01.1997
0033 // -------------------------------------------------------------------
0034 #ifndef G4RKG3_STEPPER_HH
0035 #define G4RKG3_STEPPER_HH
0036 
0037 #include "G4Types.hh"
0038 #include "G4MagIntegratorStepper.hh"
0039 #include "G4ThreeVector.hh"
0040 
0041 class G4Mag_EqRhs;
0042 
0043 class G4RKG3_Stepper : public G4MagIntegratorStepper
0044 {
0045   public:
0046 
0047     G4RKG3_Stepper(G4Mag_EqRhs* EqRhs);
0048       // Integrate over 6 variables only:  position & velocity.
0049       // Not implemented yet !
0050 
0051     ~G4RKG3_Stepper() override;
0052 
0053     void Stepper( const G4double yIn[],
0054                   const G4double dydx[],
0055                         G4double h,
0056                         G4double yOut[],
0057                         G4double yErr[] ) override;
0058       // The method which must be provided, even if less efficient.
0059 
0060     G4double  DistChord() const override ;
0061  
0062     void StepNoErr( const G4double tIn[8],
0063                     const G4double dydx[6],
0064                           G4double Step,
0065                           G4double tOut[8],
0066                           G4double B[3] );
0067       // Integrator RK Stepper from G3 with only two field evaluation per 
0068       // Step. It is used in propagation initial Step by small substeps
0069       // after solution error and delta geometry considerations. 
0070       // B[3] is magnetic field which is passed from substep to substep.
0071 
0072     void StepWithEst( const G4double  tIn[8],
0073                       const G4double dydx[6],
0074                             G4double Step,
0075                             G4double tOut[8],
0076                             G4double& alpha2,
0077                             G4double& beta2,
0078                       const G4double B1[3],
0079                             G4double B2[3] );
0080       // Integrator for RK from G3 with evaluation of error in solution and
0081       // delta geometry based on naive similarity with the case of uniform
0082       // magnetic field.
0083       // B1[3] is input  and is the first magnetic field values
0084       // B2[3] is output and is the final magnetic field values.
0085 
0086     inline G4int IntegratorOrder() const override { return 4; }
0087 
0088   private:
0089 
0090     G4ThreeVector fyInitial,
0091                   fyMidPoint,
0092                   fyFinal;
0093     G4ThreeVector fpInitial;
0094     G4ThreeVector BfldIn;
0095     G4double      hStep = 0.0;
0096 };
0097 
0098 #endif