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
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// GEANT 4 class implementation file --- Copyright CERN 1995
28
29// History:
30// 29.02.04 V.Ivanchenko create
31// 28.07.05, P.Gumplinger add G4ProcessType to constructor
32
34
35#include "G4EmProcessSubType.hh"
36#include "G4LossTableManager.hh"
37#include "G4Material.hh"
39#include "G4Region.hh"
41#include "G4VTRModel.hh"
42
43///////////////////////////////////////////////////////////////////////
45 G4ProcessType type)
46 : G4VDiscreteProcess(processName, type)
47 , region(nullptr)
48 , model(nullptr)
49 , gammaMin(100.)
50 , cosDThetaMax(std::cos(0.1))
51 , nSteps(0)
52{
54 Clear();
55 theManager = G4LossTableManager::Instance();
56 theManager->Register(this);
57}
58
59///////////////////////////////////////////////////////////////////////
61{
62 Clear();
63 theManager->DeRegister(this);
64}
65
66void G4VTransitionRadiation::ProcessDescription(std::ostream& out) const
67{
68 out << "Generic process of transition radiation.\n";
69
70 if(model)
71 model->PrintInfo();
72}
73
74///////////////////////////////////////////////////////////////////////
76{
77 materials.clear();
78 steps.clear();
79 normals.clear();
80 nSteps = 0;
81}
82
83///////////////////////////////////////////////////////////////////////
85 const G4Step& step)
86{
87 // Fill temporary vectors
88 const G4Material* material = track.GetMaterial();
89 G4double length = step.GetStepLength();
90 G4ThreeVector direction = track.GetMomentumDirection();
91
92 if(nSteps == 0)
93 {
94 nSteps = 1;
95 materials.push_back(material);
96 steps.push_back(length);
97 const G4StepPoint* point = step.GetPreStepPoint();
98 startingPosition = point->GetPosition();
99 startingDirection = point->GetMomentumDirection();
100 G4bool valid = true;
103 ->GetLocalExitNormal(&valid);
104 if(valid)
105 normals.push_back(n);
106 else
107 normals.push_back(direction);
108 }
109 else
110 {
111 if(material == materials[nSteps - 1])
112 {
113 steps[nSteps - 1] += length;
114 }
115 else
116 {
117 ++nSteps;
118 materials.push_back(material);
119 steps.push_back(length);
120 G4bool valid = true;
123 ->GetLocalExitNormal(&valid);
124 if(valid)
125 normals.push_back(n);
126 else
127 normals.push_back(direction);
128 }
129 }
130
131 // Check PostStepPoint condition
132 if(track.GetTrackStatus() == fStopAndKill ||
133 track.GetVolume()->GetLogicalVolume()->GetRegion() != region ||
134 startingDirection.x() * direction.x() +
135 startingDirection.y() * direction.y() +
136 startingDirection.z() * direction.z() <
137 cosDThetaMax)
138 {
139 if(model)
140 {
141 model->GenerateSecondaries(*pParticleChange, materials, steps, normals,
142 startingPosition, track);
143 }
144 Clear();
145 }
146
147 return pParticleChange;
148}
149
150///////////////////////////////////////////////////////////////////////
152 const G4ParticleDefinition& aParticle)
153{
154 return (aParticle.GetPDGCharge() != 0.0);
155}
156
157///////////////////////////////////////////////////////////////////////
158void G4VTransitionRadiation::SetRegion(const G4Region* reg) { region = reg; }
159
160///////////////////////////////////////////////////////////////////////
162
163///////////////////////////////////////////////////////////////////////
166{
167 if(nSteps > 0)
168 {
170 }
171 else
172 {
174 if(track.GetKineticEnergy() / track.GetDefinition()->GetPDGMass() + 1.0 >
175 gammaMin &&
176 track.GetVolume()->GetLogicalVolume()->GetRegion() == region)
177 {
179 }
180 }
181 return DBL_MAX; // so TR doesn't limit mean free path
182}
@ 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:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
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:68
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
virtual G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4VTransitionRadiation(const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
void SetRegion(const G4Region *reg)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step) override
void ProcessDescription(std::ostream &) const override
#define DBL_MAX
Definition: templates.hh:62