Geant4 11.1.1
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
 
void ProcessDescription (std::ostream &) const override
 
void DumpInfo () const override
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step) override
 
void SetRegion (const G4Region *reg)
 
void SetModel (G4VTRModel *m)
 
void Clear ()
 
G4VTransitionRadiationoperator= (const G4VTransitionRadiation &right)=delete
 
 G4VTransitionRadiation (const G4VTransitionRadiation &)=delete
 
- 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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- 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 const G4VProcessGetCreatorProcess () const
 
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 &)
 

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 44 of file G4VTransitionRadiation.cc.

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

◆ ~G4VTransitionRadiation()

G4VTransitionRadiation::~G4VTransitionRadiation ( )
virtual

Definition at line 60 of file G4VTransitionRadiation.cc.

61{
62 Clear();
63 theManager->DeRegister(this);
64}
void DeRegister(G4VEnergyLossProcess *p)

◆ G4VTransitionRadiation() [2/2]

G4VTransitionRadiation::G4VTransitionRadiation ( const G4VTransitionRadiation )
delete

Member Function Documentation

◆ Clear()

void G4VTransitionRadiation::Clear ( )

Definition at line 75 of file G4VTransitionRadiation.cc.

76{
77 materials.clear();
78 steps.clear();
79 normals.clear();
80 nSteps = 0;
81}

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

◆ DumpInfo()

void G4VTransitionRadiation::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 66 of file G4VTransitionRadiation.hh.

G4GLOB_DLL std::ostream G4cout
void ProcessDescription(std::ostream &) const override

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 164 of file G4VTransitionRadiation.cc.

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}
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 151 of file G4VTransitionRadiation.cc.

153{
154 return (aParticle.GetPDGCharge() != 0.0);
155}

◆ operator=()

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

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 84 of file G4VTransitionRadiation.cc.

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}
@ 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:325
virtual void GenerateSecondaries(G4VParticleChange &pChange, std::vector< const G4Material * > &materials, std::vector< G4double > &steps, std::vector< G4ThreeVector > &normals, G4ThreeVector &startingPosition, const G4Track &track)

◆ ProcessDescription()

void G4VTransitionRadiation::ProcessDescription ( std::ostream &  out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 66 of file G4VTransitionRadiation.cc.

67{
68 out << "Generic process of transition radiation.\n";
69
70 if(model)
71 model->PrintInfo();
72}
virtual void PrintInfo()
Definition: G4VTRModel.hh:68

Referenced by DumpInfo().

◆ SetModel()

void G4VTransitionRadiation::SetModel ( G4VTRModel m)

Definition at line 161 of file G4VTransitionRadiation.cc.

161{ model = mod; }

◆ SetRegion()

void G4VTransitionRadiation::SetRegion ( const G4Region reg)

Definition at line 158 of file G4VTransitionRadiation.cc.

158{ region = reg; }

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