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 &)

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

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().

00049   : G4VDiscreteProcess(processName, type),
00050     region(0),
00051     model(0),
00052   nSteps(0),
00053   gammaMin(100),
00054   cosDThetaMax(std::cos(0.1))
00055 {
00056   Clear();
00057 }

G4VTransitionRadiation::~G4VTransitionRadiation (  )  [virtual]

Definition at line 61 of file G4VTransitionRadiation.cc.

References Clear().

00062 {
00063   Clear();
00064 }

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().

00069 {
00070   materials.clear();
00071   steps.clear();
00072   normals.clear();
00073   nSteps = 0;
00074 }

G4double G4VTransitionRadiation::GetMeanFreePath ( const G4Track track,
G4double  ,
G4ForceCondition condition 
) [inline, virtual]

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.

00106 {
00107   if(nSteps > 0) {
00108     *condition = StronglyForced;
00109   } else {
00110     *condition = NotForced;
00111     if(track.GetKineticEnergy()/track.GetDefinition()->GetPDGMass() + 1.0 > gammaMin &&
00112        track.GetVolume()->GetLogicalVolume()->GetRegion() == region) {
00113          *condition = StronglyForced;
00114     }
00115   }
00116   return DBL_MAX;      // so TR doesn't limit mean free path
00117 }

G4bool G4VTransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 139 of file G4VTransitionRadiation.cc.

References G4ParticleDefinition::GetPDGCharge().

00141 {
00142   return ( aParticle.GetPDGCharge() != 0.0 );
00143 }

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(), materials, model, CLHEP::detail::n, normals, nSteps, G4VProcess::pParticleChange, region, startingDirection, startingPosition, and steps.

00081 {
00082 
00083   // Fill temporary vectors
00084 
00085   const G4Material* material = track.GetMaterial();
00086   G4double length = step.GetStepLength();
00087   G4ThreeVector direction = track.GetMomentumDirection();
00088 
00089   if(nSteps == 0) {
00090 
00091     nSteps = 1;
00092     materials.push_back(material);
00093     steps.push_back(length);
00094     const G4StepPoint* point = step.GetPreStepPoint();
00095     startingPosition = point->GetPosition();
00096     startingDirection = point->GetMomentumDirection();
00097     G4bool valid = true;
00098     G4ThreeVector n = G4TransportationManager::GetTransportationManager()
00099                     ->GetNavigatorForTracking()->GetLocalExitNormal(&valid);
00100     if(valid) normals.push_back(n);
00101     else      normals.push_back(direction);
00102 
00103   } else {
00104 
00105     if(material == materials[nSteps-1]) {
00106       steps[nSteps-1] += length;
00107     } else {
00108       nSteps++;
00109       materials.push_back(material);
00110       steps.push_back(length);
00111       G4bool valid = true;
00112       G4ThreeVector n = G4TransportationManager::GetTransportationManager()
00113                       ->GetNavigatorForTracking()->GetLocalExitNormal(&valid);
00114       if(valid) normals.push_back(n);
00115       else      normals.push_back(direction);
00116     }
00117   }
00118 
00119   // Check POstStepPoint condition
00120 
00121   if(track.GetTrackStatus() == fStopAndKill ||
00122      track.GetVolume()->GetLogicalVolume()->GetRegion() != region ||
00123      startingDirection.x()*direction.x() +
00124      startingDirection.y()*direction.y() +
00125      startingDirection.z()*direction.z() < cosDThetaMax)
00126   {
00127      if(model) {
00128        model->GenerateSecondaries(*pParticleChange, materials, steps,
00129                                   normals, startingPosition, track);
00130      }
00131      Clear();
00132   }
00133 
00134   return pParticleChange;
00135 }

void G4VTransitionRadiation::PrintInfoDefinition (  )  [virtual]

Definition at line 162 of file G4VTransitionRadiation.cc.

References model, and G4VTRModel::PrintInfo().

00163 {
00164   if(model) model->PrintInfo();
00165 }

void G4VTransitionRadiation::SetModel ( G4VTRModel m  ) 

Definition at line 155 of file G4VTransitionRadiation.cc.

References model.

00156 {
00157   model = mod;
00158 }

void G4VTransitionRadiation::SetRegion ( const G4Region reg  ) 

Definition at line 148 of file G4VTransitionRadiation.cc.

References region.

00149 {
00150   region = reg;
00151 }


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.

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

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:
Generated on Mon May 27 17:53:56 2013 for Geant4 by  doxygen 1.4.7