Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BorisScheme.hh
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25// G4BorisScheme
26//
27// Class description:
28//
29// Implementation of the Boris algorithm for advancing
30// charged particles in an electromagnetic field.
31
32// Author: Divyansh Tiwari (CERN, Google Summer of Code 2022), 05.11.2022
33// Supervision: John Apostolakis (CERN), Renee Fatemi, Soon Yung Jun (FNAL)
34// --------------------------------------------------------------------
35#ifndef G4BORIS_SCHEME_HH
36#define G4BORIS_SCHEME_HH
37
38#include "G4Types.hh"
39
41
43
44/**
45 * @brief The G4BorisScheme class implements of the Boris algorithm for
46 * advancing charged particles in an electromagnetic field.
47 */
48
50{
51 public:
52
53 /**
54 * Default Constructor.
55 */
56 G4BorisScheme() = default;
57
58 /**
59 * Constructor for the equation of motion.
60 * @param[in] equation Pointer to the equation of motion algorithm.
61 * @param[in] nvar The number of integration variables.
62 */
63 G4BorisScheme( G4EquationOfMotion* equation, G4int nvar = 6 );
64
65 /**
66 * Default Destructor.
67 */
68 ~G4BorisScheme() = default;
69
70 /**
71 * Does one step, updating velocity and position.
72 * @param[in] restMass Particle mass.
73 * @param[in] charge Particle charge.
74 * @param[in] yIn Initial position.
75 * @param[out] yOut Updated position.
76 * @param[in] hstep Proposed step.
77 */
78 void DoStep(G4double restMass, G4double charge, const G4double yIn[],
79 G4double yOut[], G4double hstep) const;
80
81 /**
82 * Adopts the Boris Scheme Stepping to estimate the integration error.
83 * Uses two half-steps (comparing to a full step) to obtain output and
84 * error estimate.
85 * @param[in] yIn Initial position.
86 * @param[in] restMass Particle mass.
87 * @param[in] charge Particle charge.
88 * @param[in] hstep Proposed step.
89 * @param[out] yOut Updated position.
90 * @param[out] yErr The estimated error.
91 */
92 void StepWithErrorEstimate(const G4double yIn[], G4double restMass,
93 G4double charge, G4double hstep,
94 G4double yOut[], G4double yErr[]) const;
95
96 /**
97 * Adopts the Boris Scheme Stepping to estimate the integration error.
98 * Uses two half-steps (comparing to a full step) to obtain output and
99 * error estimate. Same as above, but also returns the mid-point evaluation.
100 * @param[in] yIn Initial position.
101 * @param[in] restMass Particle mass.
102 * @param[in] charge Particle charge.
103 * @param[in] hstep Proposed step.
104 * @param[out] yMid tThe mid-point evaluation.
105 * @param[out] yOut Updated position.
106 * @param[out] yErr The estimated error.
107 */
108 void StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass,
109 G4double charge, G4double hstep,
110 G4double yMid[], G4double yOut[], G4double yErr[]) const;
111
112 /**
113 * Auxiliary methods returning a pointer to the equation of motion
114 * and the number of integration variables.
115 */
118
119 private:
120
121 /**
122 * Internal methods for updating position and velocity, used in DoStep().
123 */
124 void UpdatePosition(const G4double restMass, const G4double charge,
125 const G4double yIn[], G4double yOut[], G4double hstep) const;
126 void UpdateVelocity(const G4double restMass, const G4double charge,
127 const G4double yIn[], G4double yOut[], G4double hstep) const;
128
129 /**
130 * Utility to mem-copy 'src' array data to 'dst'.
131 */
132 void copy(G4double dst[], const G4double src[]) const;
133
134 private:
135
136 G4EquationOfMotion* fEquation = nullptr;
137 G4int fnvar = 8;
138 static constexpr G4double c_l = CLHEP::c_light/CLHEP::m*CLHEP::second;
139};
140
141#include "G4BorisScheme.icc"
142
143#endif
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
~G4BorisScheme()=default
G4int GetNumberOfVariables() const
void DoStep(G4double restMass, G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
G4BorisScheme()=default
G4EquationOfMotion * GetEquationOfMotion() const
void StepWithErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep, G4double yOut[], G4double yErr[]) const
void StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep, G4double yMid[], G4double yOut[], G4double yErr[]) const
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...