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

G4BorisDriver is a driver class using the second order Boris method to integrate the equation of motion. More...

#include <G4BorisDriver.hh>

Inheritance diagram for G4BorisDriver:

Public Member Functions

 G4BorisDriver (G4double hminimum, G4BorisScheme *Boris, G4int numberOfComponents=6, G4bool verbosity=false)
 ~G4BorisDriver () override=default
 G4BorisDriver (const G4BorisDriver &)=delete
G4BorisDriveroperator= (const G4BorisDriver &)=delete
G4bool AccurateAdvance (G4FieldTrack &track, G4double stepLen, G4double eps, G4double beginStep=0) override
G4bool QuickAdvance (G4FieldTrack &y_val, const G4double dydx[], G4double hstep, G4double &missDist, G4double &dyerr) override
void OneGoodStep (G4double yCurrentState[], G4double &curveLength, G4double htry, G4double epsilon_rel, G4double restMass, G4double charge, G4double &hdid, G4double &hnext)
G4double AdvanceChordLimited (G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance) override
void OnStartTracking () override
void OnComputeStep (const G4FieldTrack *) override
G4bool DoesReIntegrate () const override
G4double ComputeNewStepSize (G4double errMaxNorm, G4double hstepCurrent) override
G4double ShrinkStepSize2 (G4double h, G4double error2) const
G4double GrowStepSize2 (G4double h, G4double error2) const
void GetDerivatives (const G4FieldTrack &track, G4double dydx[]) const override
void GetDerivatives (const G4FieldTrack &track, G4double dydx[], G4double field[]) const override
void SetVerboseLevel (G4int level) override
G4int GetVerboseLevel () const override
G4EquationOfMotionGetEquationOfMotion () override
const G4EquationOfMotionGetEquationOfMotion () const
void SetEquationOfMotion (G4EquationOfMotion *equation) override
void StreamInfo (std::ostream &os) const override
const G4MagIntegratorStepperGetStepper () const override
G4MagIntegratorStepperGetStepper () override
Public Member Functions inherited from G4VIntegrationDriver
virtual ~G4VIntegrationDriver ()=default
virtual void RenewStepperAndAdjust (G4MagIntegratorStepper *pItsStepper)
Public Member Functions inherited from G4ChordFinderDelegate< G4BorisDriver >
virtual ~G4ChordFinderDelegate ()
G4double AdvanceChordLimitedImpl (G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance)
void ResetStepEstimate ()
G4double GetLastStepEstimateUnc ()
void SetLastStepEstimateUnc (G4double stepEst)
G4int GetNoCalls ()
G4int GetNoTrials ()
G4int GetNoMaxTrials ()
void SetFractions_Last_Next (G4double fractLast=0.90, G4double fractNext=0.95)
void SetFirstFraction (G4double fractFirst)
G4double GetFirstFraction ()
G4double GetFractionLast ()
G4double GetFractionNextEstimate ()
void StreamDelegateInfo (std::ostream &os) const
void TestChordPrint (G4int noTrials, G4int lastStepTrial, G4double dChordStep, G4double fDeltaChord, G4double nextStepTrial)

Additional Inherited Members

Static Protected Attributes inherited from G4VIntegrationDriver
static constexpr G4double max_stepping_increase = 5
static constexpr G4double max_stepping_decrease = 0.1

Detailed Description

G4BorisDriver is a driver class using the second order Boris method to integrate the equation of motion.

Definition at line 47 of file G4BorisDriver.hh.

Constructor & Destructor Documentation

◆ G4BorisDriver() [1/2]

G4BorisDriver::G4BorisDriver ( G4double hminimum,
G4BorisScheme * Boris,
G4int numberOfComponents = 6,
G4bool verbosity = false )

Constructor for G4BorisDriver.

Parameters
[in]hminimumThe minumum allowed step.
[in]BorisPointer to the Boris motion algorithm.
[in]numberOfComponentsThe number of integration variables.
[in]verbosityFlag for verbosity.

Definition at line 46 of file G4BorisDriver.cc.

49 : fMinimumStep(hminimum),
50 fVerbosity(verbosity),
51 boris(Boris)
52 // , interval_sequence{2,4}
53{
54 assert(boris->GetNumberOfVariables() == numberOfComponents);
55
56 if(boris->GetNumberOfVariables() != numberOfComponents)
57 {
58 std::ostringstream msg;
59 msg << "Disagreement in number of variables = "
60 << boris->GetNumberOfVariables()
61 << " vs no of components = " << numberOfComponents;
62 G4Exception("G4BorisDriver Constructor:",
63 "GeomField1001", FatalException, msg);
64 }
65
66}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

Referenced by G4BorisDriver(), and operator=().

◆ ~G4BorisDriver()

G4BorisDriver::~G4BorisDriver ( )
overridedefault

Default Destructor.

◆ G4BorisDriver() [2/2]

G4BorisDriver::G4BorisDriver ( const G4BorisDriver & )
delete

Copy constructor and assignment operator not allowed.

Member Function Documentation

◆ AccurateAdvance()

G4bool G4BorisDriver::AccurateAdvance ( G4FieldTrack & track,
G4double stepLen,
G4double eps,
G4double beginStep = 0 )
overridevirtual

Advances integration accurately by relative accuracy better than 'eps'.

Parameters
[in,out]trackThe current track in field.
[in]stepLenProposed step length.
[in]epsilonRequested accuracy, y_err/hstep.
[in]beginStepInitial minimum integration step.
Returns
true if integration succeeds.

Implements G4VIntegrationDriver.

Definition at line 70 of file G4BorisDriver.cc.

74{
75 // Specification: Driver with adaptive stepsize control.
76 // Integrate starting values at y_current over hstep x2 with (relative) accuracy 'eps'.
77 // On output 'track' is replaced by values at the end of the integration interval.
78
79 // Ensure that hstep > 0
80 if(hstep == 0)
81 {
82 std::ostringstream message;
83 message << "Proposed step is zero; hstep = " << hstep << " !";
84 G4Exception("G4BorisDriver::AccurateAdvance()",
85 "GeomField1001", JustWarning, message);
86 return true;
87 }
88 if(hstep < 0)
89 {
90 std::ostringstream message;
91 message << "Invalid run condition." << G4endl
92 << "Proposed step is negative; hstep = " << hstep << G4endl
93 << "Requested step cannot be negative! Aborting event.";
94 G4Exception("G4BorisDriver::AccurateAdvance()",
95 "GeomField0003", EventMustBeAborted, message);
96 return false;
97 }
98
99 if( hinitial == 0.0 ) { hinitial = hstep; }
100 if( hinitial < 0.0 ) { hinitial = std::fabs( hinitial ); }
101 // G4double htrial = std::min( hstep, hinitial );
102 G4double htrial = hstep;
103 // Decide first step size
104
105 // G4int noOfSteps = h/hstep;
106
107 // integration variables
108 //
109 track.DumpToArray(yCurrent);
110
111 const G4double restMass = track.GetRestMass();
112 const G4double charge = track.GetCharge()*e_SI;
113 const G4int nvar= GetNumberOfVariables();
114
115 // copy non-integration variables to out array
116 //
117 std::memcpy(yOut + nvar,
118 yCurrent + nvar,
119 sizeof(G4double)*(G4FieldTrack::ncompSVEC-nvar));
120
121 G4double curveLength = track.GetCurveLength(); // starting value
122 const G4double endCurveLength = curveLength + hstep;
123
124 // -- Initial version: Did it in one step -- did not account for errors !!!
125 // G4FieldTrack yFldTrk(track);
126 // yFldTrk.LoadFromArray(yCurrent, G4FieldTrack::ncompSVEC);
127 // yFldTrk.SetCurveLength(curveLength);
128 // G4double dchord_step, dyerr_len;
129 // QuickAdvance(yFldTrk, dydxCurrent, htrial, dchord_step, dyerr_len);
130
131 const G4double hThreshold =
132 std::max(epsilon * hstep, fSmallestFraction * curveLength);
133
134 G4double htry= htrial;
135
136 for (G4int nstp = 0; nstp < fMaxNoSteps; ++nstp)
137 {
138 G4double hdid= 0.0, hnext=0.0;
139
140 OneGoodStep(yCurrent, curveLength, htry, epsilon, restMass, charge, hdid, hnext);
141 //*********
142
143 // Simple check: move (distance of displacement) is smaller than length along curve!
146 CheckStep(EndPos, StartPos, hdid);
147
148 // Check 1. for finish and 2. *avoid* numerous small last steps
149 if (curveLength >= endCurveLength || htry < hThreshold)
150 {
151 break;
152 }
153
154 htry = std::max(hnext, fMinimumStep);
155 if (curveLength + htry > endCurveLength)
156 {
157 htry = endCurveLength - curveLength;
158 }
159 }
160
161 // upload new state
162 track.LoadFromArray(yCurrent, G4FieldTrack::ncompSVEC);
163 track.SetCurveLength(curveLength);
164
165 return true;
166}
G4double epsilon(G4double density, G4double temperature)
@ JustWarning
@ EventMustBeAborted
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
void OneGoodStep(G4double yCurrentState[], G4double &curveLength, G4double htry, G4double epsilon_rel, G4double restMass, G4double charge, G4double &hdid, G4double &hnext)
G4double GetCurveLength() const
G4double GetCharge() const
void SetCurveLength(G4double nCurve_s)
G4double GetRestMass() const
void DumpToArray(G4double valArr[ncompSVEC]) const
void LoadFromArray(const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
G4ThreeVector makeVector(const ArrayType &array, Value3D value)

◆ AdvanceChordLimited()

G4double G4BorisDriver::AdvanceChordLimited ( G4FieldTrack & track,
G4double hstep,
G4double eps,
G4double chordDistance )
inlineoverridevirtual

Computes the step to take, based on chord limits.

Parameters
[in,out]trackThe current track in field.
[in]hstepProposed step length.
[in]epsRequested accuracy, y_err/hstep.
[in]chordDistanceMaximum sagitta distance.
Returns
The length of step taken.

Implements G4VIntegrationDriver.

◆ ComputeNewStepSize()

G4double G4BorisDriver::ComputeNewStepSize ( G4double errMaxNorm,
G4double hstepCurrent )
inlineoverridevirtual

Computes a step size for the next step, taking the last step's normalised error 'errMaxNorm'.

Parameters
[in]errMaxNormThe normalised error on last step.
[in]hstepCurrentThe current proposed step.
Returns
The step size for the next step.

Implements G4VIntegrationDriver.

◆ DoesReIntegrate()

G4bool G4BorisDriver::DoesReIntegrate ( ) const
inlineoverridevirtual

The driver implements re-integration. Returns true. It would be false if it just used interpolation to provide a result.

Implements G4VIntegrationDriver.

◆ GetDerivatives() [1/2]

void G4BorisDriver::GetDerivatives ( const G4FieldTrack & track,
G4double dydx[] ) const
overridevirtual

Getters for derivatives.

Implements G4VIntegrationDriver.

Definition at line 269 of file G4BorisDriver.cc.

271{
272 G4double EBfieldValue[6];
273 GetDerivatives(yTrack, dydx, EBfieldValue);
274}
void GetDerivatives(const G4FieldTrack &track, G4double dydx[]) const override

Referenced by GetDerivatives().

◆ GetDerivatives() [2/2]

void G4BorisDriver::GetDerivatives ( const G4FieldTrack & track,
G4double dydx[],
G4double field[] ) const
overridevirtual

Implements G4VIntegrationDriver.

Definition at line 278 of file G4BorisDriver.cc.

281{
282 // G4Exception("G4BorisDriver::GetDerivatives()",
283 // "GeomField0003", FatalException, "This method is not implemented.");
285 yTrack.DumpToArray(ytemp);
286 GetEquationOfMotion()->EvaluateRhsReturnB(ytemp, dydx, EBfieldValue);
287}
G4EquationOfMotion * GetEquationOfMotion() override
void EvaluateRhsReturnB(const G4double y[], G4double dydx[], G4double Field[]) const

◆ GetEquationOfMotion() [1/2]

const G4EquationOfMotion * G4BorisDriver::GetEquationOfMotion ( ) const
inline

◆ GetEquationOfMotion() [2/2]

G4EquationOfMotion * G4BorisDriver::GetEquationOfMotion ( )
inlineoverridevirtual

Getters for the equation of motion.

Implements G4VIntegrationDriver.

Referenced by GetDerivatives().

◆ GetStepper() [1/2]

const G4MagIntegratorStepper * G4BorisDriver::GetStepper ( ) const
inlineoverridevirtual

Accessors for stepper. Not relevant for Boris and other non-RK methods.

Implements G4VIntegrationDriver.

◆ GetStepper() [2/2]

G4MagIntegratorStepper * G4BorisDriver::GetStepper ( )
inlineoverridevirtual

Implements G4VIntegrationDriver.

◆ GetVerboseLevel()

G4int G4BorisDriver::GetVerboseLevel ( ) const
inlineoverridevirtual

Implements G4VIntegrationDriver.

◆ GrowStepSize2()

G4double G4BorisDriver::GrowStepSize2 ( G4double h,
G4double error2 ) const

Definition at line 302 of file G4BorisDriver.cc.

304{
305 if (error2 < fErrorConstraintGrow * fErrorConstraintGrow)
306 {
307 return fMaxSteppingIncrease * h;
308 }
309 return fSafetyFactor * h * std::pow(error2, 0.5 * fPowerGrow);
310}

Referenced by OneGoodStep().

◆ OnComputeStep()

void G4BorisDriver::OnComputeStep ( const G4FieldTrack * )
inlineoverridevirtual

Dispatch interface method for computing step. Does nothing here.

Implements G4VIntegrationDriver.

◆ OneGoodStep()

void G4BorisDriver::OneGoodStep ( G4double yCurrentState[],
G4double & curveLength,
G4double htry,
G4double epsilon_rel,
G4double restMass,
G4double charge,
G4double & hdid,
G4double & hnext )

Takes one Step that is as large as possible while satisfying the accuracy criterion.

Parameters
[in,out]yCurrentStateThe current track state, y.
[in,out]curveLengthStep start, x.
[in]htryStep to attempt.
[in]epsilon_relThe relative accuracy.
[in]restMassMass value for computing velocity.
[in]chargeCharge value for computing momentum.
[out]hdidStep achieved.
[out]hnextProposed next step.
Returns
true if integration succeeds.

Definition at line 170 of file G4BorisDriver.cc.

178{
179 G4double error2 = DBL_MAX;
181
182 G4double h = htry;
183
184 const G4int max_trials = 100;
185
186 for (G4int iter = 0; iter < max_trials; ++iter)
187 {
188 boris->StepWithErrorEstimate(y, restMass, charge, h, ytemp, yerr);
189
190 error2 = field_utils::relativeError2(y, yerr, std::max(h, fMinimumStep),
191 epsilon_rel);
192 if (error2 <= 1.0)
193 {
194 break;
195 }
196
197 h = ShrinkStepSize2(h, error2);
198
199 G4double xnew = curveLength + h;
200 if(xnew == curveLength)
201 {
202 std::ostringstream message;
203 message << "Stepsize underflow in Stepper !" << G4endl
204 << "- Step's start x=" << curveLength
205 << " and end x= " << xnew
206 << " are equal !! " << G4endl
207 << " Due to step-size= " << h
208 << ". Note that input step was " << htry;
209 G4Exception("G4IntegrationDriver::OneGoodStep()",
210 "GeomField1001", JustWarning, message);
211 break;
212 }
213 }
214
215 hnext = GrowStepSize2(h, error2);
216 curveLength += (hdid = h);
217
218 field_utils::copy(y, ytemp, GetNumberOfVariables());
219}
G4double ShrinkStepSize2(G4double h, G4double error2) const
G4double GrowStepSize2(G4double h, G4double error2) const
G4double relativeError2(const G4double y[], const G4double yerr[], G4double hstep, G4double errorTolerance)
void copy(G4double dst[], const G4double src[], std::size_t size=G4FieldTrack::ncompSVEC)
#define DBL_MAX
Definition templates.hh:62

Referenced by AccurateAdvance().

◆ OnStartTracking()

void G4BorisDriver::OnStartTracking ( )
inlineoverridevirtual

Dispatch interface method for initialisation/reset of driver.

Implements G4VIntegrationDriver.

◆ operator=()

G4BorisDriver & G4BorisDriver::operator= ( const G4BorisDriver & )
delete

◆ QuickAdvance()

G4bool G4BorisDriver::QuickAdvance ( G4FieldTrack & y_val,
const G4double dydx[],
G4double hstep,
G4double & missDist,
G4double & dyerr )
overridevirtual

Attempts one integration step, and returns estimated error 'dyerr'. It does not ensure accuracy.

Parameters
[in,out]y_valThe current track in field.
[in]dydxdydx array.
[in]hstepProposed step length.
[out]missDistEstimated sagitta distance.
[out]dyerrEstimated error.
Returns
true if integration succeeds.

Reimplemented from G4VIntegrationDriver.

Definition at line 223 of file G4BorisDriver.cc.

226{
227 const auto nvar = boris->GetNumberOfVariables();
228
229 track.DumpToArray(yIn);
230 const G4double curveLength = track.GetCurveLength();
231
232 // call the boris method for step length hstep
233 G4double restMass = track.GetRestMass();
234 G4double charge = track.GetCharge()*e_SI;
235
236 // boris->DoStep(restMass, charge, yIn, yMid, hstep*0.5);
237 // boris->DoStep(restMass, charge, yMid, yOut, hstep*0.5); // Use mid-point !!
238 boris->StepWithMidAndErrorEstimate(yIn, restMass, charge, hstep,
239 yMid, yOut, yError);
240 // Same, and also return mid-point evaluation
241
242 // How to calculate chord length??
243 const auto mid = field_utils::makeVector(yMid,
245 const auto in = field_utils::makeVector(yIn,
247 const auto out = field_utils::makeVector(yOut,
249
250 missDist = G4LineSection::Distline(mid, in, out);
251
252 dyerr = field_utils::absoluteError(yOut, yError, hstep);
253
254 // copy non-integrated variables to output array
255 //
256 std::memcpy(yOut + nvar, yIn + nvar,
257 sizeof(G4double) * (G4FieldTrack::ncompSVEC - nvar));
258
259 // set new state
260 //
261 track.LoadFromArray(yOut, G4FieldTrack::ncompSVEC);
262 track.SetCurveLength(curveLength + hstep);
263
264 return true;
265}
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4double absoluteError(const G4double y[], const G4double yerr[], G4double hstep)

◆ SetEquationOfMotion()

void G4BorisDriver::SetEquationOfMotion ( G4EquationOfMotion * equation)
overridevirtual

Setter for the equation of motion. Issues an exception, as not foreseen to change equation of motion for the Boris stepper.

Implements G4VIntegrationDriver.

Definition at line 314 of file G4BorisDriver.cc.

315{
316 G4Exception("G4BorisDriver::SetEquationOfMotion()", "GeomField0003", FatalException,
317 "This method is not implemented. BorisDriver/Stepper should keep its equation");
318}

◆ SetVerboseLevel()

void G4BorisDriver::SetVerboseLevel ( G4int level)
inlineoverridevirtual

Setter and getter for verbosity.

Implements G4VIntegrationDriver.

◆ ShrinkStepSize2()

G4double G4BorisDriver::ShrinkStepSize2 ( G4double h,
G4double error2 ) const

Methods to calculate the next step size given the square of the relative error.

Definition at line 291 of file G4BorisDriver.cc.

292{
293 if (error2 > fErrorConstraintShrink * fErrorConstraintShrink)
294 {
295 return fMaxSteppingDecrease * h;
296 }
297 return fSafetyFactor * h * std::pow(error2, 0.5 * fPowerShrink);
298}

Referenced by OneGoodStep().

◆ StreamInfo()

void G4BorisDriver::StreamInfo ( std::ostream & os) const
overridevirtual

Writes out to stream the parameters/state of the driver.

Implements G4VIntegrationDriver.

Definition at line 323 of file G4BorisDriver.cc.

324{
325 os << "State of G4BorisDriver: " << std::endl;
326 os << " Method is implemented, but gives no information. " << std::endl;
327}

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