Geant4  10.00.p01
G4HelixMixedStepper.cc
Go to the documentation of this file.
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 // class G4HelixMixedStepper
27 //
28 // Class description:
29 //
30 // G4HelixMixedStepper split the Method used for Integration in two:
31 //
32 // If Stepping Angle ( h / R_curve) < pi/3
33 // use Stepper for small step(ClassicalRK4 by default)
34 // Else use HelixExplicitEuler Stepper
35 //
36 // History:
37 // Derived from ExactHelicalStepper 18/05/07
38 //
39 // -------------------------------------------------------------------------
40 
41 #include "G4HelixMixedStepper.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4ClassicalRK4.hh"
44 #include "G4CashKarpRKF45.hh"
45 #include "G4SimpleRunge.hh"
46 #include "G4HelixImplicitEuler.hh"
47 #include "G4HelixExplicitEuler.hh"
48 #include "G4HelixSimpleRunge.hh"
49 #include "G4ExactHelixStepper.hh"
50 #include "G4ExplicitEuler.hh"
51 #include "G4ImplicitEuler.hh"
52 #include "G4SimpleHeum.hh"
53 #include "G4RKG3_Stepper.hh"
54 #include "G4NystromRK4.hh"
55 
56 #include "G4ThreeVector.hh"
57 #include "G4LineSection.hh"
58 
60  : G4MagHelicalStepper(EqRhs), fNumCallsRK4(0), fNumCallsHelix(0)
61 {
62  SetVerbose(1);
63  if( angleThr < 0.0 ){
64  fAngle_threshold= 0.33*pi;
65  }else{
66  fAngle_threshold= angleThr;
67  }
68 
69  if(fStepperNumber<0)
70  fStepperNumber=4; // Default is RK4
71  fRK4Stepper = SetupStepper(EqRhs, fStepperNumber);
72 }
73 
75 {
76  delete(fRK4Stepper);
77  if (fVerbose>0){ PrintCalls();};
78 }
79 
80 void G4HelixMixedStepper::Stepper( const G4double yInput[7],
81  const G4double dydx[7],
82  G4double Step,
83  G4double yOut[7],
84  G4double yErr[])
85 
86 {
87  //Estimation of the Stepping Angle
88 
89  G4ThreeVector Bfld;
90  MagFieldEvaluate(yInput, Bfld);
91 
92  G4double Bmag = Bfld.mag();
93  const G4double *pIn = yInput+3;
94  G4ThreeVector initVelocity= G4ThreeVector( pIn[0], pIn[1], pIn[2]);
95  G4double velocityVal = initVelocity.mag();
96  G4double R_1;
97  G4double Ang_curve;
98 
99  R_1=std::abs(GetInverseCurve(velocityVal,Bmag));
100  Ang_curve=R_1*Step;
101  SetAngCurve(Ang_curve);
102  SetCurve(std::abs(1/R_1));
103 
104  if(Ang_curve< fAngle_threshold){
105  fNumCallsRK4++;
106  fRK4Stepper->Stepper(yInput,dydx,Step,yOut,yErr);
107  }
108  else
109  {
110  fNumCallsHelix++;
111  const G4int nvar = 6 ;
112  G4int i;
113  G4double yTemp[7], yIn[7], yTemp2[7];
114  G4ThreeVector Bfld_midpoint;
115 
116  // Saving yInput because yInput and yOut can be aliases for same array
117  for(i=0;i<nvar;i++) yIn[i]=yInput[i];
118 
119  G4double halfS = Step * 0.5;
120  // 1. Do first half step and full step
121  AdvanceHelix(yIn, Bfld, halfS, yTemp, yTemp2); // yTemp2 for s=2*h (halfS)
122  //**********
123 
124  MagFieldEvaluate(yTemp, Bfld_midpoint) ;
125 
126  // 2. Do second half step - with revised field
127  // NOTE: Could avoid this call if 'Bfld_midpoint == Bfld' or diff 'almost' zero
128  AdvanceHelix(yTemp, Bfld_midpoint, halfS, yOut); // Not requesting y at s=2*h (halfS)
129  //**********
130 
131  // 3. Estimate the integration error - should be (nearly) zero if Bfield= constant
132  for(i=0;i<nvar;i++) {
133  yErr[i] = yOut[i] - yTemp2[i] ;
134  }
135  }
136 }
137 
138 void
140  G4ThreeVector Bfld,
141  G4double h,
142  G4double yOut[])
143 {
144  AdvanceHelix(yIn, Bfld, h, yOut);
145 }
146 
148 {
149  // Implementation : must check whether h/R > 2 pi !!
150  // If( h/R < pi) use G4LineSection::DistLine
151  // Else DistChord=R_helix
152  //
153  G4double distChord;
154  G4double Ang_curve=GetAngCurve();
155 
156  if(Ang_curve<=pi){
157  distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
158  }
159  else
160  {
161  if(Ang_curve<twopi){
162  distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
163  }
164  else{
165  distChord=2.*GetRadHelix();
166  }
167  }
168 
169  return distChord;
170 }
171 
172 // ---------------------------------------------------------------------------
174 {
175  G4cout<<"In HelixMixedStepper::Number of calls to smallStepStepper = "<<fNumCallsRK4
176  <<" and Number of calls to Helix = "<<fNumCallsHelix<<G4endl;
177 }
178 
181 {
182  G4MagIntegratorStepper* pStepper;
183  if (fVerbose>0)
184  G4cout<<"In G4HelixMixedStepper: Choosing Stepper for small steps. Choice is ";
185  switch ( StepperNumber )
186  {
187  case 0:
188  case 4: pStepper = new G4ClassicalRK4( pE ); if (fVerbose>0)G4cout<<"G4ClassicalRK4"<<G4endl; break;
189 
190  // Steppers with embedded estimation of error
191  case 1: pStepper = new G4NystromRK4( pE ); if (fVerbose>0)G4cout<<"G4NystromRK4"<<G4endl; break;
192  case 8: pStepper = new G4CashKarpRKF45( pE ); if (fVerbose>0)G4cout<<"G4CashKarpRKF45"<<G4endl; break;
193 
194  // Lower order RK Steppers - ok overall, good for uneven fields
195  case 2: pStepper = new G4SimpleRunge( pE ); if (fVerbose>0)G4cout<<"G4SimpleRunge"<<G4endl; break;
196  case 3: pStepper = new G4SimpleHeum( pE ); if (fVerbose>0)G4cout<<"G4SimpleHeum"<<G4endl;break;
197 
198  // Helical Steppers -
199  // Since these are used for long steps, less useful
200  case 5: pStepper = new G4HelixExplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4HelixExplicitEuler"<<G4endl; break;
201  case 6: pStepper = new G4HelixImplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4HelixImplicitEuler"<<G4endl; break;
202  case 7: pStepper = new G4HelixSimpleRunge( pE ); if (fVerbose>0)G4cout<<"G4HelixSimpleRunge"<<G4endl; break;
203 
204  // Exact Helix - only for cases of i)uniform field
205  // ii) segmented uniform field (maybe)
206  case 9: pStepper = new G4ExactHelixStepper( pE ); if (fVerbose>0)G4cout<<"G4ExactHelixStepper"<<G4endl; break;
207 
208  case 10: pStepper = new G4RKG3_Stepper( pE ); if (fVerbose>0)G4cout<<"G4RKG3_Stepper"<<G4endl; break;
209 
210  // Low Order Steppers - not good except for very weak fields
211  case 11: pStepper = new G4ExplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4ExplicitEuler"<<G4endl; break;
212  case 12: pStepper = new G4ImplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4ImplicitEuler"<<G4endl; break;
213 
214  case -1:
215  default: pStepper = new G4ClassicalRK4( pE );G4cout<<"G4ClassicalRK4 (Default)."<<G4endl; break;
216 
217  }
218  return pStepper;
219 }
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
CLHEP::Hep3Vector G4ThreeVector
G4double GetRadHelix() const
virtual void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
const G4double pi
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
G4MagIntegratorStepper * SetupStepper(G4Mag_EqRhs *EqRhs, G4int StepperName)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4HelixMixedStepper(G4Mag_EqRhs *EqRhs, G4int fStepperNumber=-1, G4double Angle_threshold=-1.0)
G4double GetInverseCurve(const G4double Momentum, const G4double Bmag)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
void SetAngCurve(const G4double Ang)
G4MagIntegratorStepper * fRK4Stepper
void SetCurve(const G4double Curve)
void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])
Definition: Step.hh:41
#define G4endl
Definition: G4ios.hh:61
void SetVerbose(G4int newvalue)
double G4double
Definition: G4Types.hh:76
G4double GetAngCurve() const
G4double DistChord() const