Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4QSStepper< QSS > Class Template Reference

#include <G4QSStepper.hh>

+ Inheritance diagram for G4QSStepper< QSS >:

Public Member Functions

 G4QSStepper (G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
 
 ~G4QSStepper () override
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
 
void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])
 
void SetupInterpolation ()
 
void Interpolate (G4double tau, G4double yOut[])
 
G4double DistChord () const override
 
G4int IntegratorOrder () const override
 
void reset (const G4FieldTrack *track)
 
void SetPrecision (G4double dq_rel, G4double dq_min)
 
G4EquationOfMotionGetSpecificEquation ()
 
const field_utils::StateGetYOut () const
 
G4double GetLastStepLength ()
 
- 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)
 

Static Public Member Functions

static G4QSStepper< G4QSS2 > * build_QSS2 (G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
 
static G4QSStepper< G4QSS3 > * build_QSS3 (G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
 

Additional Inherited Members

- Protected Member Functions inherited from G4MagIntegratorStepper
void SetIntegrationOrder (G4int order)
 
void SetFSAL (G4bool flag=true)
 

Detailed Description

template<class QSS>
class G4QSStepper< QSS >

Definition at line 53 of file G4QSStepper.hh.

Constructor & Destructor Documentation

◆ G4QSStepper()

template<class QSS>
G4QSStepper< QSS >::G4QSStepper ( G4EquationOfMotion * EqRhs,
G4int numberOfVariables = 6,
G4bool primary = true )

Referenced by DistChord().

◆ ~G4QSStepper()

template<class QSS>
G4QSStepper< QSS >::~G4QSStepper ( )
inlineoverride

Definition at line 192 of file G4QSStepper.hh.

193{
194 for (auto & i : simulator->SD) { free(i); }
195
196 free(SUBSTEPS(this->simulator));
197 free(this->simulator);
198}
#define SUBSTEPS(sim)
Definition G4qss_misc.hh:86

Member Function Documentation

◆ build_QSS2()

template<class QSS>
static G4QSStepper< G4QSS2 > * G4QSStepper< QSS >::build_QSS2 ( G4EquationOfMotion * EqRhs,
G4int numberOfVariables = 6,
G4bool primary = true )
static

◆ build_QSS3()

template<class QSS>
static G4QSStepper< G4QSS3 > * G4QSStepper< QSS >::build_QSS3 ( G4EquationOfMotion * EqRhs,
G4int numberOfVariables = 6,
G4bool primary = true )
static

◆ DistChord()

template<class QSS>
G4double G4QSStepper< QSS >::DistChord ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 331 of file G4QSStepper.hh.

332{
333 G4double yMid[6];
334 const_cast<G4QSStepper<QSS>*>(this)->Interpolate(0.5, yMid);
335
339
341}
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
void Interpolate(G4double tau, G4double yOut[])
G4QSStepper(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)

◆ GetLastStepLength()

template<class QSS>
G4double G4QSStepper< QSS >::GetLastStepLength ( )
inline

Definition at line 101 of file G4QSStepper.hh.

101{ return fLastStepLength; }

◆ GetSpecificEquation()

template<class QSS>
G4EquationOfMotion * G4QSStepper< QSS >::GetSpecificEquation ( )
inline

Definition at line 97 of file G4QSStepper.hh.

97{ return GetEquationOfMotion(); }
G4EquationOfMotion * GetEquationOfMotion()

◆ GetYOut()

template<class QSS>
const field_utils::State & G4QSStepper< QSS >::GetYOut ( ) const
inline

Definition at line 99 of file G4QSStepper.hh.

99{ return fyOut; }

◆ IntegratorOrder()

template<class QSS>
G4int G4QSStepper< QSS >::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 82 of file G4QSStepper.hh.

82{ return method->order(); }

◆ Interpolate()

template<class QSS>
void G4QSStepper< QSS >::Interpolate ( G4double tau,
G4double yOut[] )
inline

Definition at line 344 of file G4QSStepper.hh.

345{
346 G4double length = tau * fLastStepLength;
347 G4int idx = 0, j = LAST_SUBSTEP(simulator);
349
350 if (j >= 15) {
351 G4int i = 0, k = j;
352 idx = j >> 1;
353 while (idx < k && i < j - 1) {
354 if (length < SUBSTEP_LEN(simulator, idx)) {
355 j = idx;
356 } else if (length >= SUBSTEP_LEN(simulator, idx + 1)) {
357 i = idx;
358 } else {
359 break;
360 }
361
362 idx = (i + j) >> 1;
363 }
364 }
365 else {
366 for (; idx < j && length >= SUBSTEP_LEN(simulator, idx + 1); idx++) {;}
367 }
368
369 auto* const substep = &SUBSTEP_STRUCT(simulator, idx);
370 end_time = substep->start_time + (length - substep->len) / fVelocity;
371
372 this->get_state_from_poly(substep->x, substep->tx, end_time, yOut);
373
374 velocity_to_momentum(yOut);
375}
#define LAST_SUBSTEP(sim)
Definition G4qss_misc.hh:84
#define SUBSTEP_STRUCT(sim, i)
Definition G4qss_misc.hh:70
#define SUBSTEP_LEN(sim, i)
Definition G4qss_misc.hh:74

Referenced by DistChord().

◆ reset()

template<class QSS>
void G4QSStepper< QSS >::reset ( const G4FieldTrack * track)
inline

Definition at line 378 of file G4QSStepper.hh.

379{
380 using Qss_misc::PXidx;
381 using Qss_misc::PYidx;
382 using Qss_misc::PZidx;
383 using Qss_misc::VXidx;
384 using Qss_misc::VYidx;
385 using Qss_misc::VZidx;
386
387 G4ThreeVector pos = track->GetPosition();
388 G4ThreeVector momentum = track->GetMomentum();
389
390 f_mass = track->GetRestMass();
391 set_charge(track->GetCharge());
392 set_gamma(momentum.mag2());
393 G4double c_mg = cLight_local / fMassGamma;
394 set_velocity(momentum.mag() * c_mg);
395
396 method->reset_state(PXidx, pos.getX());
397 method->reset_state(PYidx, pos.getY());
398 method->reset_state(PZidx, pos.getZ());
399
400 method->reset_state(VXidx, momentum.getX() * c_mg);
401 method->reset_state(VYidx, momentum.getY() * c_mg);
402 method->reset_state(VZidx, momentum.getZ() * c_mg);
403
404 this->update_field();
405 method->full_definition(get_coeff());
406
407 method->recompute_all_state_times(0);
408
409 simulator->time = 0;
410}

◆ SetPrecision()

template<class QSS>
void G4QSStepper< QSS >::SetPrecision ( G4double dq_rel,
G4double dq_min )
inline

Definition at line 413 of file G4QSStepper.hh.

414{
415 G4double* dQMin = simulator->dQMin;
416 G4double* dQRel = simulator->dQRel;
417 G4int n_vars = simulator->states;
418
419 if (dq_min <= 0) { dq_min = dq_rel * 1e-3; }
420
421 for (G4int i = 0; i < n_vars; ++i) {
422 dQRel[i] = dq_rel;
423 dQMin[i] = dq_min;
424 }
425}

◆ SetupInterpolation()

template<class QSS>
void G4QSStepper< QSS >::SetupInterpolation ( )
inline

Definition at line 77 of file G4QSStepper.hh.

77{}

◆ Stepper() [1/2]

template<class QSS>
void G4QSStepper< QSS >::Stepper ( const G4double y[],
const G4double dydx[],
G4double h,
G4double yout[],
G4double yerr[] )
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 242 of file G4QSStepper.hh.

247{
249 G4double t, prev_time = 0;
250 G4double length = 0.;
251 G4int index;
252
253 const G4int coeffs = method->order() + 1;
254 G4double* tq = simulator->tq;
255 G4double* tx = simulator->tx;
256 G4double* dQRel = simulator->dQRel;
257 G4double* dQMin = simulator->dQMin;
258 G4double* lqu = simulator->lqu;
259 G4double* x = simulator->x;
260 G4int** SD = simulator->SD;
262
263 CUR_SUBSTEP(simulator) = 0;
264
265 this->save_substep(0, length);
266
267 this->update_time();
268 t = simulator->time;
269 index = simulator->minIndex;
270
271 while (length < max_length && t < Qss_misc::INF && CUR_SUBSTEP(simulator) < QSS_MAX_SUBSTEPS) {
272 cf0 = index * coeffs;
273 elapsed = t - tx[index];
274 method->advance_time_x(cf0, elapsed);
275 tx[index] = t;
277 if (lqu[index] < dQMin[index]) {
278 lqu[index] = dQMin[index];
279 }
280 method->update_quantized_state(index);
281 tq[index] = t;
282 method->next_time(index, t);
283 for (G4int i = 0; i < 3; i++) {
284 G4int j = SD[index][i];
285 elapsed = t - tx[j];
286 infCf0 = j * coeffs;
287 if (elapsed > 0) {
288 x[infCf0] = method->evaluate_x_poly(infCf0, elapsed, x);
289 tx[j] = t;
290 }
291 }
292
293 this->update_field();
294 this->recompute_derivatives(index);
295 method->recompute_next_times(SD[index], t);
296
297 if (t > prev_time) {
298 length += fVelocity * (t - prev_time);
299 if (length <= max_length) { this->save_substep(t, length); }
300 else { break; }
301 }
302
303 this->update_time();
304 prev_time = t;
305 t = simulator->time;
306 index = simulator->minIndex;
307 }
308
309 if(CUR_SUBSTEP(simulator) >= QSS_MAX_SUBSTEPS) {
311 }
312
313 auto* const substep = &LAST_SUBSTEP_STRUCT(simulator);
314 t = substep->start_time + (max_length - substep->len) / fVelocity;
315
316 this->get_state_from_poly(substep->x, substep->tx, t, yOut);
317
318 velocity_to_momentum(yOut);
319
321 for (G4int i = 0; i < numberOfVariables; ++i) {
322 // Store Input and Final values, for possible use in calculating chord
323 fyIn[i] = yInput[i];
324 fyOut[i] = yOut[i];
325 }
326
327 fLastStepLength = max_length;
328}
#define CUR_SUBSTEP(sim)
Definition G4qss_misc.hh:83
#define LAST_SUBSTEP_STRUCT(sim)
Definition G4qss_misc.hh:76
G4int GetNumberOfVariables() const

Referenced by Stepper().

◆ Stepper() [2/2]

template<class QSS>
void G4QSStepper< QSS >::Stepper ( const G4double yInput[],
const G4double dydx[],
G4double hstep,
G4double yOutput[],
G4double yError[],
G4double dydxOutput[] )
inline

Definition at line 201 of file G4QSStepper.hh.

207{
209}
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override

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