Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QSSDriver.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4QSSDriver
27//
28// QSS Interpolator Driver
29
30// Authors: Lucio Santi, Rodrigo Castro (Univ. Buenos Aires), 2018-2021
31// --------------------------------------------------------------------
32#ifndef G4QSSDriver_HH
33#define G4QSSDriver_HH
34
36#include "G4QSSMessenger.hh"
37
38/**
39 * @brief G4QSSDriver is a templated driver class defining the QSS
40 * (Quantum State Simulation) Interpolator Driver.
41 */
42
43template <class T>
44class G4QSSDriver : public G4InterpolationDriver<T, true>
45{
46 public:
47
48 /**
49 * Constructor for G4QSSDriver.
50 * @param[in] T Pointer to the stepper algorithm.
51 */
52 inline G4QSSDriver(T* stepper);
53
54 /**
55 * Copy constructor and assignment operator not allowed.
56 */
57 G4QSSDriver(const G4QSSDriver&) = delete;
58 const G4QSSDriver& operator=(const G4QSSDriver&) = delete;
59
60 /**
61 * Dispatch interface method for initialisation/reset of driver.
62 */
63 void OnStartTracking() override;
64
65 /**
66 * Computes the step to take, based on chord limits.
67 * @param[in,out] track The current track in field.
68 * @param[in] hstep Proposed step length.
69 * @param[in] eps Requested accuracy, y_err/hstep.
70 * @param[in] chordDistance Maximum sagitta distance.
71 * @returns The length of step taken.
72 */
74 G4double hstep,
75 G4double eps,
76 G4double chordDistance) override;
77
78 /**
79 * Dispatch interface method for computing step.
80 */
81 inline void OnComputeStep(const G4FieldTrack* track) override;
82
83 /**
84 * Setter for driver precision parameters.
85 */
86 inline void SetPrecision(G4double dq_rel, G4double dq_min);
87
88 /**
89 * Takes one Step that is as large as possible while satisfying the
90 * accuracy criterion.
91 * @param[in] it Stepper iterator.
92 * @param[in,out] y The current track state, y.
93 * @param[in] dydx dydx array.
94 * @param[in,out] hstep Step to attempt.
95 * @param[in] epsStep The relative accuracy.
96 * @param[in] curveLength Step start, x.
97 * @param[in,out] track Pointer to the Field track. Not used.
98 * @returns The step achieved.
99 */
102 field_utils::State& dydx,
103 G4double& hstep,
104 G4double epsStep,
105 G4double curveLength,
106 G4FieldTrack* track) override;
107
108 private:
109
111
112 G4bool initializedOnFirstRun = false;
113};
114
115#include "G4QSSDriver.icc"
116
117#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
G4InterpolationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
typename std::vector< InterpStepper >::iterator StepperIterator
void SetPrecision(G4double dq_rel, G4double dq_min)
G4QSSDriver(T *stepper)
G4double AdvanceChordLimited(G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance) override
void OnStartTracking() override
void OnComputeStep(const G4FieldTrack *track) override
G4QSSDriver(const G4QSSDriver &)=delete
G4double OneGoodStep(typename G4InterpolationDriver< T, true >::StepperIterator it, field_utils::State &y, field_utils::State &dydx, G4double &hstep, G4double epsStep, G4double curveLength, G4FieldTrack *track) override
const G4QSSDriver & operator=(const G4QSSDriver &)=delete
G4double[G4FieldTrack::ncompSVEC] State