G4VAdjointReverseReaction Class Reference

#include <G4VAdjointReverseReaction.hh>

Inheritance diagram for G4VAdjointReverseReaction:

G4VDiscreteProcess G4VProcess G4eInverseBremsstrahlung G4eInverseCompton G4eInverseIonisation G4hInverseIonisation G4InversePEEffect G4IonInverseIonisation

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)

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)

Protected Attributes

G4VEmAdjointModeltheAdjointEMModel
G4ParticleChangefParticleChange
G4AdjointCSManagertheAdjointCSManager
G4bool IsScatProjToProjCase

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.

References fParticleChange, G4AdjointCSManager::GetAdjointCSManager(), IsScatProjToProjCase, and theAdjointCSManager.

00045                                                                               :
00046                         G4VDiscreteProcess(process_name)
00047 {theAdjointCSManager = G4AdjointCSManager::GetAdjointCSManager();
00048  IsScatProjToProjCase=whichScatCase;
00049  fParticleChange=new G4ParticleChange();
00050  IsFwdCSUsed=false;
00051  IsIntegralModeUsed=false;
00052  lastCS=0.;
00053 }

G4VAdjointReverseReaction::~G4VAdjointReverseReaction (  )  [virtual]

Definition at line 57 of file G4VAdjointReverseReaction.cc.

References fParticleChange.

00058 { if (fParticleChange) delete fParticleChange;
00059 }                       


Member Function Documentation

void G4VAdjointReverseReaction::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 67 of file G4VAdjointReverseReaction.cc.

References G4AdjointCSManager::BuildCrossSectionMatrices(), G4AdjointCSManager::BuildTotalSigmaTables(), and theAdjointCSManager.

00068 {
00069 
00070  theAdjointCSManager->BuildCrossSectionMatrices(); //do not worry it will be done just once
00071  theAdjointCSManager->BuildTotalSigmaTables();
00072 
00073 }

G4double G4VAdjointReverseReaction::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, virtual]

Implements G4VDiscreteProcess.

Definition at line 107 of file G4VAdjointReverseReaction.cc.

References G4VEmAdjointModel::GetAdjointCrossSection(), G4AdjointCSManager::GetCrossSectionCorrection(), G4Track::GetDefinition(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), IsScatProjToProjCase, NotForced, theAdjointCSManager, and theAdjointEMModel.

00110 { *condition = NotForced;
00111   G4double preStepKinEnergy = track.GetKineticEnergy();
00112   
00113   /*G4double Sigma =
00114                 theAdjointEMModel->AdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase);*/
00115                 
00116   G4double Sigma =
00117                 theAdjointEMModel->GetAdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase); 
00118   G4double fwd_TotCS;
00119   Sigma *=  theAdjointCSManager->GetCrossSectionCorrection(track.GetDefinition(),preStepKinEnergy,track.GetMaterialCutsCouple(),IsFwdCSUsed, fwd_TotCS);
00120   //G4cout<<fwd_TotCS<<G4endl;
00121   /*if (IsFwdCSUsed && IsIntegralModeUsed){ //take the maximum cross section only for charged particle          
00122         G4double e_sigma_max, sigma_max;
00123         theAdjointCSManager->GetMaxFwdTotalCS(track.GetDefinition(),
00124                                      track.GetMaterialCutsCouple(), e_sigma_max, sigma_max);
00125         if (e_sigma_max > preStepKinEnergy){
00126                 Sigma*=sigma_max/fwd_TotCS;
00127         }                    
00128   }
00129   */            
00130 
00131   G4double mean_free_path = 1.e60 *mm; 
00132   if (Sigma>0) mean_free_path = 1./Sigma;
00133   lastCS=Sigma;
00134   
00135   /*G4cout<<"Sigma  "<<Sigma<<G4endl;
00136   G4cout<<"mean_free_path [mm] "<<mean_free_path/mm<<G4endl;
00137   */
00138   
00139 
00140   return mean_free_path;
00141 }                                        

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

Reimplemented from G4VDiscreteProcess.

Definition at line 76 of file G4VAdjointReverseReaction.cc.

References G4VProcess::ClearNumberOfInteractionLengthLeft(), fParticleChange, G4ParticleChange::Initialize(), IsScatProjToProjCase, G4VEmAdjointModel::SampleSecondaries(), and theAdjointEMModel.

00077 { 
00078   
00079   
00080   
00081   fParticleChange->Initialize(track);
00082  
00083  /* if (IsFwdCSUsed && IsIntegralModeUsed){ //INtegral mode still unstable
00084          G4double Tkin = step.GetPostStepPoint()->GetKineticEnergy();
00085          G4double fwdCS = theAdjointCSManager->GetTotalForwardCS(track.GetDefinition(), Tkin, track.GetMaterialCutsCouple());
00086          //G4cout<<"lastCS "<<lastCS<<G4endl;
00087          if (fwdCS<lastCS*G4UniformRand()) { // the reaction does not take place, same integral method as the one used for forward ionisation in  G4 
00088                 ClearNumberOfInteractionLengthLeft();
00089                 return fParticleChange;
00090          } 
00091          
00092   }
00093  */
00094  
00095   theAdjointEMModel->SampleSecondaries(track,
00096                                        IsScatProjToProjCase,
00097                                         fParticleChange);
00098   
00099   ClearNumberOfInteractionLengthLeft();
00100   return fParticleChange;
00101                         
00102    
00103   
00104 }

void G4VAdjointReverseReaction::PreparePhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 62 of file G4VAdjointReverseReaction.cc.

00063 {;
00064 }

void G4VAdjointReverseReaction::SetIntegralMode ( G4bool  aBool  )  [inline]

Definition at line 89 of file G4VAdjointReverseReaction.hh.

Referenced by G4eInverseBremsstrahlung::G4eInverseBremsstrahlung(), G4eInverseCompton::G4eInverseCompton(), G4eInverseIonisation::G4eInverseIonisation(), G4hInverseIonisation::G4hInverseIonisation(), G4InversePEEffect::G4InversePEEffect(), and G4IonInverseIonisation::G4IonInverseIonisation().

00089 {IsIntegralModeUsed = aBool;}


Field Documentation

G4ParticleChange* G4VAdjointReverseReaction::fParticleChange [protected]

Definition at line 99 of file G4VAdjointReverseReaction.hh.

Referenced by G4VAdjointReverseReaction(), PostStepDoIt(), and ~G4VAdjointReverseReaction().

G4bool G4VAdjointReverseReaction::IsScatProjToProjCase [protected]

Definition at line 101 of file G4VAdjointReverseReaction.hh.

Referenced by G4eInverseBremsstrahlung::G4eInverseBremsstrahlung(), G4VAdjointReverseReaction(), GetMeanFreePath(), and PostStepDoIt().

G4AdjointCSManager* G4VAdjointReverseReaction::theAdjointCSManager [protected]

Definition at line 100 of file G4VAdjointReverseReaction.hh.

Referenced by BuildPhysicsTable(), G4VAdjointReverseReaction(), and GetMeanFreePath().

G4VEmAdjointModel* G4VAdjointReverseReaction::theAdjointEMModel [protected]

Definition at line 98 of file G4VAdjointReverseReaction.hh.

Referenced by G4eInverseBremsstrahlung::G4eInverseBremsstrahlung(), G4eInverseCompton::G4eInverseCompton(), G4eInverseIonisation::G4eInverseIonisation(), G4hInverseIonisation::G4hInverseIonisation(), G4InversePEEffect::G4InversePEEffect(), G4IonInverseIonisation::G4IonInverseIonisation(), GetMeanFreePath(), and PostStepDoIt().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:39 2013 for Geant4 by  doxygen 1.4.7