G4ErrorEnergyLoss Class Reference

#include <G4ErrorEnergyLoss.hh>

Inheritance diagram for G4ErrorEnergyLoss:

G4VContinuousProcess G4VProcess

Public Member Functions

 G4ErrorEnergyLoss (const G4String &processName="G4ErrorEnergyLoss", G4ProcessType type=fElectromagnetic)
virtual ~G4ErrorEnergyLoss ()
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
G4double GetContinuousStepLimit (const G4Track &aTrack, G4double, G4double currentMinimumStep, G4double &)
G4VParticleChangeAlongStepDoIt (const G4Track &aTrack, const G4Step &aStep)
G4double GetStepLimit () const
void SetStepLimit (G4double val)

Detailed Description

Definition at line 47 of file G4ErrorEnergyLoss.hh.


Constructor & Destructor Documentation

G4ErrorEnergyLoss::G4ErrorEnergyLoss ( const G4String processName = "G4ErrorEnergyLoss",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 32 of file G4ErrorEnergyLoss.cc.

References G4cout, G4endl, G4VProcess::GetProcessName(), and G4VProcess::verboseLevel.

00034            : G4VContinuousProcess(processName, type)
00035 {
00036   if (verboseLevel>2) {
00037     G4cout << GetProcessName() << " is created " << G4endl;
00038   }
00039 
00040   theELossForExtrapolator = new G4EnergyLossForExtrapolator;
00041 
00042   theStepLimit = 1.;
00043 }

G4ErrorEnergyLoss::~G4ErrorEnergyLoss (  )  [virtual]

Definition at line 51 of file G4ErrorEnergyLoss.cc.

00052 {
00053   delete theELossForExtrapolator;
00054 }


Member Function Documentation

G4VParticleChange * G4ErrorEnergyLoss::AlongStepDoIt ( const G4Track aTrack,
const G4Step aStep 
) [virtual]

Reimplemented from G4VContinuousProcess.

Definition at line 58 of file G4ErrorEnergyLoss.cc.

References G4VProcess::aParticleChange, G4VParticleChange::ClearDebugFlag(), G4EnergyLossForExtrapolator::EnergyAfterStep(), G4EnergyLossForExtrapolator::EnergyBeforeStep(), G4cout, G4endl, G4ErrorMode_PropBackwards, G4ErrorMode_PropForwards, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4ErrorPropagatorData::GetErrorPropagatorData(), G4Track::GetKineticEnergy(), G4Track::GetMaterial(), G4ErrorPropagatorData::GetMode(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), G4Step::GetStepLength(), G4ParticleChange::Initialize(), G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::SetNumberOfSecondaries(), and G4ErrorPropagatorData::verbose().

00059 {
00060   aParticleChange.Initialize(aTrack);
00061 
00062   G4ErrorPropagatorData* g4edata =  G4ErrorPropagatorData::GetErrorPropagatorData();
00063 
00064   G4double kinEnergyStart = aTrack.GetKineticEnergy();
00065   G4double step_length  = aStep.GetStepLength();
00066 
00067   const G4Material* aMaterial = aTrack.GetMaterial();
00068   const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
00069   G4double kinEnergyEnd = kinEnergyStart;
00070 
00071   if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
00072     kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart, 
00073                                                               step_length, 
00074                                                               aMaterial, 
00075                                                               aParticleDef );
00076     G4double kinEnergyHalfStep = kinEnergyStart - (kinEnergyStart-kinEnergyEnd)/2.;
00077 
00078 #ifdef G4VERBOSE
00079   if(G4ErrorPropagatorData::verbose() >= 3 ) 
00080     G4cout << " G4ErrorEnergyLoss FWD  end " << kinEnergyEnd 
00081            << " halfstep " << kinEnergyHalfStep << G4endl;
00082 #endif
00083 
00084     //--- rescale to energy lost at 1/2 step
00085     kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyHalfStep, 
00086                                                               step_length, 
00087                                                               aMaterial, 
00088                                                               aParticleDef );
00089     kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
00090   }else if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropForwards) ) {
00091     kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyStart, 
00092                                                              step_length, 
00093                                                              aMaterial, 
00094                                                              aParticleDef );
00095     G4double kinEnergyHalfStep = kinEnergyStart - (kinEnergyStart-kinEnergyEnd)/2.;
00096 #ifdef G4VERBOSE
00097   if(G4ErrorPropagatorData::verbose() >= 3 ) 
00098     G4cout << " G4ErrorEnergyLoss BCKD  end " << kinEnergyEnd 
00099            << " halfstep " << kinEnergyHalfStep << G4endl;
00100 #endif
00101 
00102     //--- rescale to energy lost at 1/2 step
00103     kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyHalfStep, 
00104                                                              step_length, 
00105                                                              aMaterial, 
00106                                                              aParticleDef );
00107     kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
00108   }
00109 
00110   G4double edepo = kinEnergyEnd - kinEnergyStart;
00111 
00112 #ifdef G4VERBOSE
00113   if( G4ErrorPropagatorData::verbose() >= 2 ) 
00114     G4cout << "AlongStepDoIt Estart= " << kinEnergyStart << " Eend " << kinEnergyEnd 
00115            << " Ediff " << kinEnergyStart-kinEnergyEnd << " step= " << step_length 
00116            << " mate= " << aMaterial->GetName() 
00117            << " particle= " << aParticleDef->GetParticleName() << G4endl;
00118 #endif
00119 
00120   aParticleChange.ClearDebugFlag();
00121   aParticleChange.ProposeLocalEnergyDeposit( edepo );
00122   aParticleChange.SetNumberOfSecondaries(0);
00123  
00124   aParticleChange.ProposeEnergy( kinEnergyEnd );
00125   
00126   return &aParticleChange;
00127 }

G4double G4ErrorEnergyLoss::GetContinuousStepLimit ( const G4Track aTrack,
G4double  ,
G4double  currentMinimumStep,
G4double  
) [virtual]

Implements G4VContinuousProcess.

Definition at line 131 of file G4ErrorEnergyLoss.cc.

References DBL_MAX, G4EnergyLossForExtrapolator::EnergyAfterStep(), G4EnergyLossForExtrapolator::EnergyBeforeStep(), G4cout, G4endl, G4ErrorMode_PropBackwards, G4ErrorMode_PropForwards, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4ErrorPropagatorData::GetErrorPropagatorData(), G4Track::GetKineticEnergy(), G4Track::GetMaterial(), G4ErrorPropagatorData::GetMode(), and G4ErrorPropagatorData::verbose().

00135 {
00136   G4double Step = DBL_MAX;
00137   if( theStepLimit != 1. ) { 
00138     G4double kinEnergyStart = aTrack.GetKineticEnergy();
00139     G4double kinEnergyLoss = kinEnergyStart;
00140     const G4Material* aMaterial = aTrack.GetMaterial();
00141     const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
00142     G4ErrorPropagatorData* g4edata =  G4ErrorPropagatorData::GetErrorPropagatorData();
00143     if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
00144       kinEnergyLoss = - kinEnergyStart + 
00145         theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart, currentMinimumStep, 
00146                                                    aMaterial, aParticleDef );
00147     }else if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropForwards) ) {
00148       kinEnergyLoss = kinEnergyStart - 
00149         theELossForExtrapolator->EnergyAfterStep( kinEnergyStart, currentMinimumStep, 
00150                                                   aMaterial, aParticleDef );
00151     }
00152 #ifdef G4VERBOSE
00153   if(G4ErrorPropagatorData::verbose() >= 3 ) 
00154     G4cout << " G4ErrorEnergyLoss: currentMinimumStep " <<currentMinimumStep 
00155            << "  kinEnergyLoss " << kinEnergyLoss 
00156            << " kinEnergyStart " << kinEnergyStart << G4endl;
00157 #endif
00158     if( kinEnergyLoss / kinEnergyStart > theStepLimit ) {
00159       Step = theStepLimit / (kinEnergyLoss / kinEnergyStart)  * currentMinimumStep;
00160 #ifdef G4VERBOSE
00161   if(G4ErrorPropagatorData::verbose() >= 2 ) 
00162     G4cout << " G4ErrorEnergyLoss: limiting Step " << Step 
00163            << " energy loss fraction " << kinEnergyLoss / kinEnergyStart 
00164            << " > " << theStepLimit << G4endl;
00165 #endif
00166     }
00167   }
00168   
00169   return Step;
00170 
00171 }

G4double G4ErrorEnergyLoss::GetStepLimit (  )  const [inline]

Definition at line 72 of file G4ErrorEnergyLoss.hh.

00072 { return theStepLimit; }

G4bool G4ErrorEnergyLoss::IsApplicable ( const G4ParticleDefinition aParticleType  )  [inline, virtual]

Reimplemented from G4VProcess.

Definition at line 95 of file G4ErrorEnergyLoss.hh.

References G4ParticleDefinition::GetPDGCharge().

00096 {
00097    return (aParticleType.GetPDGCharge() != 0);
00098 }

void G4ErrorEnergyLoss::SetStepLimit ( G4double  val  )  [inline]

Definition at line 73 of file G4ErrorEnergyLoss.hh.

Referenced by G4ErrorMessenger::SetNewValue().

00073 { theStepLimit = val; }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:54 2013 for Geant4 by  doxygen 1.4.7