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

#include <G4AdjointAlongStepWeightCorrection.hh>

Inheritance diagram for G4AdjointAlongStepWeightCorrection:
Collaboration diagram for G4AdjointAlongStepWeightCorrection:

Public Member Functions

 G4AdjointAlongStepWeightCorrection (const G4String &name="ContinuousWeightCorrection", G4ProcessType type=fElectromagnetic)
 
virtual ~G4AdjointAlongStepWeightCorrection ()
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- 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 G4bool IsApplicable (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 &)
 

Protected Member Functions

virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
- 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

G4ParticleChangefParticleChange
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 70 of file G4AdjointAlongStepWeightCorrection.hh.

Constructor & Destructor Documentation

G4AdjointAlongStepWeightCorrection::G4AdjointAlongStepWeightCorrection ( const G4String name = "ContinuousWeightCorrection",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 37 of file G4AdjointAlongStepWeightCorrection.cc.

38  : G4VContinuousProcess(name, type)
40  currentMaterialIndex=0;
41  preStepKinEnergy=1.;
42  currentCouple=0;
43 }
G4AdjointAlongStepWeightCorrection::~G4AdjointAlongStepWeightCorrection ( )
virtual

Definition at line 48 of file G4AdjointAlongStepWeightCorrection.cc.

Member Function Documentation

G4VParticleChange * G4AdjointAlongStepWeightCorrection::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4VContinuousProcess.

Definition at line 66 of file G4AdjointAlongStepWeightCorrection.cc.

68 {
69 
71 
72  // Get the actual (true) Step length
73  //----------------------------------
74  G4double length = step.GetStepLength();
75 
76 
78  G4ParticleDefinition* thePartDef= const_cast<G4ParticleDefinition*> (track.GetDynamicParticle()->GetDefinition());
80  preStepKinEnergy,Tkin, currentCouple,length);
81 
82 
83 
84 
85  //Caution!!!
86  // It is important to select the weight of the post_step_point
87  // as the current weight and not the weight of the track, as t
88  // the weight of the track is changed after having applied all
89  // the along_step_do_it.
90 
91  // G4double new_weight=weight_correction*track.GetWeight(); //old
92  G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight();
93 
94 
95 
96  //if (weight_correction >2.) new_weight=1.e-300;
97 
98 
99  //The following test check for zero weight.
100  //This happens after weight correction of gamma for photo electric effect.
101  //When the new weight is 0 it will be later on consider as nan by G4.
102  //Therefore we do put a lower limit of 1.e-300. for new_weight
103  //Correction by L.Desorgher on 15 July 2009
104  if (new_weight==0 || (new_weight<=0 && new_weight>0)){
105  //G4cout<<new_weight<<'\t'<<weight_correction<<'\t'<<track.GetWeight()<<G4endl;
106  new_weight=1.e-300;
107  }
108 
109  //G4cout<<new_weight<<'\t'<<weight_correction<<'\t'<<track.GetWeight()<<G4endl;
113 
114 
115  return fParticleChange;
116 
117 }
G4double GetWeight() const
G4double GetStepLength() const
const G4DynamicParticle * GetDynamicParticle() const
void ProposeParentWeight(G4double finalWeight)
G4ParticleDefinition * GetDefinition() const
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
virtual void Initialize(const G4Track &)
G4double GetContinuousWeightCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, G4double AfterStepEkin, const G4MaterialCutsCouple *aCouple, G4double step_length)
G4StepPoint * GetPostStepPoint() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
static G4AdjointCSManager * GetAdjointCSManager()

Here is the call graph for this function:

void G4AdjointAlongStepWeightCorrection::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 61 of file G4AdjointAlongStepWeightCorrection.cc.

62 {;
63 }
G4double G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
protectedvirtual

Implements G4VContinuousProcess.

Definition at line 120 of file G4AdjointAlongStepWeightCorrection.cc.

122 {
123  G4double x = DBL_MAX;
124  DefineMaterial(track.GetMaterialCutsCouple());
125  preStepKinEnergy = track.GetKineticEnergy();
126  return x;
127 }
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

void G4AdjointAlongStepWeightCorrection::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 53 of file G4AdjointAlongStepWeightCorrection.cc.

55 {
56 ;
57 }

Member Data Documentation

G4ParticleChange* G4AdjointAlongStepWeightCorrection::fParticleChange
protected

Definition at line 112 of file G4AdjointAlongStepWeightCorrection.hh.


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