Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Data Fields
G4VTransitionRadiation Class Reference

#include <G4VTransitionRadiation.hh>

Inheritance diagram for G4VTransitionRadiation:
G4VDiscreteProcess G4VProcess

Public Member Functions

 G4VTransitionRadiation (const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VTransitionRadiation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step)
 
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 &, 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 &)
 

Data Fields

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 
)

Definition at line 47 of file G4VTransitionRadiation.cc.

References Clear().

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 }
G4VTransitionRadiation::~G4VTransitionRadiation ( )
virtual

Definition at line 61 of file G4VTransitionRadiation.cc.

References Clear().

62 {
63  Clear();
64 }
G4VTransitionRadiation::G4VTransitionRadiation ( const G4VTransitionRadiation )

Member Function Documentation

void G4VTransitionRadiation::Clear ( )

Definition at line 68 of file G4VTransitionRadiation.cc.

References materials, normals, nSteps, and steps.

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

69 {
70  materials.clear();
71  steps.clear();
72  normals.clear();
73  nSteps = 0;
74 }
std::vector< G4ThreeVector > normals
std::vector< const G4Material * > materials
std::vector< G4double > steps
G4double G4VTransitionRadiation::GetMeanFreePath ( const G4Track track,
G4double  ,
G4ForceCondition condition 
)
inlinevirtual

Implements G4VDiscreteProcess.

Definition at line 103 of file G4VTransitionRadiation.hh.

References DBL_MAX, gammaMin, G4Track::GetDefinition(), G4Track::GetKineticEnergy(), G4VPhysicalVolume::GetLogicalVolume(), G4ParticleDefinition::GetPDGMass(), G4LogicalVolume::GetRegion(), G4Track::GetVolume(), NotForced, nSteps, region, and StronglyForced.

106 {
107  if(nSteps > 0) {
109  } else {
110  *condition = NotForced;
111  if(track.GetKineticEnergy()/track.GetDefinition()->GetPDGMass() + 1.0 > gammaMin &&
112  track.GetVolume()->GetLogicalVolume()->GetRegion() == region) {
114  }
115  }
116  return DBL_MAX; // so TR doesn't limit mean free path
117 }
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
G4bool G4VTransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 139 of file G4VTransitionRadiation.cc.

References G4ParticleDefinition::GetPDGCharge().

141 {
142  return ( aParticle.GetPDGCharge() != 0.0 );
143 }
G4VTransitionRadiation& G4VTransitionRadiation::operator= ( const G4VTransitionRadiation right)
G4VParticleChange * G4VTransitionRadiation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 78 of file G4VTransitionRadiation.cc.

References Clear(), cosDThetaMax, fStopAndKill, G4Navigator::GetLocalExitNormal(), G4VPhysicalVolume::GetLogicalVolume(), G4Track::GetMaterial(), G4StepPoint::GetMomentumDirection(), G4Track::GetMomentumDirection(), G4TransportationManager::GetNavigatorForTracking(), G4StepPoint::GetPosition(), G4Step::GetPreStepPoint(), G4LogicalVolume::GetRegion(), G4Step::GetStepLength(), G4Track::GetTrackStatus(), G4TransportationManager::GetTransportationManager(), G4Track::GetVolume(), eplot::material, materials, n, normals, nSteps, G4VProcess::pParticleChange, region, startingDirection, startingPosition, steps, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

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();
95  startingPosition = point->GetPosition();
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) {
128  model->GenerateSecondaries(*pParticleChange, materials, steps,
129  normals, startingPosition, track);
130  }
131  Clear();
132  }
133 
134  return pParticleChange;
135 }
G4double GetStepLength() const
double x() const
G4TrackStatus GetTrackStatus() const
G4Navigator * GetNavigatorForTracking() const
G4Region * GetRegion() const
double z() const
string material
Definition: eplot.py:19
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
const XML_Char XML_Content * model
const G4int n
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
void G4VTransitionRadiation::PrintInfoDefinition ( )
virtual

Definition at line 162 of file G4VTransitionRadiation.cc.

163 {
164  if(model) model->PrintInfo();
165 }
const XML_Char XML_Content * model
void G4VTransitionRadiation::SetModel ( G4VTRModel m)

Definition at line 155 of file G4VTransitionRadiation.cc.

156 {
157  model = mod;
158 }
const XML_Char XML_Content * model
void G4VTransitionRadiation::SetRegion ( const G4Region reg)

Definition at line 148 of file G4VTransitionRadiation.cc.

References region.

149 {
150  region = reg;
151 }

Field Documentation

G4double G4VTransitionRadiation::cosDThetaMax

Definition at line 99 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt().

G4double G4VTransitionRadiation::gammaMin

Definition at line 98 of file G4VTransitionRadiation.hh.

Referenced by GetMeanFreePath().

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

Definition at line 87 of file G4VTransitionRadiation.hh.

Referenced by Clear(), and PostStepDoIt().

G4VTRModel* G4VTransitionRadiation::model

Definition at line 94 of file G4VTransitionRadiation.hh.

std::vector<G4ThreeVector> G4VTransitionRadiation::normals

Definition at line 89 of file G4VTransitionRadiation.hh.

Referenced by Clear(), and PostStepDoIt().

G4int G4VTransitionRadiation::nSteps

Definition at line 96 of file G4VTransitionRadiation.hh.

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

const G4Region* G4VTransitionRadiation::region

Definition at line 93 of file G4VTransitionRadiation.hh.

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

G4ThreeVector G4VTransitionRadiation::startingDirection

Definition at line 92 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt().

G4ThreeVector G4VTransitionRadiation::startingPosition

Definition at line 91 of file G4VTransitionRadiation.hh.

Referenced by PostStepDoIt().

std::vector<G4double> G4VTransitionRadiation::steps

Definition at line 88 of file G4VTransitionRadiation.hh.

Referenced by Clear(), and PostStepDoIt().


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