Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4HelixMixedStepper Class Reference

#include <G4HelixMixedStepper.hh>

Inheritance diagram for G4HelixMixedStepper:
G4MagHelicalStepper G4MagIntegratorStepper

Public Member Functions

 G4HelixMixedStepper (G4Mag_EqRhs *EqRhs, G4int fStepperNumber=-1, G4double Angle_threshold=-1.0)
 
 ~G4HelixMixedStepper ()
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
 
void DumbStepper (const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])
 
G4double DistChord () const
 
void SetVerbose (G4int newvalue)
 
void PrintCalls ()
 
G4MagIntegratorStepperSetupStepper (G4Mag_EqRhs *EqRhs, G4int StepperName)
 
void SetAngleThreshold (G4double val)
 
G4double GetAngleThreshold ()
 
G4int IntegratorOrder () const
 
- Public Member Functions inherited from G4MagHelicalStepper
 G4MagHelicalStepper (G4Mag_EqRhs *EqRhs)
 
virtual ~G4MagHelicalStepper ()
 
G4double DistChord () const
 
- Public Member Functions inherited from G4MagIntegratorStepper
 G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12)
 
virtual ~G4MagIntegratorStepper ()
 
virtual void ComputeRightHandSide (const G4double y[], G4double dydx[])
 
void NormaliseTangentVector (G4double vec[6])
 
void NormalisePolarizationVector (G4double vec[12])
 
void RightHandSide (const double y[], double dydx[])
 
G4int GetNumberOfVariables () const
 
G4int GetNumberOfStateVariables () const
 
G4EquationOfMotionGetEquationOfMotion ()
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 

Additional Inherited Members

- Protected Member Functions inherited from G4MagHelicalStepper
void LinearStep (const G4double yIn[], G4double h, G4double yHelix[]) const
 
void AdvanceHelix (const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
 
void MagFieldEvaluate (const G4double y[], G4ThreeVector &Bfield)
 
G4double GetInverseCurve (const G4double Momentum, const G4double Bmag)
 
void SetAngCurve (const G4double Ang)
 
G4double GetAngCurve () const
 
void SetCurve (const G4double Curve)
 
G4double GetCurve () const
 
void SetRadHelix (const G4double Rad)
 
G4double GetRadHelix () const
 

Detailed Description

Definition at line 58 of file G4HelixMixedStepper.hh.

Constructor & Destructor Documentation

G4HelixMixedStepper::G4HelixMixedStepper ( G4Mag_EqRhs EqRhs,
G4int  fStepperNumber = -1,
G4double  Angle_threshold = -1.0 
)

Definition at line 59 of file G4HelixMixedStepper.cc.

References python.hepunit::pi, SetupStepper(), and SetVerbose().

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 }
G4MagIntegratorStepper * SetupStepper(G4Mag_EqRhs *EqRhs, G4int StepperName)
G4MagHelicalStepper(G4Mag_EqRhs *EqRhs)
void SetVerbose(G4int newvalue)
G4HelixMixedStepper::~G4HelixMixedStepper ( )

Definition at line 74 of file G4HelixMixedStepper.cc.

References PrintCalls().

75 {
76  delete(fRK4Stepper);
77  if (fVerbose>0){ PrintCalls();};
78 }

Member Function Documentation

G4double G4HelixMixedStepper::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 147 of file G4HelixMixedStepper.cc.

References G4MagHelicalStepper::GetAngCurve(), G4MagHelicalStepper::GetRadHelix(), python.hepunit::pi, and python.hepunit::twopi.

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 }
G4double GetRadHelix() const
double G4double
Definition: G4Types.hh:76
G4double GetAngCurve() const
void G4HelixMixedStepper::DumbStepper ( const G4double  y[],
G4ThreeVector  Bfld,
G4double  h,
G4double  yout[] 
)
virtual

Implements G4MagHelicalStepper.

Definition at line 139 of file G4HelixMixedStepper.cc.

References G4MagHelicalStepper::AdvanceHelix().

143 {
144  AdvanceHelix(yIn, Bfld, h, yOut);
145 }
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
G4double G4HelixMixedStepper::GetAngleThreshold ( )
inline

Definition at line 93 of file G4HelixMixedStepper.hh.

93 { return fAngle_threshold; }
G4int G4HelixMixedStepper::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 95 of file G4HelixMixedStepper.hh.

95 { return 4; }
void G4HelixMixedStepper::PrintCalls ( )

Definition at line 173 of file G4HelixMixedStepper.cc.

References G4cout, and G4endl.

Referenced by ~G4HelixMixedStepper().

174 {
175  G4cout<<"In HelixMixedStepper::Number of calls to smallStepStepper = "<<fNumCallsRK4
176  <<" and Number of calls to Helix = "<<fNumCallsHelix<<G4endl;
177 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4HelixMixedStepper::SetAngleThreshold ( G4double  val)
inline

Definition at line 92 of file G4HelixMixedStepper.hh.

92 { fAngle_threshold= val;}
G4MagIntegratorStepper * G4HelixMixedStepper::SetupStepper ( G4Mag_EqRhs EqRhs,
G4int  StepperName 
)

Definition at line 180 of file G4HelixMixedStepper.cc.

References G4cout, and G4endl.

Referenced by G4HelixMixedStepper().

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 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4HelixMixedStepper::SetVerbose ( G4int  newvalue)
inline

Definition at line 86 of file G4HelixMixedStepper.hh.

Referenced by G4HelixMixedStepper().

86 {fVerbose=newvalue;}
void G4HelixMixedStepper::Stepper ( const G4double  y[],
const G4double  dydx[],
G4double  h,
G4double  yout[],
G4double  yerr[] 
)
virtual

Reimplemented from G4MagHelicalStepper.

Definition at line 80 of file G4HelixMixedStepper.cc.

References G4MagHelicalStepper::AdvanceHelix(), G4MagHelicalStepper::GetInverseCurve(), CLHEP::Hep3Vector::mag(), G4MagHelicalStepper::MagFieldEvaluate(), G4MagHelicalStepper::SetAngCurve(), G4MagHelicalStepper::SetCurve(), and G4MagIntegratorStepper::Stepper().

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 }
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
CLHEP::Hep3Vector G4ThreeVector
virtual void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
int G4int
Definition: G4Types.hh:78
G4double GetInverseCurve(const G4double Momentum, const G4double Bmag)
void SetAngCurve(const G4double Ang)
void SetCurve(const G4double Curve)
Definition: Step.hh:41
double G4double
Definition: G4Types.hh:76
double mag() const

The documentation for this class was generated from the following files: