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

#include <G4VAdjointReverseReaction.hh>

Inheritance diagram for G4VAdjointReverseReaction:
Collaboration 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 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
 
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 GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
- 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 ( G4String  process_name,
G4bool  whichScatCase 
)

Definition at line 45 of file G4VAdjointReverseReaction.cc.

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

Here is the call graph for this function:

G4VAdjointReverseReaction::~G4VAdjointReverseReaction ( )
virtual

Definition at line 58 of file G4VAdjointReverseReaction.cc.

59 { if (fParticleChange) delete fParticleChange;
60 }

Member Function Documentation

void G4VAdjointReverseReaction::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 68 of file G4VAdjointReverseReaction.cc.

69 {
70 
71  theAdjointCSManager->BuildCrossSectionMatrices(); //do not worry it will be done just once
73 
74 }
G4AdjointCSManager * theAdjointCSManager

Here is the call graph for this function:

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

Implements G4VDiscreteProcess.

Definition at line 106 of file G4VAdjointReverseReaction.cc.

109 { *condition = NotForced;
110  G4double preStepKinEnergy = track.GetKineticEnergy();
111 
112  if(track.GetTrackID() != trackid) {
113  trackid = track.GetTrackID();
114  nstep = 0;
115  }
116  ++nstep;
117 
118 
119 
120  /*G4double Sigma =
121  theAdjointEMModel->AdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase);*/
122 
123  G4double Sigma =
125 
126  //G4double sig = Sigma;
127 
128  G4double fwd_TotCS;
129  G4double corr = theAdjointCSManager->GetCrossSectionCorrection(track.GetDefinition(),preStepKinEnergy,track.GetMaterialCutsCouple(),IsFwdCSUsed, fwd_TotCS);
130 
131  if(std::fabs(corr) > 100.) { Sigma = 0.0; }
132  else { Sigma *= corr; }
133 
134  //G4cout<<fwd_TotCS<<G4endl;
135  /*if (IsFwdCSUsed && IsIntegralModeUsed){ //take the maximum cross section only for charged particle
136  G4double e_sigma_max, sigma_max;
137  theAdjointCSManager->GetMaxFwdTotalCS(track.GetDefinition(),
138  track.GetMaterialCutsCouple(), e_sigma_max, sigma_max);
139  if (e_sigma_max > preStepKinEnergy){
140  Sigma*=sigma_max/fwd_TotCS;
141  }
142  }
143  */
144 
145  G4double mean_free_path = 1.e60 *mm;
146  if (Sigma>0) mean_free_path = 1./Sigma;
147  lastCS=Sigma;
148  /*
149  if(nstep > 100) {
150 
151  G4cout << "#* " << track.GetDefinition()->GetParticleName()
152  << " " << GetProcessName()
153  << " Nstep " << nstep
154  << " E(MeV)= " << preStepKinEnergy << " Sig0= " << sig
155  << " sig1= " << Sigma << " mfp= " << mean_free_path << G4endl;
156 
157  }
158  if (nstep > 20000) {
159  exit(1);
160  }
161  */
162  /*G4cout<<"Sigma "<<Sigma<<G4endl;
163  G4cout<<"mean_free_path [mm] "<<mean_free_path/mm<<G4endl;
164  */
165 
166 
167  return mean_free_path;
168 }
G4double condition(const G4ErrorSymMatrix &m)
G4ParticleDefinition * GetDefinition() const
G4double GetCrossSectionCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, const G4MaterialCutsCouple *aCouple, G4bool &fwd_is_used, G4double &fwd_TotCS)
static constexpr double mm
Definition: G4SIunits.hh:115
virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetKineticEnergy() const
G4AdjointCSManager * theAdjointCSManager
G4int GetTrackID() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 77 of file G4VAdjointReverseReaction.cc.

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

Here is the call graph for this function:

void G4VAdjointReverseReaction::PreparePhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 63 of file G4VAdjointReverseReaction.cc.

64 {;
65 }
void G4VAdjointReverseReaction::SetIntegralMode ( G4bool  aBool)
inline

Definition at line 89 of file G4VAdjointReverseReaction.hh.

89 {IsIntegralModeUsed = aBool;}

Here is the caller graph for this function:

Member Data Documentation

G4ParticleChange* G4VAdjointReverseReaction::fParticleChange
protected

Definition at line 99 of file G4VAdjointReverseReaction.hh.

G4bool G4VAdjointReverseReaction::IsScatProjToProjCase
protected

Definition at line 101 of file G4VAdjointReverseReaction.hh.

G4AdjointCSManager* G4VAdjointReverseReaction::theAdjointCSManager
protected

Definition at line 100 of file G4VAdjointReverseReaction.hh.

G4VEmAdjointModel* G4VAdjointReverseReaction::theAdjointEMModel
protected

Definition at line 98 of file G4VAdjointReverseReaction.hh.


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