Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4RKPropagation.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
27#ifndef G4RKPropagation_h
28#define G4RKPropagation_h 1
29
31#include "G4VNuclearField.hh"
32#include "G4V3DNucleus.hh"
33#include "G4KM_DummyField.hh"
34#include "G4Mag_EqRhs.hh"
35#include <map>
36
37
39{
40
41public:
43 virtual ~G4RKPropagation();
44
45private:
46 G4RKPropagation(const G4RKPropagation &right);
47 const G4RKPropagation & operator=(const G4RKPropagation & right);
48 G4int operator==(const G4RKPropagation & right) const;
49 G4int operator!=(const G4RKPropagation & right) const;
50
51public:
52
53 virtual void Init(G4V3DNucleus * nucleus);
54 virtual void Transport(G4KineticTrackVector &theActive,
55 const G4KineticTrackVector &theSpectators,
56 G4double theTimeStep);
58 G4double & t1, G4double & t2);
60private:
61 G4double theOuterRadius;
62 G4V3DNucleus * theNucleus;
63 std::map <G4int, G4VNuclearField *, std::less<G4int> > * theFieldMap;
64 std::map <G4int, G4Mag_EqRhs *, std::less<G4int> > * theEquationMap;
65 G4KM_DummyField * theField;
66
67 G4ThreeVector theMomentumTranfer;
68
70 const G4ThreeVector & currentPos,
71 const G4LorentzVector & momentum,
72 G4double & t1, G4double & t2);
73// implementation
74
75 G4bool FieldTransport(G4KineticTrack * track, const G4double timestep);
76 G4bool FreeTransport(G4KineticTrack * track, const G4double timestep);
77
78 void delete_FieldsAndMap(
79 std::map <G4int, G4VNuclearField *, std::less<G4int> > * aMap);
80 void delete_EquationsAndMap(
81 std::map <G4int, G4Mag_EqRhs *, std::less<G4int> > * aMap);
82
83 public:
85 {
86 std::map <G4int, G4VNuclearField *, std::less<G4int> >::iterator iter;
87 iter = theFieldMap->find(encoding);
88 if(iter == theFieldMap->end()) return 0;
89 return (*theFieldMap)[encoding]->GetBarrier();
90 }
91
93 {
94 std::map <G4int, G4VNuclearField *, std::less<G4int> >::iterator iter;
95 iter = theFieldMap->find(encoding);
96 if(iter == theFieldMap->end()) return 0;
97 return (*theFieldMap)[encoding]->GetField(pos);
98 }
99
100};
101
102#endif
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetField(G4int encoding, G4ThreeVector pos)
virtual ~G4RKPropagation()
G4double GetBarrier(G4int encoding)
virtual void Transport(G4KineticTrackVector &theActive, const G4KineticTrackVector &theSpectators, G4double theTimeStep)
virtual void Init(G4V3DNucleus *nucleus)
G4ThreeVector GetMomentumTransfer() const
G4bool GetSphereIntersectionTimes(const G4KineticTrack *track, G4double &t1, G4double &t2)
#define encoding
Definition: xmlparse.cc:588