Geant4 9.6.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// $Id$
28//
29// G4VTransitionRadiation class -- implementation file
30
31// GEANT 4 class implementation file --- Copyright CERN 1995
32// CERN Geneva Switzerland
33
34// History:
35// 29.02.04 V.Ivanchenko create
36// 28.07.05, P.Gumplinger add G4ProcessType to constructor
37
40#include "G4VTRModel.hh"
41#include "G4Material.hh"
42#include "G4Region.hh"
44
45///////////////////////////////////////////////////////////////////////
46
48 G4ProcessType type )
49 : G4VDiscreteProcess(processName, type),
50 region(0),
51 model(0),
52 nSteps(0),
53 gammaMin(100),
54 cosDThetaMax(std::cos(0.1))
55{
56 Clear();
57}
58
59///////////////////////////////////////////////////////////////////////
60
62{
63 Clear();
64}
65
66///////////////////////////////////////////////////////////////////////
67
69{
70 materials.clear();
71 steps.clear();
72 normals.clear();
73 nSteps = 0;
74}
75
76///////////////////////////////////////////////////////////////////////
77
79 const G4Track& track,
80 const G4Step& step)
81{
82
83 // Fill temporary vectors
84
85 const G4Material* material = track.GetMaterial();
86 G4double length = step.GetStepLength();
87 G4ThreeVector direction = track.GetMomentumDirection();
88
89 if(nSteps == 0) {
90
91 nSteps = 1;
92 materials.push_back(material);
93 steps.push_back(length);
94 const G4StepPoint* point = step.GetPreStepPoint();
97 G4bool valid = true;
100 if(valid) normals.push_back(n);
101 else normals.push_back(direction);
102
103 } else {
104
105 if(material == materials[nSteps-1]) {
106 steps[nSteps-1] += length;
107 } else {
108 nSteps++;
109 materials.push_back(material);
110 steps.push_back(length);
111 G4bool valid = true;
114 if(valid) normals.push_back(n);
115 else normals.push_back(direction);
116 }
117 }
118
119 // Check POstStepPoint condition
120
121 if(track.GetTrackStatus() == fStopAndKill ||
122 track.GetVolume()->GetLogicalVolume()->GetRegion() != region ||
123 startingDirection.x()*direction.x() +
124 startingDirection.y()*direction.y() +
125 startingDirection.z()*direction.z() < cosDThetaMax)
126 {
127 if(model) {
129 normals, startingPosition, track);
130 }
131 Clear();
132 }
133
134 return pParticleChange;
135}
136
137///////////////////////////////////////////////////////////////////////
138
140 const G4ParticleDefinition& aParticle)
141{
142 return ( aParticle.GetPDGCharge() != 0.0 );
143}
144
145///////////////////////////////////////////////////////////////////////
146
147
149{
150 region = reg;
151}
152
153///////////////////////////////////////////////////////////////////////
154
156{
157 model = mod;
158}
159
160///////////////////////////////////////////////////////////////////////
161
163{
164 if(model) model->PrintInfo();
165}
166
167///////////////////////////////////////////////////////////////////////
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
double z() const
double x() const
double y() const
G4Region * GetRegion() const
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
G4double GetPDGCharge() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
Definition: G4Step.hh:78
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() const
G4Material * GetMaterial() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4LogicalVolume * GetLogicalVolume() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
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:71
std::vector< G4ThreeVector > normals
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)
std::vector< const G4Material * > materials
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType)