Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AdjointAlongStepWeightCorrection Class Reference

#include <G4AdjointAlongStepWeightCorrection.hh>

+ Inheritance diagram for G4AdjointAlongStepWeightCorrection:

Public Member Functions

 G4AdjointAlongStepWeightCorrection (const G4String &name="ContinuousWeightCorrection", G4ProcessType type=fElectromagnetic)
 
 ~G4AdjointAlongStepWeightCorrection () override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
void ProcessDescription (std::ostream &) const override
 
void DumpInfo () const override
 
 G4AdjointAlongStepWeightCorrection (G4AdjointAlongStepWeightCorrection &)=delete
 
G4AdjointAlongStepWeightCorrectionoperator= (const G4AdjointAlongStepWeightCorrection &right)=delete
 
- Public Member Functions inherited from G4VContinuousProcess
 G4VContinuousProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VContinuousProcess (G4VContinuousProcess &)
 
virtual ~G4VContinuousProcess ()
 
G4VContinuousProcessoperator= (const G4VContinuousProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 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 const G4VProcessGetCreatorProcess () const
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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

G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
 
- Protected Member Functions inherited from G4VContinuousProcess
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 59 of file G4AdjointAlongStepWeightCorrection.hh.

Constructor & Destructor Documentation

◆ G4AdjointAlongStepWeightCorrection() [1/2]

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

Definition at line 36 of file G4AdjointAlongStepWeightCorrection.cc.

38 : G4VContinuousProcess(name, type)
39{
40 fParticleChange = new G4ParticleChange();
42}
static G4AdjointCSManager * GetAdjointCSManager()

◆ ~G4AdjointAlongStepWeightCorrection()

G4AdjointAlongStepWeightCorrection::~G4AdjointAlongStepWeightCorrection ( )
override

Definition at line 45 of file G4AdjointAlongStepWeightCorrection.cc.

46{
47 delete fParticleChange;
48}

◆ G4AdjointAlongStepWeightCorrection() [2/2]

G4AdjointAlongStepWeightCorrection::G4AdjointAlongStepWeightCorrection ( G4AdjointAlongStepWeightCorrection )
delete

Member Function Documentation

◆ AlongStepDoIt()

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

Reimplemented from G4VContinuousProcess.

Definition at line 68 of file G4AdjointAlongStepWeightCorrection.cc.

70{
71 fParticleChange->Initialize(track);
72
73 // Get the actual (true) Step length
74 G4double length = step.GetStepLength();
75
77 G4ParticleDefinition* thePartDef = const_cast<G4ParticleDefinition*>(
79 G4double weight_correction = fCSManager->GetContinuousWeightCorrection(
80 thePartDef, fPreStepKinEnergy, Tkin, fCurrentCouple, length);
81
82 // Caution!!!
83 // It is important to select the weight of the post_step_point as the current
84 // weight and not the weight of the track, as the weight of the track is
85 // changed after having applied all the along_step_do_it.
86 G4double new_weight =
87 weight_correction * step.GetPostStepPoint()->GetWeight();
88
89 // The following test check for zero weight.
90 // This happens after weight correction of gamma for photo electric effect.
91 // When the new weight is 0 it will be later on considered as NaN by G4.
92 // Therefore we put a lower limit of 1.e-300. for new_weight
93 if(new_weight == 0. || (new_weight <= 0. && new_weight > 0.))
94 {
95 new_weight = 1.e-300;
96 }
97
98 fParticleChange->SetParentWeightByProcess(false);
99 fParticleChange->SetSecondaryWeightByProcess(false);
100 fParticleChange->ProposeParentWeight(new_weight);
101
102 return fParticleChange;
103}
double G4double
Definition: G4Types.hh:83
G4double GetContinuousWeightCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, G4double AfterStepEkin, const G4MaterialCutsCouple *aCouple, G4double step_length)
G4ParticleDefinition * GetDefinition() const
void Initialize(const G4Track &) override
G4double GetKineticEnergy() const
G4double GetWeight() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)

◆ DumpInfo()

void G4AdjointAlongStepWeightCorrection::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 71 of file G4AdjointAlongStepWeightCorrection.hh.

G4GLOB_DLL std::ostream G4cout
void ProcessDescription(std::ostream &) const override

◆ GetContinuousStepLimit()

G4double G4AdjointAlongStepWeightCorrection::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
)
overrideprotectedvirtual

Implements G4VContinuousProcess.

Definition at line 106 of file G4AdjointAlongStepWeightCorrection.cc.

108{
109 DefineMaterial(track.GetMaterialCutsCouple());
110 fPreStepKinEnergy = track.GetKineticEnergy();
111 return DBL_MAX;
112}
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
#define DBL_MAX
Definition: templates.hh:62

◆ operator=()

G4AdjointAlongStepWeightCorrection & G4AdjointAlongStepWeightCorrection::operator= ( const G4AdjointAlongStepWeightCorrection right)
delete

◆ ProcessDescription()

void G4AdjointAlongStepWeightCorrection::ProcessDescription ( std::ostream &  out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 51 of file G4AdjointAlongStepWeightCorrection.cc.

53{
54 out <<
55 "Continuous processes act on adjoint particles to continuously correct their "
56 "weight during the adjoint reverse tracking. This process is needed when "
57 "the adjoint cross sections are not scaled such that the total adjoint cross "
58 "section matches the total forward cross section. By default the mode where "
59 "the total adjoint cross section is equal to the total forward cross section "
60 "is used and therefore this along step weightcorrection factor is 1. However "
61 "in some cases (some energy ranges) the total forward cross section or the "
62 "total adjoint cross section can be zero. In this case the along step weight "
63 "correction is needed and is given by exp(-(Sigma_tot_adj-Sigma_tot_fwd).dx)"
64 "\n";
65}

Referenced by DumpInfo().


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