Geant4 Cross Reference |
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 // G4HelixExplicitEuler implementation 27 // 28 // Helix Explicit Euler: x_1 = x_0 + helix(h) 29 // with helix(h) being a helix piece of length h. 30 // Most simple approach for solving linear differential equations. 31 // Take the current derivative and add it to the current position. 32 // 33 // Author: W.Wander <wwc@mit.edu>, 12.09.1997 34 // ------------------------------------------------------------------- 35 36 #include "G4HelixExplicitEuler.hh" 37 #include "G4PhysicalConstants.hh" 38 #include "G4ThreeVector.hh" 39 40 G4HelixExplicitEuler::G4HelixExplicitEuler(G4Mag_EqRhs* EqRhs) 41 : G4MagHelicalStepper(EqRhs) 42 { 43 } 44 45 G4HelixExplicitEuler::~G4HelixExplicitEuler() = default; 46 47 void G4HelixExplicitEuler::Stepper( const G4double yInput[], // [7] 48 const G4double*, 49 G4double Step, 50 G4double yOut[], // [7] 51 G4double yErr[] ) 52 { 53 // Estimation of the Stepping Angle 54 // 55 G4ThreeVector Bfld; 56 MagFieldEvaluate(yInput, Bfld); 57 58 const G4int nvar = 6 ; 59 G4double yTemp[8], yIn[8] ; 60 G4ThreeVector Bfld_midpoint; 61 62 // Saving yInput because yInput and yOut can be aliases for same array 63 // 64 for(G4int i=0; i<nvar; ++i) 65 { 66 yIn[i] = yInput[i]; 67 } 68 69 G4double h = Step * 0.5; 70 71 // Do full step and two half steps 72 // 73 G4double yTemp2[7]; 74 AdvanceHelix(yIn, Bfld, h, yTemp2,yTemp); 75 MagFieldEvaluate(yTemp2, Bfld_midpoint) ; 76 AdvanceHelix(yTemp2, Bfld_midpoint, h, yOut); 77 SetAngCurve(GetAngCurve() * 2); 78 79 // Error estimation 80 // 81 for(G4int i=0; i<nvar; ++i) 82 { 83 yErr[i] = yOut[i] - yTemp[i]; 84 } 85 } 86 87 G4double G4HelixExplicitEuler::DistChord() const 88 { 89 // Implementation : must check whether h/R > 2 pi !! 90 // If( h/R < pi) use G4LineSection::DistLine 91 // Else DistChord=R_helix 92 // 93 G4double distChord; 94 G4double Ang_curve=GetAngCurve(); 95 96 97 if(Ang_curve<=pi) 98 { 99 distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve)); 100 } 101 else if(Ang_curve<twopi) 102 { 103 distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve))); 104 } 105 else 106 { 107 distChord=2.*GetRadHelix(); 108 } 109 110 return distChord; 111 } 112 113 void G4HelixExplicitEuler::DumbStepper( const G4double yIn[], 114 G4ThreeVector Bfld, 115 G4double h, 116 G4double yOut[] ) 117 { 118 AdvanceHelix(yIn, Bfld, h, yOut); 119 } 120