Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
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)
 
- Public Member Functions inherited from G4VContinuousProcess
 G4VContinuousProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousProcess (G4VContinuousProcess &)
 
virtual ~G4VContinuousProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangePostStepDoIt (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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VContinuousProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- 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 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.

34  : G4VContinuousProcess(processName, type)
35 {
36  if (verboseLevel>2) {
37  G4cout << GetProcessName() << " is created " << G4endl;
38  }
39 
40  theELossForExtrapolator = new G4EnergyLossForExtrapolator;
41 
42  theStepLimit = 1.;
43 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
G4ErrorEnergyLoss::~G4ErrorEnergyLoss ( )
virtual

Definition at line 51 of file G4ErrorEnergyLoss.cc.

52 {
53  delete theELossForExtrapolator;
54 }

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

59 {
61 
63 
64  G4double kinEnergyStart = aTrack.GetKineticEnergy();
65  G4double step_length = aStep.GetStepLength();
66 
67  const G4Material* aMaterial = aTrack.GetMaterial();
68  const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
69  G4double kinEnergyEnd = kinEnergyStart;
70 
71  if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
72  kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart,
73  step_length,
74  aMaterial,
75  aParticleDef );
76  G4double kinEnergyHalfStep = kinEnergyStart - (kinEnergyStart-kinEnergyEnd)/2.;
77 
78 #ifdef G4VERBOSE
80  G4cout << " G4ErrorEnergyLoss FWD end " << kinEnergyEnd
81  << " halfstep " << kinEnergyHalfStep << G4endl;
82 #endif
83 
84  //--- rescale to energy lost at 1/2 step
85  kinEnergyEnd = theELossForExtrapolator->EnergyBeforeStep( kinEnergyHalfStep,
86  step_length,
87  aMaterial,
88  aParticleDef );
89  kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
90  }else if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropForwards) ) {
91  kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyStart,
92  step_length,
93  aMaterial,
94  aParticleDef );
95  G4double kinEnergyHalfStep = kinEnergyStart - (kinEnergyStart-kinEnergyEnd)/2.;
96 #ifdef G4VERBOSE
98  G4cout << " G4ErrorEnergyLoss BCKD end " << kinEnergyEnd
99  << " halfstep " << kinEnergyHalfStep << G4endl;
100 #endif
101 
102  //--- rescale to energy lost at 1/2 step
103  kinEnergyEnd = theELossForExtrapolator->EnergyAfterStep( kinEnergyHalfStep,
104  step_length,
105  aMaterial,
106  aParticleDef );
107  kinEnergyEnd = kinEnergyStart - (kinEnergyHalfStep - kinEnergyEnd );
108  }
109 
110  G4double edepo = kinEnergyEnd - kinEnergyStart;
111 
112 #ifdef G4VERBOSE
113  if( G4ErrorPropagatorData::verbose() >= 2 )
114  G4cout << "AlongStepDoIt Estart= " << kinEnergyStart << " Eend " << kinEnergyEnd
115  << " Ediff " << kinEnergyStart-kinEnergyEnd << " step= " << step_length
116  << " mate= " << aMaterial->GetName()
117  << " particle= " << aParticleDef->GetParticleName() << G4endl;
118 #endif
119 
123 
124  aParticleChange.ProposeEnergy( kinEnergyEnd );
125 
126  return &aParticleChange;
127 }
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
const G4String & GetName() const
Definition: G4Material.hh:176
G4ParticleDefinition * GetDefinition() const
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ErrorMode GetMode() const
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4Material * GetMaterial() const
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
virtual void Initialize(const G4Track &)
void SetNumberOfSecondaries(G4int totSecondaries)
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4ErrorPropagatorData * GetErrorPropagatorData()
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().

135 {
137  if( theStepLimit != 1. ) {
138  G4double kinEnergyStart = aTrack.GetKineticEnergy();
139  G4double kinEnergyLoss = kinEnergyStart;
140  const G4Material* aMaterial = aTrack.GetMaterial();
141  const G4ParticleDefinition* aParticleDef = aTrack.GetDynamicParticle()->GetDefinition();
143  if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropBackwards) ) {
144  kinEnergyLoss = - kinEnergyStart +
145  theELossForExtrapolator->EnergyBeforeStep( kinEnergyStart, currentMinimumStep,
146  aMaterial, aParticleDef );
147  }else if( g4edata->GetMode() == G4ErrorMode(G4ErrorMode_PropForwards) ) {
148  kinEnergyLoss = kinEnergyStart -
149  theELossForExtrapolator->EnergyAfterStep( kinEnergyStart, currentMinimumStep,
150  aMaterial, aParticleDef );
151  }
152 #ifdef G4VERBOSE
153  if(G4ErrorPropagatorData::verbose() >= 3 )
154  G4cout << " G4ErrorEnergyLoss: currentMinimumStep " <<currentMinimumStep
155  << " kinEnergyLoss " << kinEnergyLoss
156  << " kinEnergyStart " << kinEnergyStart << G4endl;
157 #endif
158  if( kinEnergyLoss / kinEnergyStart > theStepLimit ) {
159  Step = theStepLimit / (kinEnergyLoss / kinEnergyStart) * currentMinimumStep;
160 #ifdef G4VERBOSE
161  if(G4ErrorPropagatorData::verbose() >= 2 )
162  G4cout << " G4ErrorEnergyLoss: limiting Step " << Step
163  << " energy loss fraction " << kinEnergyLoss / kinEnergyStart
164  << " > " << theStepLimit << G4endl;
165 #endif
166  }
167  }
168 
169  return Step;
170 
171 }
const G4DynamicParticle * GetDynamicParticle() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ErrorMode GetMode() const
G4double EnergyBeforeStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
G4Material * GetMaterial() const
G4double EnergyAfterStep(G4double kinEnergy, G4double step, const G4Material *, const G4ParticleDefinition *)
Definition: Step.hh:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4ErrorPropagatorData * GetErrorPropagatorData()
#define DBL_MAX
Definition: templates.hh:83
G4double G4ErrorEnergyLoss::GetStepLimit ( ) const
inline

Definition at line 72 of file G4ErrorEnergyLoss.hh.

72 { return theStepLimit; }
G4bool G4ErrorEnergyLoss::IsApplicable ( const G4ParticleDefinition aParticleType)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 95 of file G4ErrorEnergyLoss.hh.

References G4ParticleDefinition::GetPDGCharge().

96 {
97  return (aParticleType.GetPDGCharge() != 0);
98 }
G4double GetPDGCharge() const
void G4ErrorEnergyLoss::SetStepLimit ( G4double  val)
inline

Definition at line 73 of file G4ErrorEnergyLoss.hh.

Referenced by G4ErrorMessenger::SetNewValue().

73 { theStepLimit = val; }

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