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 // G4TMagErrorStepper 27 // 28 // Class description: 29 // 30 // Templated version of G4MagErrorStepper 31 // 32 // 33 // Created: Josh Xie (supported by Google Summer of Code 2014 ) 34 // Supervisors: Sandro Wenzel, John Apostolakis (CERN) 35 // Adapted from G4G4TMagErrorStepper class 36 // -------------------------------------------------------------------- 37 #ifndef G4TMAG_ERROR_STEPPER_HH 38 #define G4TMAG_ERROR_STEPPER_HH 39 40 #include "G4Types.hh" 41 #include "G4MagIntegratorStepper.hh" 42 #include "G4ThreeVector.hh" 43 #include "G4LineSection.hh" 44 45 template <class T_Stepper, class T_Equation, unsigned int N> 46 class G4TMagErrorStepper : public G4MagIntegratorStepper 47 { 48 public: // with description 49 G4TMagErrorStepper(T_Equation* EqRhs, G4int numberOfVariables, 50 G4int numStateVariables = 12) 51 : G4MagIntegratorStepper(EqRhs, numberOfVariables, numStateVariables) 52 , fEquation_Rhs(EqRhs) 53 { 54 // G4int nvar = std::max(this->GetNumberOfVariables(), 8); 55 } 56 57 virtual ~G4TMagErrorStepper() { ; } 58 59 inline void RightHandSide(G4double y[], G4double dydx[]) 60 { 61 fEquation_Rhs->T_Equation::RightHandSide(y, dydx); 62 } 63 64 inline void Stepper(const G4double yInput[], const G4double dydx[], 65 G4double hstep, G4double yOutput[], G4double yError[]) override final; 66 67 inline G4double DistChord() const override final; 68 69 private: 70 G4TMagErrorStepper(const G4TMagErrorStepper&); 71 G4TMagErrorStepper& operator=(const G4TMagErrorStepper&); 72 // Private copy constructor and assignment operator. 73 74 private: 75 // STATE 76 G4ThreeVector fInitialPoint, fMidPoint, fFinalPoint; 77 // Data stored in order to find the chord 78 79 // Dependent Objects, owned --- part of the STATE 80 G4double yInitial[N < 8 ? 8 : N]; 81 G4double yMiddle[N < 8 ? 8 : N]; 82 G4double dydxMid[N < 8 ? 8 : N]; 83 G4double yOneStep[N < 8 ? 8 : N]; 84 // The following arrays are used only for temporary storage 85 // they are allocated at the class level only for efficiency - 86 // so that calls to new and delete are not made in Stepper(). 87 88 T_Equation* fEquation_Rhs; 89 }; 90 91 // ------------ Implementation ----------------------- 92 93 template <class T_Stepper, class T_Equation, unsigned int N > 94 void G4TMagErrorStepper<T_Stepper,T_Equation,N>:: 95 Stepper(const G4double yInput[], 96 const G4double dydx[], 97 G4double hstep, 98 G4double yOutput[], 99 G4double yError[]) 100 // The stepper for the Runge Kutta integration. The stepsize 101 // is fixed, with the Step size given by hstep. 102 // Integrates ODE starting values y[0 to N]. 103 // Outputs yout[] and its estimated error yerr[]. 104 { 105 const unsigned int maxvar = GetNumberOfStateVariables(); 106 107 // Saving yInput because yInput and yOutput can be aliases for same array 108 for(unsigned int i = 0; i < N; ++i) 109 yInitial[i] = yInput[i]; 110 yInitial[7] = 111 yInput[7]; // Copy the time in case ... even if not really needed 112 yMiddle[7] = yInput[7]; // Copy the time from initial value 113 yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ? 114 // yOutput[7] = yInput[7]; // -> dumb stepper does it too for RK4 115 for(unsigned int i = N; i < maxvar; ++i) 116 yOutput[i] = yInput[i]; 117 118 G4double halfStep = hstep * 0.5; 119 120 // Do two half steps 121 122 static_cast<T_Stepper*>(this)->DumbStepper(yInitial, dydx, halfStep, 123 yMiddle); 124 this->RightHandSide(yMiddle, dydxMid); 125 static_cast<T_Stepper*>(this)->DumbStepper(yMiddle, dydxMid, halfStep, 126 yOutput); 127 128 // Store midpoint, chord calculation 129 130 fMidPoint = G4ThreeVector(yMiddle[0], yMiddle[1], yMiddle[2]); 131 132 // Do a full Step 133 static_cast<T_Stepper*>(this)->DumbStepper(yInitial, dydx, hstep, yOneStep); 134 for(unsigned int i = 0; i < N; ++i) 135 { 136 yError[i] = yOutput[i] - yOneStep[i]; 137 yOutput[i] += 138 yError[i] * 139 T_Stepper::IntegratorCorrection; // Provides accuracy increased 140 // by 1 order via the 141 // Richardson Extrapolation 142 } 143 144 fInitialPoint = G4ThreeVector(yInitial[0], yInitial[1], yInitial[2]); 145 fFinalPoint = G4ThreeVector(yOutput[0], yOutput[1], yOutput[2]); 146 147 return; 148 } 149 150 151 template <class T_Stepper, class T_Equation, unsigned int N > 152 inline G4double 153 G4TMagErrorStepper<T_Stepper,T_Equation,N>::DistChord() const 154 { 155 // Estimate the maximum distance from the curve to the chord 156 // 157 // We estimate this using the distance of the midpoint to 158 // chord (the line between 159 // 160 // Method below is good only for angle deviations < 2 pi, 161 // This restriction should not a problem for the Runge cutta methods, 162 // which generally cannot integrate accurately for large angle deviations. 163 G4double distLine, distChord; 164 165 if(fInitialPoint != fFinalPoint) 166 { 167 distLine = G4LineSection::Distline(fMidPoint, fInitialPoint, fFinalPoint); 168 // This is a class method that gives distance of Mid 169 // from the Chord between the Initial and Final points. 170 171 distChord = distLine; 172 } 173 else 174 { 175 distChord = (fMidPoint - fInitialPoint).mag(); 176 } 177 178 return distChord; 179 } 180 181 #endif /* G4TMagErrorStepper_HH */ 182