Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VTransitionRadiation.cc
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// G4VTransitionRadiation class -- implementation file
29
30// GEANT 4 class implementation file --- Copyright CERN 1995
31// CERN Geneva Switzerland
32
33// History:
34// 29.02.04 V.Ivanchenko create
35// 28.07.05, P.Gumplinger add G4ProcessType to constructor
36
39#include "G4VTRModel.hh"
40#include "G4Material.hh"
41#include "G4Region.hh"
43#include "G4EmProcessSubType.hh"
44#include "G4LossTableManager.hh"
45
46///////////////////////////////////////////////////////////////////////
47
49 G4ProcessType type )
50 : G4VDiscreteProcess(processName, type),
51 region(nullptr),
52 model(nullptr),
53 nSteps(0),
54 gammaMin(100.),
55 cosDThetaMax(std::cos(0.1))
56{
58 Clear();
60 theManager->Register(this);
61}
62
63///////////////////////////////////////////////////////////////////////
64
66{
67 Clear();
69}
70
71///////////////////////////////////////////////////////////////////////
72
74{
75 materials.clear();
76 steps.clear();
77 normals.clear();
78 nSteps = 0;
79}
80
81///////////////////////////////////////////////////////////////////////
82
84 const G4Track& track,
85 const G4Step& step)
86{
87
88 // Fill temporary vectors
89
90 const G4Material* material = track.GetMaterial();
91 G4double length = step.GetStepLength();
92 G4ThreeVector direction = track.GetMomentumDirection();
93
94 if(nSteps == 0) {
95
96 nSteps = 1;
97 materials.push_back(material);
98 steps.push_back(length);
99 const G4StepPoint* point = step.GetPreStepPoint();
100 startingPosition = point->GetPosition();
102 G4bool valid = true;
105 if(valid) normals.push_back(n);
106 else normals.push_back(direction);
107
108 } else {
109
110 if(material == materials[nSteps-1]) {
111 steps[nSteps-1] += length;
112 } else {
113 nSteps++;
114 materials.push_back(material);
115 steps.push_back(length);
116 G4bool valid = true;
119 if(valid) normals.push_back(n);
120 else normals.push_back(direction);
121 }
122 }
123
124 // Check POstStepPoint condition
125
126 if(track.GetTrackStatus() == fStopAndKill ||
127 track.GetVolume()->GetLogicalVolume()->GetRegion() != region ||
128 startingDirection.x()*direction.x() +
129 startingDirection.y()*direction.y() +
130 startingDirection.z()*direction.z() < cosDThetaMax)
131 {
132 if(model) {
134 normals, startingPosition, track);
135 }
136 Clear();
137 }
138
139 return pParticleChange;
140}
141
142///////////////////////////////////////////////////////////////////////
143
145 const G4ParticleDefinition& aParticle)
146{
147 return ( aParticle.GetPDGCharge() != 0.0 );
148}
149
150///////////////////////////////////////////////////////////////////////
151
152
154{
155 region = reg;
156}
157
158///////////////////////////////////////////////////////////////////////
159
161{
162 model = mod;
163}
164
165///////////////////////////////////////////////////////////////////////
166
168{
169 if(model) model->PrintInfo();
170}
171
172///////////////////////////////////////////////////////////////////////
173
175 const G4Track& track, G4double,
177{
178 if(nSteps > 0) {
180 } else {
182 if(track.GetKineticEnergy()/track.GetDefinition()->GetPDGMass() + 1.0 > gammaMin &&
183 track.GetVolume()->GetLogicalVolume()->GetRegion() == region) {
185 }
186 }
187 return DBL_MAX; // so TR doesn't limit mean free path
188}
189
190///////////////////////////////////////////////////////////////////////
@ fTransitionRadiation
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ StronglyForced
@ NotForced
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
double z() const
double x() const
double y() const
G4Region * GetRegion() const
static G4LossTableManager * Instance()
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
G4double GetPDGCharge() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() const
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
virtual void GenerateSecondaries(G4VParticleChange &pChange, std::vector< const G4Material * > &materials, std::vector< G4double > &steps, std::vector< G4ThreeVector > &normals, G4ThreeVector &startingPosition, const G4Track &track)
virtual void PrintInfo()
Definition: G4VTRModel.hh:70
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
std::vector< G4ThreeVector > normals
virtual G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4VTransitionRadiation(const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
std::vector< G4double > steps
void SetRegion(const G4Region *reg)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step) override
std::vector< const G4Material * > materials
G4LossTableManager * theManager
#define DBL_MAX
Definition: templates.hh:62