Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChordFinder.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// G4ChordFinder
27//
28// Class description:
29//
30// A class that provides RK integration of motion ODE (as does g4magtr)
31// and also has a method that returns an Approximate point on the curve
32// near to a (chord) point.
33
34// Author: John Apostolakis (CERN), 25.02.1997 - Design and implementation
35// -------------------------------------------------------------------
36#ifndef G4CHORDFINDER_HH
37#define G4CHORDFINDER_HH
38
40#include "G4FieldParameters.hh"
42
43#include <memory>
44
46
47class G4MagneticField;
49class G4HelixHeum;
50class G4QSStepper;
51
52/**
53 * @brief G4ChordFinder is a class that provides Runge-Kutta integration of
54 * motion ODE and also has a method that returns an approximate point on the
55 * curve near to a (chord) point.
56 */
57
59{
60 public:
61
65
66 /**
67 * The most flexible constructor, which allows the user to specify
68 * any type of field, equation, stepper and integration driver.
69 * @param[in] pIntegrationDriver Pointer to the integrator driver to use.
70 */
71 explicit G4ChordFinder( G4VIntegrationDriver* pIntegrationDriver );
72
73 /**
74 * Constructor that creates defaults for all "children" classes.
75 * The type of equation of motion is fixed.
76 * A default type of stepper (Dormand Prince since release 10.4) is used,
77 * and the corresponding integration driver.
78 * @param[in] itsMagField Pointer to the magnetic field.
79 * @param[in] stepMinimum Pointer to the magnetic field.
80 * @param[in] pItsStepper Optional pointer to the stepper algorithm.
81 * @param[in] stepperDriverChoice Type of stepper driver.
82 */
83 G4ChordFinder( G4MagneticField* itsMagField,
85 G4MagIntegratorStepper* pItsStepper = nullptr,
86 G4int stepperDriverChoice = kTemplatedStepperType );
87
88 /**
89 * Destructor.
90 */
92
93 /**
94 * Copy constructor and assignment operator not allowed.
95 */
96 G4ChordFinder(const G4ChordFinder&) = delete;
98
99 /**
100 * Computes the step to take, based on chord limits.
101 * Uses ODE solver's driver to find the endpoint that satisfies
102 * the chord criterion that: d_chord < delta_chord.
103 * @param[in,out] yCurrent The current track in field.
104 * @param[in] stepInitial Proposed initial step length.
105 * @param[in] epsStep_Relative Requested accuracy.
106 * @param[in] latestSafetyOrigin Last safety origin point. Unused.
107 * @param[in] lasestSafetyRadius Last safety distance. Unused.
108 * @returns The length of step taken.
109 */
111 G4double stepInitial,
112 G4double epsStep_Relative,
113 const G4ThreeVector& latestSafetyOrigin,
114 G4double lasestSafetyRadius );
115
116 /**
117 * Uses the Brent algorithm when possible, to determine the closest point
118 * on the curve. Given a starting curve point A (CurveA_PointVelocity),
119 * curve point B (CurveB_PointVelocity), a point E which is (generally)
120 * not on the curve and a point F which is on the curve (first
121 * approximation), find new point S on the curve closer to point E.
122 * While advancing towards S utilise 'eps_step' as a measure of the
123 * relative accuracy of each Step.
124 * @returns The end point on the curve closer to the given point E.
125 */
126 G4FieldTrack ApproxCurvePointS( const G4FieldTrack& curveAPointVelocity,
127 const G4FieldTrack& curveBPointVelocity,
128 const G4FieldTrack& ApproxCurveV,
129 const G4ThreeVector& currentEPoint,
130 const G4ThreeVector& currentFPoint,
131 const G4ThreeVector& PointG,
132 G4bool first, G4double epsStep );
133
134 /**
135 * If r=|AE|/|AB|, and s=true path lenght (AB)
136 * returns the point that is r*s along the curve.
137 */
138 G4FieldTrack ApproxCurvePointV( const G4FieldTrack& curveAPointVelocity,
139 const G4FieldTrack& curveBPointVelocity,
140 const G4ThreeVector& currentEPoint,
141 G4double epsStep);
142
143 /**
144 * Calculates the inverse parabolic through the three points (x,y) and
145 * returns the value x that, for the inverse parabolic, corresponds to y=0.
146 */
147 inline G4double InvParabolic( const G4double xa, const G4double ya,
148 const G4double xb, const G4double yb,
149 const G4double xc, const G4double yc );
150
151 /**
152 * Accessors and modifiers.
153 */
154 inline G4double GetDeltaChord() const;
155 inline void SetDeltaChord(G4double newval);
156 inline void SetIntegrationDriver(G4VIntegrationDriver* IntegrationDriver);
158
159 /**
160 * Clears the internal state (last step estimate).
161 */
162 inline void ResetStepEstimate();
163
164 /**
165 * Sets the verbosity.
166 * @returns The old verbosity value.
167 */
168 inline G4int SetVerbose( G4int newvalue=1 );
169
170 /**
171 * Dispatch interface method for computing step.
172 */
173 inline void OnComputeStep(const G4FieldTrack* track);
174
175 /**
176 * Writes out to stream the parameters/state of the driver.
177 */
178 friend std::ostream& operator<<( std::ostream& os, const G4ChordFinder& cf);
179
180 /**
181 * Sets verbosity for constructor.
182 */
183 static void SetVerboseConstruction(G4bool v = true);
184
185 private: // ............................................................
186
187 static G4bool gVerboseCtor; // Verbosity for contructor
188
189 // Constants
190 // ---------------------
191 const G4double fDefaultDeltaChord = G4FieldDefaults::kDeltaChord;
192
193 // PARAMETERS
194 // ---------------------
195 G4double fDeltaChord; // Maximum miss distance
196
197 G4int fStatsVerbose = 0; // if > 0, print Statistics in destructor
198
199 // DEPENDENT Objects
200 // ---------------------
201 G4VIntegrationDriver* fIntgrDriver = nullptr;
202 G4MagIntegratorStepper* fRegularStepperOwned = nullptr;
203 G4MagIntegratorStepper* fNewFSALStepperOwned = nullptr;
204 std::unique_ptr<G4HelixHeum> fLongStepper;
205 G4CachedMagneticField* fCachedField = nullptr;
206 G4QSStepper* fQssStepperOwned = nullptr;
207 G4EquationOfMotion* fEquation = nullptr;
208};
209
210// Inline function implementation:
211
212#include "G4ChordFinder.icc"
213
214#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4CachedMagneticField is a specialisation of G4MagneticField and is used to cache the Magnetic Field ...
G4double InvParabolic(const G4double xa, const G4double ya, const G4double xb, const G4double yb, const G4double xc, const G4double yc)
static void SetVerboseConstruction(G4bool v=true)
void SetIntegrationDriver(G4VIntegrationDriver *IntegrationDriver)
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
G4FieldTrack ApproxCurvePointS(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4FieldTrack &ApproxCurveV, const G4ThreeVector &currentEPoint, const G4ThreeVector &currentFPoint, const G4ThreeVector &PointG, G4bool first, G4double epsStep)
G4ChordFinder(const G4ChordFinder &)=delete
G4double GetDeltaChord() const
G4VIntegrationDriver * GetIntegrationDriver()
G4int SetVerbose(G4int newvalue=1)
void ResetStepEstimate()
void OnComputeStep(const G4FieldTrack *track)
void SetDeltaChord(G4double newval)
G4double AdvanceChordLimited(G4FieldTrack &yCurrent, G4double stepInitial, G4double epsStep_Relative, const G4ThreeVector &latestSafetyOrigin, G4double lasestSafetyRadius)
G4ChordFinder & operator=(const G4ChordFinder &)=delete
friend std::ostream & operator<<(std::ostream &os, const G4ChordFinder &cf)
G4ChordFinder(G4VIntegrationDriver *pIntegrationDriver)
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
G4HelixHeum implements a simple Heum stepper for magnetic field with 3rd order solver.
G4MagIntegratorStepper is an abstract base class for integrator of particle's equation of motion,...
G4QSStepper is an integrator of particle's equation of motion based on the QSS implementation.
G4VFSALIntegrationStepper is a class similar to G4VMagIntegratorStepper, but for steppers which estim...
constexpr G4double kDeltaChord
Default delta chord in G4ChordFinder.
constexpr G4double kMinimumStep
Default minimum step in G4ChordFinder.