Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VTransitionRadiation Class Reference

#include <G4VTransitionRadiation.hh>

Inheritance diagram for G4VTransitionRadiation:
Collaboration 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)=delete
 
 G4VTransitionRadiation (const G4VTransitionRadiation &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
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 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
 
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

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)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 52 of file G4VTransitionRadiation.hh.

Constructor & Destructor Documentation

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

Definition at line 48 of file G4VTransitionRadiation.cc.

Here is the call graph for this function:

G4VTransitionRadiation::~G4VTransitionRadiation ( )
virtual

Definition at line 63 of file G4VTransitionRadiation.cc.

64 {
65  Clear();
66 }

Here is the call graph for this function:

G4VTransitionRadiation::G4VTransitionRadiation ( const G4VTransitionRadiation )
delete

Member Function Documentation

void G4VTransitionRadiation::Clear ( )

Definition at line 70 of file G4VTransitionRadiation.cc.

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

Here is the caller graph for this function:

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

Implements G4VDiscreteProcess.

Definition at line 105 of file G4VTransitionRadiation.hh.

108 {
109  if(nSteps > 0) {
111  } else {
112  *condition = NotForced;
113  if(track.GetKineticEnergy()/track.GetDefinition()->GetPDGMass() + 1.0 > gammaMin &&
114  track.GetVolume()->GetLogicalVolume()->GetRegion() == region) {
116  }
117  }
118  return DBL_MAX; // so TR doesn't limit mean free path
119 }
G4double condition(const G4ErrorSymMatrix &m)
G4ParticleDefinition * GetDefinition() const
G4Region * GetRegion() const
G4double GetKineticEnergy() const
G4LogicalVolume * GetLogicalVolume() const
G4double GetPDGMass() const
G4VPhysicalVolume * GetVolume() const
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

G4bool G4VTransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 141 of file G4VTransitionRadiation.cc.

143 {
144  return ( aParticle.GetPDGCharge() != 0.0 );
145 }

Here is the call graph for this function:

G4VTransitionRadiation& G4VTransitionRadiation::operator= ( const G4VTransitionRadiation right)
delete
G4VParticleChange * G4VTransitionRadiation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 80 of file G4VTransitionRadiation.cc.

83 {
84 
85  // Fill temporary vectors
86 
87  const G4Material* material = track.GetMaterial();
88  G4double length = step.GetStepLength();
89  G4ThreeVector direction = track.GetMomentumDirection();
90 
91  if(nSteps == 0) {
92 
93  nSteps = 1;
94  materials.push_back(material);
95  steps.push_back(length);
96  const G4StepPoint* point = step.GetPreStepPoint();
97  startingPosition = point->GetPosition();
99  G4bool valid = true;
102  if(valid) normals.push_back(n);
103  else normals.push_back(direction);
104 
105  } else {
106 
107  if(material == materials[nSteps-1]) {
108  steps[nSteps-1] += length;
109  } else {
110  nSteps++;
111  materials.push_back(material);
112  steps.push_back(length);
113  G4bool valid = true;
116  if(valid) normals.push_back(n);
117  else normals.push_back(direction);
118  }
119  }
120 
121  // Check POstStepPoint condition
122 
123  if(track.GetTrackStatus() == fStopAndKill ||
124  track.GetVolume()->GetLogicalVolume()->GetRegion() != region ||
125  startingDirection.x()*direction.x() +
126  startingDirection.y()*direction.y() +
127  startingDirection.z()*direction.z() < cosDThetaMax)
128  {
129  if(model) {
130  model->GenerateSecondaries(*pParticleChange, materials, steps,
131  normals, startingPosition, track);
132  }
133  Clear();
134  }
135 
136  return pParticleChange;
137 }
G4double GetStepLength() const
double x() const
G4TrackStatus GetTrackStatus() const
G4Navigator * GetNavigatorForTracking() const
G4Region * GetRegion() const
double z() const
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
std::vector< G4ThreeVector > normals
virtual G4ThreeVector GetLocalExitNormal(G4bool *valid)
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
G4Material * GetMaterial() const
static G4TransportationManager * GetTransportationManager()
const G4ThreeVector & GetMomentumDirection() const
G4LogicalVolume * GetLogicalVolume() const
double y() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VPhysicalVolume * GetVolume() const
double G4double
Definition: G4Types.hh:76
std::vector< const G4Material * > materials
std::vector< G4double > steps
const XML_Char XML_Content * model
Definition: expat.h:151

Here is the call graph for this function:

void G4VTransitionRadiation::PrintInfoDefinition ( )
virtual

Definition at line 164 of file G4VTransitionRadiation.cc.

165 {
166  if(model) model->PrintInfo();
167 }
const XML_Char XML_Content * model
Definition: expat.h:151
void G4VTransitionRadiation::SetModel ( G4VTRModel m)

Definition at line 157 of file G4VTransitionRadiation.cc.

158 {
159  model = mod;
160 }
const XML_Char XML_Content * model
Definition: expat.h:151
void G4VTransitionRadiation::SetRegion ( const G4Region reg)

Definition at line 150 of file G4VTransitionRadiation.cc.

151 {
152  region = reg;
153 }
static const G4double reg

Member Data Documentation

G4double G4VTransitionRadiation::cosDThetaMax

Definition at line 101 of file G4VTransitionRadiation.hh.

G4double G4VTransitionRadiation::gammaMin

Definition at line 100 of file G4VTransitionRadiation.hh.

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

Definition at line 89 of file G4VTransitionRadiation.hh.

G4VTRModel* G4VTransitionRadiation::model

Definition at line 96 of file G4VTransitionRadiation.hh.

std::vector<G4ThreeVector> G4VTransitionRadiation::normals

Definition at line 91 of file G4VTransitionRadiation.hh.

G4int G4VTransitionRadiation::nSteps

Definition at line 98 of file G4VTransitionRadiation.hh.

const G4Region* G4VTransitionRadiation::region

Definition at line 95 of file G4VTransitionRadiation.hh.

G4ThreeVector G4VTransitionRadiation::startingDirection

Definition at line 94 of file G4VTransitionRadiation.hh.

G4ThreeVector G4VTransitionRadiation::startingPosition

Definition at line 93 of file G4VTransitionRadiation.hh.

std::vector<G4double> G4VTransitionRadiation::steps

Definition at line 90 of file G4VTransitionRadiation.hh.


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