Geant4 11.1.1
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, Google Summer of Code 2022
33// Supervision: John Apostolakis,Renee Fatemi, Soon Yung Jun
34// --------------------------------------------------------------------
35#ifndef G4BORIS_SCHEME_HH
36#define G4BORIS_SCHEME_HH
37
39
40#include "G4Types.hh"
41
42
43// class G4EqMagElectricField;
44
45// #include "G4FieldTrack.hh"
46
48
50{
51 public:
52
53 G4BorisScheme() = default;
54 G4BorisScheme( // G4EqMagElectricField
55 G4EquationOfMotion* equation,
56 G4int nvar = 6);
57 ~G4BorisScheme() = default;
58
59 void DoStep( G4double restMass, G4double charge, const G4double yIn[],
60 G4double yOut[], G4double hstep) const;
61
62 protected:
63 // Used to implement the 'DoStep' method above
64 void UpdatePosition(const G4double restMass, const G4double charge, const G4double yIn[],
65 G4double yOut[], G4double hstep) const;
66
67 void UpdateVelocity(const G4double restMass, const G4double charge, const G4double yIn[],
68 G4double yOut[], G4double hstep) const;
69
70 public:
71 // - Methods using the Boris Scheme Stepping to estimate integration error
72 void StepWithErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep,
73 G4double yOut[], G4double yErr[]) const;
74 // Use two half-steps (comparing to a full step) to obtain output and error estimate
75
76 void StepWithMidAndErrorEstimate(const G4double yIn[], G4double restMass, G4double charge, G4double hstep,
77 G4double yMid[], G4double yOut[], G4double yErr[]) const;
78 // Same, and also return mid-point evaluation
79
80 // Auxiliary method
82 // inline void SetEquationOfMotion(G4EquationOfMotion* equation); // Un-needed, dangerous
83
85
86 private:
87
88 void copy(G4double dst[], const G4double src[]) const;
89
90 private:
91
92 G4EquationOfMotion* fEquation = nullptr;
93 G4int fnvar = 8;
94 static constexpr G4double c_l = CLHEP::c_light/CLHEP::m*CLHEP::second;
95};
96
97#include "G4BorisScheme.icc"
98#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
~G4BorisScheme()=default
G4int GetNumberOfVariables() const
void UpdateVelocity(const G4double restMass, const G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
void DoStep(G4double restMass, G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) const
G4BorisScheme()=default
G4EquationOfMotion * GetEquationOfMotion()
void UpdatePosition(const G4double restMass, const G4double charge, const G4double yIn[], G4double yOut[], G4double hstep) 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