Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Channeling.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
28#ifndef G4Channeling_h
29#define G4Channeling_h 1
30
31#include "G4ios.hh"
32#include "globals.hh"
33#include "G4VDiscreteProcess.hh"
35#include "G4ExtendedMaterial.hh"
37
39
41{
42public:
43
45 virtual ~G4Channeling();
47 const G4Step&);
49 return(aPD.GetPDGCharge() != 0.);
50 };
51 virtual void BuildPhysicsTable(const G4ParticleDefinition&){;};
52
53protected:
54 virtual G4double GetMeanFreePath(const G4Track&,
57
58private:
59 G4ParticleDefinition* GetParticleDefinition(const G4Track& aTrack){
60 return const_cast<G4ParticleDefinition*>(aTrack.GetParticleDefinition());
61 }
62private:
63 G4StepPoint* GetPre(const G4Track& aTrack){return aTrack.GetStep()->GetPreStepPoint();}
64 G4StepPoint* GetPost(const G4Track& aTrack){return aTrack.GetStep()->GetPostStepPoint();}
65
66
67private:
68 G4ChannelingMaterialData* GetMatData(const G4Track& aTrack){
70 if(aLV->IsExtended() == true){
72 return (G4ChannelingMaterialData*) aEM->RetrieveExtension("channeling");
73 }
74 else{
75 return nullptr;
76 }
77 }
78
79 //----------------------------------------
80 // Functions for the calculations of
81 // parameters related to channeling
82 //----------------------------------------
83public:
85 return std::sqrt(2.0*GetMatData(aTrack)->GetPot()->GetMaxMin()
86 /GetPre(aTrack)->GetTotalEnergy());}
88 return (CLHEP::pi * GetMatData(aTrack)->GetPot()->GetIntSp(0)
89 / GetCriticalAngle(aTrack));
90 }
91 //----------------------------------------
92 // Channeling Auxiliary Track Information
93 //----------------------------------------
94private:
95 G4int fChannelingID;
96 G4ChannelingTrackData* GetTrackData(const G4Track&);
97
98 //----------------------------------------
99 // Variables for the integration
100 // of the particle trajectory
101 //----------------------------------------
102private:
103 G4bool UpdateIntegrationStep(const G4Track&,
105 G4double&);
106 G4bool UpdateParameters(const G4Track&);
107
108 void GetEF(const G4Track&,G4ThreeVector&,G4ThreeVector&);
109
110public:
112
113public:
114 G4double GetTransverseVariationMax() {return fTransverseVariationMax;};
115 void SetTransverseVariationMax(G4double aDouble) {fTransverseVariationMax = aDouble;};
116
117 G4double GetTimeStepMin() {return fTimeStepMin;};
118 void SetTimeStepMin(G4double aDouble) {fTimeStepMin = aDouble;};
119
120private:
121 G4double fTimeStepMin;
122 G4double fTimeStepMax;
123
124 G4double fTransverseVariationMax;
125
126 const G4ThreeVector k010;
127 G4ThreeVector fSpin;
128};
129
130#endif
131
132
133
134
135
136
137
138
139
140
G4ForceCondition
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double GetTimeStepMin()
G4double GetCriticalAngle(const G4Track &aTrack)
Definition: G4Channeling.hh:84
G4double GetTransverseVariationMax()
void PosToLattice(G4StepPoint *step, G4ThreeVector &)
Definition: G4Channeling.cc:75
virtual ~G4Channeling()
Definition: G4Channeling.cc:49
void SetTransverseVariationMax(G4double aDouble)
virtual G4bool IsApplicable(const G4ParticleDefinition &aPD)
Definition: G4Channeling.hh:48
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4Channeling.hh:51
G4double GetOscillationPeriod(const G4Track &aTrack)
Definition: G4Channeling.hh:87
void SetTimeStepMin(G4double aDouble)
G4VMaterialExtension * RetrieveExtension(const G4String &name)
virtual G4bool IsExtended() const
G4Material * GetMaterial() const
G4double GetPDGCharge() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
const G4Step * GetStep() const
G4LogicalVolume * GetLogicalVolume() const