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

#include <G4VAdjointReverseReaction.hh>

+ Inheritance diagram for G4VAdjointReverseReaction:

Public Member Functions

 G4VAdjointReverseReaction (G4String process_name, G4bool whichScatCase)
 
virtual ~G4VAdjointReverseReaction ()
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
void SetIntegralMode (G4bool aBool)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (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
 
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 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
 

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4VEmAdjointModeltheAdjointEMModel
 
G4ParticleChangefParticleChange
 
G4AdjointCSManagertheAdjointCSManager
 
G4bool IsScatProjToProjCase
 
- 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 74 of file G4VAdjointReverseReaction.hh.

Constructor & Destructor Documentation

◆ G4VAdjointReverseReaction()

G4VAdjointReverseReaction::G4VAdjointReverseReaction ( G4String  process_name,
G4bool  whichScatCase 
)

Definition at line 44 of file G4VAdjointReverseReaction.cc.

45 :
46 G4VDiscreteProcess(process_name)
48 IsScatProjToProjCase=whichScatCase;
50 IsFwdCSUsed=false;
51 IsIntegralModeUsed=false;
52 currentMaterialIndex=100000;
53 currentTcut=0.;
54 lastCS=0.;
55}
static G4AdjointCSManager * GetAdjointCSManager()
G4AdjointCSManager * theAdjointCSManager

◆ ~G4VAdjointReverseReaction()

G4VAdjointReverseReaction::~G4VAdjointReverseReaction ( )
virtual

Definition at line 58 of file G4VAdjointReverseReaction.cc.

61}

Member Function Documentation

◆ BuildPhysicsTable()

void G4VAdjointReverseReaction::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 69 of file G4VAdjointReverseReaction.cc.

70{
71
72 theAdjointCSManager->BuildCrossSectionMatrices(); //do not worry it will be done just once
74
75}

◆ GetMeanFreePath()

G4double G4VAdjointReverseReaction::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 109 of file G4VAdjointReverseReaction.cc.

113 G4double preStepKinEnergy = track.GetKineticEnergy();
114
115 /*G4double Sigma =
116 theAdjointEMModel->AdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase);*/
117
118 G4double Sigma =
120 G4double fwd_TotCS;
121 Sigma *= theAdjointCSManager->GetCrossSectionCorrection(track.GetDefinition(),preStepKinEnergy,track.GetMaterialCutsCouple(),IsFwdCSUsed, fwd_TotCS);
122 //G4cout<<fwd_TotCS<<G4endl;
123 /*if (IsFwdCSUsed && IsIntegralModeUsed){ //take the maximum cross section only for charged particle
124 G4double e_sigma_max, sigma_max;
125 theAdjointCSManager->GetMaxFwdTotalCS(track.GetDefinition(),
126 track.GetMaterialCutsCouple(), e_sigma_max, sigma_max);
127 if (e_sigma_max > preStepKinEnergy){
128 Sigma*=sigma_max/fwd_TotCS;
129 }
130 }
131 */
132
133 G4double mean_free_path = 1.e60 *mm;
134 if (Sigma>0) mean_free_path = 1./Sigma;
135 lastCS=Sigma;
136
137 /*G4cout<<"Sigma "<<Sigma<<G4endl;
138 G4cout<<"mean_free_path [mm] "<<mean_free_path/mm<<G4endl;
139 */
140
141
142 return mean_free_path;
143}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition: G4Types.hh:64
G4double GetCrossSectionCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, const G4MaterialCutsCouple *aCouple, G4bool &fwd_is_used, G4double &fwd_TotCS)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)

◆ PostStepDoIt()

G4VParticleChange * G4VAdjointReverseReaction::PostStepDoIt ( const G4Track track,
const G4Step  
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 78 of file G4VAdjointReverseReaction.cc.

79{
80
81
82
84
85 /* if (IsFwdCSUsed && IsIntegralModeUsed){ //INtegral mode still unstable
86 G4double Tkin = step.GetPostStepPoint()->GetKineticEnergy();
87 G4double fwdCS = theAdjointCSManager->GetTotalForwardCS(track.GetDefinition(), Tkin, track.GetMaterialCutsCouple());
88 //G4cout<<"lastCS "<<lastCS<<G4endl;
89 if (fwdCS<lastCS*G4UniformRand()) { // the reaction does not take place, same integral method as the one used for forward ionisation in G4
90 ClearNumberOfInteractionLengthLeft();
91 return fParticleChange;
92 }
93
94 }
95 */
96
100
102 return fParticleChange;
103
104
105
106}
virtual void Initialize(const G4Track &)
virtual void SampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418

◆ PreparePhysicsTable()

void G4VAdjointReverseReaction::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 64 of file G4VAdjointReverseReaction.cc.

65{;
66}

◆ SetIntegralMode()

Member Data Documentation

◆ fParticleChange

G4ParticleChange* G4VAdjointReverseReaction::fParticleChange
protected

◆ IsScatProjToProjCase

G4bool G4VAdjointReverseReaction::IsScatProjToProjCase
protected

◆ theAdjointCSManager

G4AdjointCSManager* G4VAdjointReverseReaction::theAdjointCSManager
protected

◆ theAdjointEMModel


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