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

#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 epsilon, 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 ()
 
void TestChordPrint (G4int noTrials, G4int lastStepTrial, G4double dChordStep, G4double fDeltaChord, G4double nextStepTrial)
 
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 ()
 
G4double GetLastStepEstimateUnc ()
 
void SetLastStepEstimateUnc (G4double stepEst)
 
void StreamDelegateInfo (std::ostream &os) const
 

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

Definition at line 44 of file G4BorisDriver.hh.

Constructor & Destructor Documentation

◆ G4BorisDriver() [1/2]

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

Definition at line 52 of file G4BorisDriver.cc.

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

◆ ~G4BorisDriver()

G4BorisDriver::~G4BorisDriver ( )
inlineoverridedefault

◆ G4BorisDriver() [2/2]

G4BorisDriver::G4BorisDriver ( const G4BorisDriver & )
inlinedelete

Member Function Documentation

◆ AccurateAdvance()

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

Implements G4VIntegrationDriver.

Definition at line 76 of file G4BorisDriver.cc.

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

Implements G4VIntegrationDriver.

Definition at line 87 of file G4BorisDriver.hh.

91 {
93 AdvanceChordLimitedImpl(track, hstep, eps, chordDistance);
94 }
G4double AdvanceChordLimitedImpl(G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance)

◆ ComputeNewStepSize()

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

Implements G4VIntegrationDriver.

◆ DoesReIntegrate()

G4bool G4BorisDriver::DoesReIntegrate ( ) const
inlineoverridevirtual

Implements G4VIntegrationDriver.

Definition at line 106 of file G4BorisDriver.hh.

106{ return true; }

◆ GetDerivatives() [1/2]

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

Implements G4VIntegrationDriver.

Definition at line 275 of file G4BorisDriver.cc.

277{
278 G4double EBfieldValue[6];
279 GetDerivatives(yTrack, dydx, EBfieldValue);
280}
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 284 of file G4BorisDriver.cc.

287{
288 // G4Exception("G4BorisDriver::GetDerivatives()",
289 // "GeomField0003", FatalException, "This method is not implemented.");
291 yTrack.DumpToArray(ytemp);
292 GetEquationOfMotion()->EvaluateRhsReturnB(ytemp, dydx, EBfieldValue);
293}
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

Implements G4VIntegrationDriver.

Referenced by GetDerivatives().

◆ GetStepper() [1/2]

const G4MagIntegratorStepper * G4BorisDriver::GetStepper ( ) const
inlineoverridevirtual

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 308 of file G4BorisDriver.cc.

310{
311 if (error2 < fErrorConstraintGrow * fErrorConstraintGrow)
312 {
313 return fMaxSteppingIncrease * h;
314 }
315 return fSafetyFactor * h * std::pow(error2, 0.5 * fPowerGrow);
316}

Referenced by OneGoodStep().

◆ OnComputeStep()

void G4BorisDriver::OnComputeStep ( const G4FieldTrack * )
inlineoverridevirtual

Implements G4VIntegrationDriver.

Definition at line 101 of file G4BorisDriver.hh.

101{}

◆ OneGoodStep()

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

Definition at line 176 of file G4BorisDriver.cc.

184{
185 G4double error2 = DBL_MAX;
187
188 G4double h = htry;
189
190 const G4int max_trials = 100;
191
192 for (G4int iter = 0; iter < max_trials; ++iter)
193 {
194 boris->StepWithErrorEstimate(y, restMass, charge, h, ytemp, yerr);
195
196 error2 = field_utils::relativeError2(y, yerr, std::max(h, fMinimumStep),
197 epsilon_rel);
198 if (error2 <= 1.0)
199 {
200 break;
201 }
202
203 h = ShrinkStepSize2(h, error2);
204
205 G4double xnew = curveLength + h;
206 if(xnew == curveLength)
207 {
208 std::ostringstream message;
209 message << "Stepsize underflow in Stepper !" << G4endl
210 << "- Step's start x=" << curveLength
211 << " and end x= " << xnew
212 << " are equal !! " << G4endl
213 << " Due to step-size= " << h
214 << ". Note that input step was " << htry;
215 G4Exception("G4IntegrationDriver::OneGoodStep()",
216 "GeomField1001", JustWarning, message);
217 break;
218 }
219 }
220
221 hnext = GrowStepSize2(h, error2);
222 curveLength += (hdid = h);
223
224 field_utils::copy(y, ytemp, GetNumberOfVariables());
225}
G4double ShrinkStepSize2(G4double h, G4double error2) const
G4double GrowStepSize2(G4double h, G4double error2) const
void StepWithErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep, G4double yOut[], G4double yErr[]) 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

◆ operator=()

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

◆ QuickAdvance()

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

Reimplemented from G4VIntegrationDriver.

Definition at line 229 of file G4BorisDriver.cc.

232{
233 const auto nvar = boris->GetNumberOfVariables();
234
235 track.DumpToArray(yIn);
236 const G4double curveLength = track.GetCurveLength();
237
238 // call the boris method for step length hstep
239 G4double restMass = track.GetRestMass();
240 G4double charge = track.GetCharge()*e_SI;
241
242 // boris->DoStep(restMass, charge, yIn, yMid, hstep*0.5);
243 // boris->DoStep(restMass, charge, yMid, yOut, hstep*0.5); // Use mid-point !!
244 boris->StepWithMidAndErrorEstimate(yIn, restMass, charge, hstep,
245 yMid, yOut, yError);
246 // Same, and also return mid-point evaluation
247
248 // How to calculate chord length??
249 const auto mid = field_utils::makeVector(yMid,
251 const auto in = field_utils::makeVector(yIn,
253 const auto out = field_utils::makeVector(yOut,
255
256 missDist = G4LineSection::Distline(mid, in, out);
257
258 dyerr = field_utils::absoluteError(yOut, yError, hstep);
259
260 // copy non-integrated variables to output array
261 //
262 std::memcpy(yOut + nvar, yIn + nvar,
263 sizeof(G4double) * (G4FieldTrack::ncompSVEC - nvar));
264
265 // set new state
266 //
267 track.LoadFromArray(yOut, G4FieldTrack::ncompSVEC);
268 track.SetCurveLength(curveLength + hstep);
269
270 return true;
271}
void StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep, G4double yMid[], G4double yOut[], G4double yErr[]) const
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

Implements G4VIntegrationDriver.

Definition at line 320 of file G4BorisDriver.cc.

321{
322 G4Exception("G4BorisDriver::SetEquationOfMotion()", "GeomField0003", FatalException,
323 "This method is not implemented. BorisDriver/Stepper should keep its equation");
324}

◆ SetVerboseLevel()

void G4BorisDriver::SetVerboseLevel ( G4int level)
inlineoverridevirtual

Implements G4VIntegrationDriver.

◆ ShrinkStepSize2()

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

Definition at line 297 of file G4BorisDriver.cc.

298{
299 if (error2 > fErrorConstraintShrink * fErrorConstraintShrink)
300 {
301 return fMaxSteppingDecrease * h;
302 }
303 return fSafetyFactor * h * std::pow(error2, 0.5 * fPowerShrink);
304}

Referenced by OneGoodStep().

◆ StreamInfo()

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

Implements G4VIntegrationDriver.

Definition at line 329 of file G4BorisDriver.cc.

330{
331 os << "State of G4BorisDriver: " << std::endl;
332 os << " Method is implemented, but gives no information. " << std::endl;
333}

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