Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HelixMixedStepper Class Reference

#include <G4HelixMixedStepper.hh>

+ Inheritance diagram for G4HelixMixedStepper:

Public Member Functions

 G4HelixMixedStepper (G4Mag_EqRhs *EqRhs, G4int fStepperNumber=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)
 
G4int IntegratorOrder () const
 
- Public Member Functions inherited from G4MagHelicalStepper
 G4MagHelicalStepper (G4Mag_EqRhs *EqRhs)
 
virtual ~G4MagHelicalStepper ()
 
virtual void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
 
virtual void DumbStepper (const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])=0
 
G4double DistChord () const
 
- Public Member Functions inherited from G4MagIntegratorStepper
 G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12)
 
virtual ~G4MagIntegratorStepper ()
 
virtual void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
 
virtual G4double DistChord () const =0
 
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
 
virtual G4int IntegratorOrder () const =0
 
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::G4HelixMixedStepper ( G4Mag_EqRhs EqRhs,
G4int  fStepperNumber = 0 
)

Definition at line 57 of file G4HelixMixedStepper.cc.

58 : G4MagHelicalStepper(EqRhs)
59
60{
61 SetVerbose(1); fNumCallsRK4=0; fNumCallsHelix=0;
62 if(!fStepperNumber) fStepperNumber=4;
63 fRK4Stepper = SetupStepper(EqRhs, fStepperNumber);
64}
void SetVerbose(G4int newvalue)
G4MagIntegratorStepper * SetupStepper(G4Mag_EqRhs *EqRhs, G4int StepperName)

◆ ~G4HelixMixedStepper()

G4HelixMixedStepper::~G4HelixMixedStepper ( )

Definition at line 67 of file G4HelixMixedStepper.cc.

67 {
68
69 delete(fRK4Stepper);
70 if (fVerbose>0){ PrintCalls();};
71}

Member Function Documentation

◆ DistChord()

G4double G4HelixMixedStepper::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 145 of file G4HelixMixedStepper.cc.

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}
double G4double
Definition: G4Types.hh:64
G4double GetRadHelix() const
G4double GetAngCurve() const

◆ DumbStepper()

void G4HelixMixedStepper::DumbStepper ( const G4double  y[],
G4ThreeVector  Bfld,
G4double  h,
G4double  yout[] 
)
virtual

Implements G4MagHelicalStepper.

Definition at line 132 of file G4HelixMixedStepper.cc.

136{
137
138
139 AdvanceHelix(yIn, Bfld, h, yOut);
140
141
142
143}
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)

◆ IntegratorOrder()

G4int G4HelixMixedStepper::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 93 of file G4HelixMixedStepper.hh.

93{ return 4; }

◆ PrintCalls()

void G4HelixMixedStepper::PrintCalls ( )

Definition at line 172 of file G4HelixMixedStepper.cc.

173{
174 G4cout<<"In HelixMixedStepper::Number of calls to smallStepStepper = "<<fNumCallsRK4
175 <<" and Number of calls to Helix = "<<fNumCallsHelix<<G4endl;
176}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

Referenced by ~G4HelixMixedStepper().

◆ SetupStepper()

G4MagIntegratorStepper * G4HelixMixedStepper::SetupStepper ( G4Mag_EqRhs EqRhs,
G4int  StepperName 
)

Definition at line 180 of file G4HelixMixedStepper.cc.

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}

Referenced by G4HelixMixedStepper().

◆ SetVerbose()

void G4HelixMixedStepper::SetVerbose ( G4int  newvalue)
inline

Definition at line 86 of file G4HelixMixedStepper.hh.

86{fVerbose=newvalue;}

Referenced by G4HelixMixedStepper().

◆ Stepper()

void G4HelixMixedStepper::Stepper ( const G4double  y[],
const G4double  dydx[],
G4double  h,
G4double  yout[],
G4double  yerr[] 
)
virtual

Reimplemented from G4MagHelicalStepper.

Definition at line 72 of file G4HelixMixedStepper.cc.

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

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