Geant4 11.2.2
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 StepperNumber=-1, G4double Angle_threshold=-1.0)
 
 ~G4HelixMixedStepper () override
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
 
void DumbStepper (const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[]) override
 
G4double DistChord () const override
 
void SetVerbose (G4int newvalue)
 
void PrintCalls ()
 
G4MagIntegratorStepperSetupStepper (G4Mag_EqRhs *EqRhs, G4int StepperName)
 
void SetAngleThreshold (G4double val)
 
G4double GetAngleThreshold ()
 
G4int IntegratorOrder () const override
 
- Public Member Functions inherited from G4MagHelicalStepper
 G4MagHelicalStepper (G4Mag_EqRhs *EqRhs)
 
 ~G4MagHelicalStepper () override
 
 G4MagHelicalStepper (const G4MagHelicalStepper &)=delete
 
G4MagHelicalStepperoperator= (const G4MagHelicalStepper &)=delete
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
 
G4double DistChord () const override
 
- Public Member Functions inherited from G4MagIntegratorStepper
 G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
 
virtual ~G4MagIntegratorStepper ()=default
 
 G4MagIntegratorStepper (const G4MagIntegratorStepper &)=delete
 
G4MagIntegratorStepperoperator= (const G4MagIntegratorStepper &)=delete
 
void NormaliseTangentVector (G4double vec[6])
 
void NormalisePolarizationVector (G4double vec[12])
 
void RightHandSide (const G4double y[], G4double dydx[]) const
 
void RightHandSide (const G4double y[], G4double dydx[], G4double field[]) const
 
G4int GetNumberOfVariables () const
 
G4int GetNumberOfStateVariables () const
 
G4int IntegrationOrder ()
 
G4EquationOfMotionGetEquationOfMotion ()
 
const G4EquationOfMotionGetEquationOfMotion () const
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 
unsigned long GetfNoRHSCalls ()
 
void ResetfNORHSCalls ()
 
G4bool IsFSAL () const
 
G4bool isQSS () const
 
void SetIsQSS (G4bool val)
 

Additional Inherited Members

- Protected Member Functions inherited from G4MagHelicalStepper
void LinearStep (const G4double yIn[], G4double h, G4double yHelix[]) const
 
void AdvanceHelix (const G4double yIn[], const G4ThreeVector &Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=nullptr)
 
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
 
- Protected Member Functions inherited from G4MagIntegratorStepper
void SetIntegrationOrder (G4int order)
 
void SetFSAL (G4bool flag=true)
 

Detailed Description

Definition at line 62 of file G4HelixMixedStepper.hh.

Constructor & Destructor Documentation

◆ G4HelixMixedStepper()

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

Definition at line 64 of file G4HelixMixedStepper.cc.

68 : G4MagHelicalStepper(EqRhs)
69{
70 if( angleThreshold < 0.0 )
71 {
72 fAngle_threshold = (1.0/3.0)*pi;
73 }
74 else
75 {
76 fAngle_threshold = angleThreshold;
77 }
78
79 if(stepperNumber<0)
80 {
81 // stepperNumber = 4; // Default is RK4 (original)
82 stepperNumber = 745; // Default is DormandPrince745 (ie DoPri5)
83 // stepperNumber = 8; // Default is CashKarp
84 }
85
86 fStepperNumber = stepperNumber; // Store the choice
87 fRK4Stepper = SetupStepper(EqRhs, fStepperNumber);
88}
G4MagIntegratorStepper * SetupStepper(G4Mag_EqRhs *EqRhs, G4int StepperName)
G4MagHelicalStepper(G4Mag_EqRhs *EqRhs)

◆ ~G4HelixMixedStepper()

G4HelixMixedStepper::~G4HelixMixedStepper ( )
override

Definition at line 91 of file G4HelixMixedStepper.cc.

92{
93 delete fRK4Stepper;
94 if (fVerbose>0) { PrintCalls(); }
95}

Member Function Documentation

◆ DistChord()

G4double G4HelixMixedStepper::DistChord ( ) const
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 178 of file G4HelixMixedStepper.cc.

179{
180 // Implementation : must check whether h/R > 2 pi !!
181 // If( h/R < pi) use G4LineSection::DistLine
182 // Else DistChord=R_helix
183 //
184 G4double distChord;
185 G4double Ang_curve=GetAngCurve();
186
187 if(Ang_curve<=pi)
188 {
189 distChord=GetRadHelix()*(1-std::cos(0.5*Ang_curve));
190 }
191 else
192 {
193 if(Ang_curve<twopi)
194 {
195 distChord=GetRadHelix()*(1+std::cos(0.5*(twopi-Ang_curve)));
196 }
197 else
198 {
199 distChord=2.*GetRadHelix();
200 }
201 }
202
203 return distChord;
204}
double G4double
Definition G4Types.hh:83
G4double GetRadHelix() const
G4double GetAngCurve() const

◆ DumbStepper()

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

Implements G4MagHelicalStepper.

Definition at line 169 of file G4HelixMixedStepper.cc.

173{
174 AdvanceHelix(yIn, Bfld, h, yOut);
175}
void AdvanceHelix(const G4double yIn[], const G4ThreeVector &Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=nullptr)

◆ GetAngleThreshold()

G4double G4HelixMixedStepper::GetAngleThreshold ( )
inline

Definition at line 94 of file G4HelixMixedStepper.hh.

94{ return fAngle_threshold; }

◆ IntegratorOrder()

G4int G4HelixMixedStepper::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 95 of file G4HelixMixedStepper.hh.

95{ return 4; }

◆ PrintCalls()

void G4HelixMixedStepper::PrintCalls ( )

Definition at line 207 of file G4HelixMixedStepper.cc.

208{
209 G4cout << "In HelixMixedStepper::Number of calls to smallStepStepper = "
210 << fNumCallsRK4
211 << " and Number of calls to Helix = " << fNumCallsHelix << G4endl;
212}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

Referenced by ~G4HelixMixedStepper().

◆ SetAngleThreshold()

void G4HelixMixedStepper::SetAngleThreshold ( G4double val)
inline

Definition at line 93 of file G4HelixMixedStepper.hh.

93{ fAngle_threshold = val; }

◆ SetupStepper()

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

Definition at line 216 of file G4HelixMixedStepper.cc.

217{
218 G4MagIntegratorStepper* pStepper;
219 if (fVerbose>0) { G4cout << " G4HelixMixedStepper: ";
220}
221 switch ( StepperNumber )
222 {
223 // Robust, classic method
224 case 4:
225 pStepper = new G4ClassicalRK4( pE );
226 if (fVerbose>0) { G4cout << "G4ClassicalRK4"; }
227 break;
228
229 // Steppers with embedded estimation of error
230 case 8:
231 pStepper = new G4CashKarpRKF45( pE );
232 if (fVerbose>0) { G4cout << "G4CashKarpRKF45"; }
233 break;
234 case 13:
235 pStepper = new G4NystromRK4( pE );
236 if (fVerbose>0) { G4cout << "G4NystromRK4"; }
237 break;
238
239 // Lowest order RK Stepper - experimental
240 case 1:
241 pStepper = new G4ImplicitEuler( pE );
242 if (fVerbose>0) { G4cout << "G4ImplicitEuler"; }
243 break;
244
245 // Lower order RK Steppers - ok overall, good for uneven fields
246 case 2:
247 pStepper = new G4SimpleRunge( pE );
248 if (fVerbose>0) { G4cout << "G4SimpleRunge"; }
249 break;
250 case 3:
251 pStepper = new G4SimpleHeum( pE );
252 if (fVerbose>0) { G4cout << "G4SimpleHeum"; }
253 break;
254 case 23:
255 pStepper = new G4BogackiShampine23( pE );
256 if (fVerbose>0) { G4cout << "G4BogackiShampine23"; }
257 break;
258
259 // Higher order RK Steppers
260 // for smoother fields and high accuracy requirements
261 case 45:
262 pStepper = new G4BogackiShampine45( pE );
263 if (fVerbose>0) { G4cout << "G4BogackiShampine45"; }
264 break;
265 case 145:
266 pStepper = new G4TsitourasRK45( pE );
267 if (fVerbose>0) { G4cout << "G4TsitourasRK45"; }
268 break;
269 case 745:
270 pStepper = new G4DormandPrince745( pE );
271 if (fVerbose>0) { G4cout << "G4DormandPrince745"; }
272 break;
273
274 // Helical Steppers
275 case 6:
276 pStepper = new G4HelixImplicitEuler( pE );
277 if (fVerbose>0) { G4cout << "G4HelixImplicitEuler"; }
278 break;
279 case 7:
280 pStepper = new G4HelixSimpleRunge( pE );
281 if (fVerbose>0) { G4cout << "G4HelixSimpleRunge"; }
282 break;
283 case 5:
284 pStepper = new G4HelixExplicitEuler( pE );
285 if (fVerbose>0) { G4cout << "G4HelixExplicitEuler"; }
286 break; // Since Helix Explicit is used for long steps,
287 // this is useful only to measure overhead.
288 // Exact Helix - likely good only for cases of
289 // i) uniform field (potentially over small distances)
290 // ii) segmented uniform field (maybe)
291 case 9:
292 pStepper = new G4ExactHelixStepper( pE );
293 if (fVerbose>0) { G4cout << "G4ExactHelixStepper"; }
294 break;
295 case 10:
296 pStepper = new G4RKG3_Stepper( pE );
297 if (fVerbose>0) { G4cout << "G4RKG3_Stepper"; }
298 break;
299
300 // Low Order Steppers - not good except for very weak fields
301 case 11:
302 pStepper = new G4ExplicitEuler( pE );
303 if (fVerbose>0) { G4cout << "G4ExplicitEuler"; }
304 break;
305 case 12:
306 pStepper = new G4ImplicitEuler( pE );
307 if (fVerbose>0) { G4cout << "G4ImplicitEuler"; }
308 break;
309
310 case 0:
311 case -1:
312 default:
313 pStepper = new G4DormandPrince745( pE ); // Was G4ClassicalRK4( pE );
314 if (fVerbose>0) { G4cout << "G4DormandPrince745 (Default)"; }
315 break;
316 }
317
318 if(fVerbose>0)
319 {
320 G4cout << " chosen as stepper for small steps in G4HelixMixedStepper."
321 << G4endl;
322 }
323
324 return pStepper;
325}

Referenced by G4HelixMixedStepper().

◆ SetVerbose()

void G4HelixMixedStepper::SetVerbose ( G4int newvalue)
inline

Definition at line 88 of file G4HelixMixedStepper.hh.

88{ fVerbose = newvalue; }

◆ Stepper()

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

Implements G4MagIntegratorStepper.

Definition at line 98 of file G4HelixMixedStepper.cc.

103{
104 // Estimation of the Stepping Angle
105 //
106 G4ThreeVector Bfld;
107 MagFieldEvaluate(yInput, Bfld);
108
109 G4double Bmag = Bfld.mag();
110 const G4double* pIn = yInput+3;
111 G4ThreeVector initVelocity = G4ThreeVector( pIn[0], pIn[1], pIn[2] );
112 G4double velocityVal = initVelocity.mag();
113
114 const G4double R_1 = std::abs(GetInverseCurve(velocityVal,Bmag));
115 // curv = inverse Radius
116 G4double Ang_curve = R_1 * Step;
117 // SetAngCurve(Ang_curve);
118 // SetCurve(std::abs(1/R_1));
119
120 if(Ang_curve < fAngle_threshold)
121 {
122 ++fNumCallsRK4;
123 fRK4Stepper->Stepper(yInput,dydx,Step,yOut,yErr);
124 }
125 else
126 {
127 constexpr G4int nvar = 6 ;
128 constexpr G4int nvarMax = 8 ;
129 G4double yTemp[nvarMax], yIn[nvarMax], yTemp2[nvarMax];
130 G4ThreeVector Bfld_midpoint;
131
132 SetAngCurve(Ang_curve);
133 SetCurve(std::abs(1.0/R_1));
134 ++fNumCallsHelix;
135
136 // Saving yInput because yInput and yOut can be aliases for same array
137 //
138 for(G4int i=0; i<nvar; ++i)
139 {
140 yIn[i]=yInput[i];
141 }
142
143 G4double halfS = Step * 0.5;
144
145 // 1. Do first half step and full step
146 //
147 AdvanceHelix(yIn, Bfld, halfS, yTemp, yTemp2); // yTemp2 for s=2*h (halfS)
148
149 MagFieldEvaluate(yTemp, Bfld_midpoint) ;
150
151 // 2. Do second half step - with revised field
152 // NOTE: Could avoid this call if 'Bfld_midpoint == Bfld'
153 // or diff 'almost' zero
154 //
155 AdvanceHelix(yTemp, Bfld_midpoint, halfS, yOut);
156 // Not requesting y at s=2*h (halfS)
157
158 // 3. Estimate the integration error
159 // should be (nearly) zero if Bfield= constant
160 //
161 for(G4int i=0; i<nvar; ++i)
162 {
163 yErr[i] = yOut[i] - yTemp2[i];
164 }
165 }
166}
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition G4Types.hh:85
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: