Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
55 #include "G4ThreeVector.hh"
56 #include "G4LineSection.hh"
58  : G4MagHelicalStepper(EqRhs)
59 
60 {
61  SetVerbose(1); fNumCallsRK4=0; fNumCallsHelix=0;
62  if(!fStepperNumber) fStepperNumber=4;
63  fRK4Stepper = SetupStepper(EqRhs, fStepperNumber);
64 }
65 
66 
68 
69  delete(fRK4Stepper);
70  if (fVerbose>0){ PrintCalls();};
71 }
72 void G4HelixMixedStepper::Stepper( const G4double yInput[7],
73  const G4double dydx[7],
74  G4double Step,
75  G4double yOut[7],
76  G4double yErr[])
77 
78 {
79 
80  //Estimation of the Stepping Angle
81 
82  G4ThreeVector Bfld;
83  MagFieldEvaluate(yInput, Bfld);
84 
85  G4double Bmag = Bfld.mag();
86  const G4double *pIn = yInput+3;
87  G4ThreeVector initVelocity= G4ThreeVector( pIn[0], pIn[1], pIn[2]);
88  G4double velocityVal = initVelocity.mag();
89  G4double R_1;
90  G4double Ang_curve;
91 
92  R_1=std::abs(GetInverseCurve(velocityVal,Bmag));
93  Ang_curve=R_1*Step;
94  SetAngCurve(Ang_curve);
95  SetCurve(std::abs(1/R_1));
96 
97 
98  if(Ang_curve<0.33*pi){
99  fNumCallsRK4++;
100  fRK4Stepper->Stepper(yInput,dydx,Step,yOut,yErr);
101 
102 
103  }
104  else{
105  fNumCallsHelix++;
106  const G4int nvar = 6 ;
107  G4int i;
108  G4double yTemp[7], yIn[7] ;
109  G4double yTemp2[7];
110  G4ThreeVector Bfld_midpoint;
111  // Saving yInput because yInput and yOut can be aliases for same array
112  for(i=0;i<nvar;i++) yIn[i]=yInput[i];
113 
114  G4double h = Step * 0.5;
115  // Do two half steps and full step
116  AdvanceHelix(yIn, Bfld, h, yTemp,yTemp2);
117  MagFieldEvaluate(yTemp, Bfld_midpoint) ;
118  AdvanceHelix(yTemp, Bfld_midpoint, h, yOut);
119  // Error estimation
120  for(i=0;i<nvar;i++) {
121  yErr[i] = yOut[i] - yTemp2[i] ;
122 
123  }
124  }
125 
126 
127 
128 
129 }
130 
131 void
133  G4ThreeVector Bfld,
134  G4double h,
135  G4double yOut[])
136 {
137 
138 
139  AdvanceHelix(yIn, Bfld, h, yOut);
140 
141 
142 
143 }
144 
146 {
147  // Implementation : must check whether h/R > 2 pi !!
148  // If( h/R < pi) use G4LineSection::DistLine
149  // Else DistChord=R_helix
150  //
151  G4double distChord;
152  G4double Ang_curve=GetAngCurve();
153 
154 
155  if(Ang_curve<=pi){
156  distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
157  }
158  else
159  if(Ang_curve<twopi){
160  distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
161  }
162  else{
163  distChord=2.*GetRadHelix();
164  }
165 
166 
167 
168  return distChord;
169 
170 }
171 // ---------------------------------------------------------------------------
173 {
174  G4cout<<"In HelixMixedStepper::Number of calls to smallStepStepper = "<<fNumCallsRK4
175  <<" and Number of calls to Helix = "<<fNumCallsHelix<<G4endl;
176 }
177 
178 
179 
181 {
182  G4MagIntegratorStepper* pStepper;
183  if (fVerbose>0)G4cout<<"In G4HelixMixedStepper Stepper for small steps is ";
184  switch ( StepperNumber )
185  {
186  case 0: pStepper = new G4ExplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4ExplicitEuler"<<G4endl; break;
187  case 1: pStepper = new G4ImplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4ImplicitEuler"<<G4endl; break;
188  case 2: pStepper = new G4SimpleRunge( pE ); if (fVerbose>0)G4cout<<"G4SimpleRunge"<<G4endl; break;
189  case 3: pStepper = new G4SimpleHeum( pE ); if (fVerbose>0)G4cout<<"G4SimpleHeum"<<G4endl;break;
190  case 4: pStepper = new G4ClassicalRK4( pE ); if (fVerbose>0)G4cout<<"G4ClassicalRK4"<<G4endl; break;
191  case 5: pStepper = new G4HelixExplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4HelixExplicitEuler"<<G4endl; break;
192  case 6: pStepper = new G4HelixImplicitEuler( pE ); if (fVerbose>0)G4cout<<"G4HelixImplicitEuler"<<G4endl; break;
193  case 7: pStepper = new G4HelixSimpleRunge( pE ); if (fVerbose>0)G4cout<<"G4HelixSimpleRunge"<<G4endl; break;
194  case 8: pStepper = new G4CashKarpRKF45( pE ); if (fVerbose>0)G4cout<<"G4CashKarpRKF45"<<G4endl; break;
195  case 9: pStepper = new G4ExactHelixStepper( pE ); if (fVerbose>0)G4cout<<"G4ExactHelixStepper"<<G4endl; break;
196  case 10: pStepper = new G4RKG3_Stepper( pE ); if (fVerbose>0)G4cout<<"G4RKG3_Stepper"<<G4endl; break;
197 
198  default: pStepper = new G4ClassicalRK4( pE );G4cout<<"Default G4ClassicalRK4"<<G4endl; break;
199 
200  }
201  return pStepper;
202 }