Geant4 Cross Reference |
1 // 1 2 // ******************************************* 3 // * License and Disclaimer 4 // * 5 // * The Geant4 software is copyright of th 6 // * the Geant4 Collaboration. It is provided 7 // * conditions of the Geant4 Software License 8 // * LICENSE and available at http://cern.ch/ 9 // * include a list of copyright holders. 10 // * 11 // * Neither the authors of this software syst 12 // * institutes,nor the agencies providing fin 13 // * work make any representation or warran 14 // * regarding this software system or assum 15 // * use. Please see the license in the file 16 // * for the full disclaimer and the limitatio 17 // * 18 // * This code implementation is the result 19 // * technical work of the GEANT4 collaboratio 20 // * By using, copying, modifying or distri 21 // * any work based on the software) you ag 22 // * use in resulting scientific publicati 23 // * acceptance of all terms of the Geant4 Sof 24 // ******************************************* 25 // 26 // G4TMagErrorStepper 27 // 28 // Class description: 29 // 30 // Templated version of G4MagErrorStepper 31 // 32 // 33 // Created: Josh Xie (supported by Google Sum 34 // Supervisors: Sandro Wenzel, John Apostolak 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, u 46 class G4TMagErrorStepper : public G4MagIntegra 47 { 48 public: // with description 49 G4TMagErrorStepper(T_Equation* EqRhs, G4int 50 G4int numStateVariables = 51 : G4MagIntegratorStepper(EqRhs, numberOfVa 52 , fEquation_Rhs(EqRhs) 53 { 54 // G4int nvar = std::max(this->GetNumberOf 55 } 56 57 virtual ~G4TMagErrorStepper() { ; } 58 59 inline void RightHandSide(G4double y[], G4do 60 { 61 fEquation_Rhs->T_Equation::RightHandSide(y 62 } 63 64 inline void Stepper(const G4double yInput[], 65 G4double hstep, G4double 66 67 inline G4double DistChord() const override f 68 69 private: 70 G4TMagErrorStepper(const G4TMagErrorStepper& 71 G4TMagErrorStepper& operator=(const G4TMagEr 72 // Private copy constructor and assignment o 73 74 private: 75 // STATE 76 G4ThreeVector fInitialPoint, fMidPoint, fFin 77 // Data stored in order to find the chord 78 79 // Dependent Objects, owned --- part of the 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 te 85 // they are allocated at the class level onl 86 // so that calls to new and delete are not m 87 88 T_Equation* fEquation_Rhs; 89 }; 90 91 // ------------ Implementation ------------- 92 93 template <class T_Stepper, class T_Equation, u 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 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 = GetNumberOfState 106 107 // Saving yInput because yInput and yOutput 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 ... 112 yMiddle[7] = yInput[7]; // Copy the time f 113 yOneStep[7] = yInput[7]; // As it contribut 114 // yOutput[7] = yInput[7]; // -> dumb stepp 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(y 123 124 this->RightHandSide(yMiddle, dydxMid); 125 static_cast<T_Stepper*>(this)->DumbStepper(y 126 y 127 128 // Store midpoint, chord calculation 129 130 fMidPoint = G4ThreeVector(yMiddle[0], yMiddl 131 132 // Do a full Step 133 static_cast<T_Stepper*>(this)->DumbStepper(y 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; // Pr 140 // by 1 order via the 141 // Richardson Extrapolation 142 } 143 144 fInitialPoint = G4ThreeVector(yInitial[0], y 145 fFinalPoint = G4ThreeVector(yOutput[0], yO 146 147 return; 148 } 149 150 151 template <class T_Stepper, class T_Equation, u 152 inline G4double 153 G4TMagErrorStepper<T_Stepper,T_Equation,N>::Di 154 { 155 // Estimate the maximum distance from the cu 156 // 157 // We estimate this using the distance of t 158 // chord (the line between 159 // 160 // Method below is good only for angle devi 161 // This restriction should not a problem f 162 // which generally cannot integrate accura 163 G4double distLine, distChord; 164 165 if(fInitialPoint != fFinalPoint) 166 { 167 distLine = G4LineSection::Distline(fMidPoi 168 // This is a class method that gives dista 169 // from the Chord between the Initial and 170 171 distChord = distLine; 172 } 173 else 174 { 175 distChord = (fMidPoint - fInitialPoint).m 176 } 177 178 return distChord; 179 } 180 181 #endif /* G4TMagErrorStepper_HH */ 182