// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // class G4HelixMixedStepper // // Class description: // // G4HelixMixedStepper split the Method used for Integration in two: // // If Stepping Angle ( h / R_curve) < pi/3 // use Stepper for small step(ClassicalRK4 by default) // Else use HelixExplicitEuler Stepper // // History: // Derived from ExactHelicalStepper 18/05/07 // // ------------------------------------------------------------------------- #include "G4HelixMixedStepper.hh" #include "G4ClassicalRK4.hh" #include "G4CashKarpRKF45.hh" #include "G4SimpleRunge.hh" #include "G4HelixImplicitEuler.hh" #include "G4HelixExplicitEuler.hh" #include "G4HelixSimpleRunge.hh" #include "G4ExactHelixStepper.hh" #include "G4ExplicitEuler.hh" #include "G4ImplicitEuler.hh" #include "G4SimpleHeum.hh" #include "G4RKG3_Stepper.hh" #include "G4ThreeVector.hh" #include "G4LineSection.hh" G4HelixMixedStepper::G4HelixMixedStepper(G4Mag_EqRhs *EqRhs,G4int fStepperNumber) : G4MagHelicalStepper(EqRhs) { SetVerbose(1); fNumCallsRK4=0; fNumCallsHelix=0; if(!fStepperNumber) fStepperNumber=4; fRK4Stepper = SetupStepper(EqRhs, fStepperNumber); } G4HelixMixedStepper::~G4HelixMixedStepper() { delete(fRK4Stepper); if (fVerbose>0){ PrintCalls();}; } void G4HelixMixedStepper::Stepper( const G4double yInput[7], const G4double dydx[7], G4double Step, G4double yOut[7], G4double yErr[]) { //Estimation of the Stepping Angle G4ThreeVector Bfld; MagFieldEvaluate(yInput, Bfld); G4double Bmag = Bfld.mag(); const G4double *pIn = yInput+3; G4ThreeVector initVelocity= G4ThreeVector( pIn[0], pIn[1], pIn[2]); G4double velocityVal = initVelocity.mag(); G4double R_1; G4double Ang_curve; R_1=std::abs(GetInverseCurve(velocityVal,Bmag)); Ang_curve=R_1*Step; SetAngCurve(Ang_curve); SetCurve(std::abs(1/R_1)); if(Ang_curve<0.33*pi){ fNumCallsRK4++; fRK4Stepper->Stepper(yInput,dydx,Step,yOut,yErr); } else{ fNumCallsHelix++; const G4int nvar = 6 ; G4int i; G4double yTemp[7], yIn[7] ; G4double yTemp2[7]; G4ThreeVector Bfld_midpoint; // Saving yInput because yInput and yOut can be aliases for same array for(i=0;i 2 pi !! // If( h/R < pi) use G4LineSection::DistLine // Else DistChord=R_helix // G4double distChord; G4double Ang_curve=GetAngCurve(); if(Ang_curve<=pi){ distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve)); } else if(Ang_curve0)G4cout<<"G4ExplicitEuler"<0)G4cout<<"G4ImplicitEuler"<0)G4cout<<"G4SimpleRunge"<0)G4cout<<"G4SimpleHeum"<0)G4cout<<"G4ClassicalRK4"<0)G4cout<<"G4HelixExplicitEuler"<0)G4cout<<"G4HelixImplicitEuler"<0)G4cout<<"G4HelixSimpleRunge"<0)G4cout<<"G4CashKarpRKF45"<0)G4cout<<"G4ExactHelixStepper"<0)G4cout<<"G4RKG3_Stepper"<