Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VTransitionRadiation Class Reference

#include <G4VTransitionRadiation.hh>

+ Inheritance diagram for G4VTransitionRadiation:

Public Member Functions

 G4VTransitionRadiation (const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VTransitionRadiation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step) override
 
virtual void PrintInfoDefinition ()
 
void SetRegion (const G4Region *reg)
 
void SetModel (G4VTRModel *m)
 
void Clear ()
 
G4VTransitionRadiationoperator= (const G4VTransitionRadiation &right)
 
 G4VTransitionRadiation (const G4VTransitionRadiation &)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Public Attributes

G4LossTableManagertheManager
 
std::vector< const G4Material * > materials
 
std::vector< G4doublesteps
 
std::vector< G4ThreeVectornormals
 
G4ThreeVector startingPosition
 
G4ThreeVector startingDirection
 
const G4Regionregion
 
G4VTRModelmodel
 
G4int nSteps
 
G4double gammaMin
 
G4double cosDThetaMax
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 52 of file G4VTransitionRadiation.hh.

Constructor & Destructor Documentation

◆ G4VTransitionRadiation() [1/2]

G4VTransitionRadiation::G4VTransitionRadiation ( const G4String processName = "TR",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 48 of file G4VTransitionRadiation.cc.

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}
@ fTransitionRadiation
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4LossTableManager * theManager

◆ ~G4VTransitionRadiation()

G4VTransitionRadiation::~G4VTransitionRadiation ( )
virtual

Definition at line 65 of file G4VTransitionRadiation.cc.

66{
67 Clear();
69}
void DeRegister(G4VEnergyLossProcess *p)

◆ G4VTransitionRadiation() [2/2]

G4VTransitionRadiation::G4VTransitionRadiation ( const G4VTransitionRadiation )

Member Function Documentation

◆ Clear()

void G4VTransitionRadiation::Clear ( )

Definition at line 73 of file G4VTransitionRadiation.cc.

74{
75 materials.clear();
76 steps.clear();
77 normals.clear();
78 nSteps = 0;
79}
std::vector< G4ThreeVector > normals
std::vector< G4double > steps
std::vector< const G4Material * > materials

Referenced by G4VTransitionRadiation(), PostStepDoIt(), and ~G4VTransitionRadiation().

◆ GetMeanFreePath()

G4double G4VTransitionRadiation::GetMeanFreePath ( const G4Track track,
G4double  ,
G4ForceCondition condition 
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 174 of file G4VTransitionRadiation.cc.

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}
G4double condition(const G4ErrorSymMatrix &m)
@ StronglyForced
@ NotForced
G4Region * GetRegion() const
G4VPhysicalVolume * GetVolume() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4LogicalVolume * GetLogicalVolume() const
#define DBL_MAX
Definition: templates.hh:62

◆ IsApplicable()

G4bool G4VTransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 144 of file G4VTransitionRadiation.cc.

146{
147 return ( aParticle.GetPDGCharge() != 0.0 );
148}

◆ operator=()

G4VTransitionRadiation & G4VTransitionRadiation::operator= ( const G4VTransitionRadiation right)

◆ PostStepDoIt()

G4VParticleChange * G4VTransitionRadiation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 83 of file G4VTransitionRadiation.cc.

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}
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
double z() const
double x() const
double y() const
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4TrackStatus GetTrackStatus() const
G4Material * GetMaterial() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
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)

◆ PrintInfoDefinition()

void G4VTransitionRadiation::PrintInfoDefinition ( )
virtual

Definition at line 167 of file G4VTransitionRadiation.cc.

168{
169 if(model) model->PrintInfo();
170}
virtual void PrintInfo()
Definition: G4VTRModel.hh:70

◆ SetModel()

void G4VTransitionRadiation::SetModel ( G4VTRModel m)

Definition at line 160 of file G4VTransitionRadiation.cc.

161{
162 model = mod;
163}

◆ SetRegion()

void G4VTransitionRadiation::SetRegion ( const G4Region reg)

Definition at line 153 of file G4VTransitionRadiation.cc.

154{
155 region = reg;
156}

Member Data Documentation

◆ cosDThetaMax

G4double G4VTransitionRadiation::cosDThetaMax

Definition at line 103 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt().

◆ gammaMin

G4double G4VTransitionRadiation::gammaMin

Definition at line 102 of file G4VTransitionRadiation.hh.

Referenced by GetMeanFreePath().

◆ materials

std::vector<const G4Material*> G4VTransitionRadiation::materials

Definition at line 91 of file G4VTransitionRadiation.hh.

Referenced by Clear(), and PostStepDoIt().

◆ model

G4VTRModel* G4VTransitionRadiation::model

Definition at line 98 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt(), PrintInfoDefinition(), and SetModel().

◆ normals

std::vector<G4ThreeVector> G4VTransitionRadiation::normals

Definition at line 93 of file G4VTransitionRadiation.hh.

Referenced by Clear(), and PostStepDoIt().

◆ nSteps

G4int G4VTransitionRadiation::nSteps

Definition at line 100 of file G4VTransitionRadiation.hh.

Referenced by Clear(), GetMeanFreePath(), and PostStepDoIt().

◆ region

const G4Region* G4VTransitionRadiation::region

Definition at line 97 of file G4VTransitionRadiation.hh.

Referenced by GetMeanFreePath(), PostStepDoIt(), and SetRegion().

◆ startingDirection

G4ThreeVector G4VTransitionRadiation::startingDirection

Definition at line 96 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt().

◆ startingPosition

G4ThreeVector G4VTransitionRadiation::startingPosition

Definition at line 95 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt().

◆ steps

std::vector<G4double> G4VTransitionRadiation::steps

Definition at line 92 of file G4VTransitionRadiation.hh.

Referenced by Clear(), and PostStepDoIt().

◆ theManager

G4LossTableManager* G4VTransitionRadiation::theManager

Definition at line 89 of file G4VTransitionRadiation.hh.

Referenced by G4VTransitionRadiation(), and ~G4VTransitionRadiation().


The documentation for this class was generated from the following files: