Geant4_10
G4VAdjointReverseReaction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4VAdjointReverseReaction.cc 68044 2013-03-13 14:29:07Z gcosmo $
27 //
29 #include "G4SystemOfUnits.hh"
30 #include "G4AdjointCSManager.hh"
31 #include "G4AdjointCSMatrix.hh"
32 #include "G4AdjointInterpolator.hh"
33 #include "G4AdjointCSMatrix.hh"
34 #include "G4VEmAdjointModel.hh"
35 #include "G4ElementTable.hh"
36 #include "G4Element.hh"
37 #include "G4Material.hh"
38 #include "G4MaterialCutsCouple.hh"
39 #include "G4AdjointCSManager.hh"
40 #include "G4ParticleChange.hh"
41 #include "G4AdjointElectron.hh"
42 
43 
45  G4VAdjointReverseReaction(G4String process_name, G4bool whichScatCase):
46  G4VDiscreteProcess(process_name)
48  IsScatProjToProjCase=whichScatCase;
50  IsFwdCSUsed=false;
51  IsIntegralModeUsed=false;
52  lastCS=0.;
53 }
55 //
58 { if (fParticleChange) delete fParticleChange;
59 }
61 //
63 {;
64 }
66 //
68 {
69 
70  theAdjointCSManager->BuildCrossSectionMatrices(); //do not worry it will be done just once
72 
73 }
75 //
77 {
78 
79 
80 
82 
83  /* if (IsFwdCSUsed && IsIntegralModeUsed){ //INtegral mode still unstable
84  G4double Tkin = step.GetPostStepPoint()->GetKineticEnergy();
85  G4double fwdCS = theAdjointCSManager->GetTotalForwardCS(track.GetDefinition(), Tkin, track.GetMaterialCutsCouple());
86  //G4cout<<"lastCS "<<lastCS<<G4endl;
87  if (fwdCS<lastCS*G4UniformRand()) { // the reaction does not take place, same integral method as the one used for forward ionisation in G4
88  ClearNumberOfInteractionLengthLeft();
89  return fParticleChange;
90  }
91 
92  }
93  */
94 
98 
100  return fParticleChange;
101 
102 
103 
104 }
106 //
108  G4double ,
110 { *condition = NotForced;
111  G4double preStepKinEnergy = track.GetKineticEnergy();
112 
113  /*G4double Sigma =
114  theAdjointEMModel->AdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase);*/
115 
116  G4double Sigma =
118  G4double fwd_TotCS;
119  Sigma *= theAdjointCSManager->GetCrossSectionCorrection(track.GetDefinition(),preStepKinEnergy,track.GetMaterialCutsCouple(),IsFwdCSUsed, fwd_TotCS);
120  //G4cout<<fwd_TotCS<<G4endl;
121  /*if (IsFwdCSUsed && IsIntegralModeUsed){ //take the maximum cross section only for charged particle
122  G4double e_sigma_max, sigma_max;
123  theAdjointCSManager->GetMaxFwdTotalCS(track.GetDefinition(),
124  track.GetMaterialCutsCouple(), e_sigma_max, sigma_max);
125  if (e_sigma_max > preStepKinEnergy){
126  Sigma*=sigma_max/fwd_TotCS;
127  }
128  }
129  */
130 
131  G4double mean_free_path = 1.e60 *mm;
132  if (Sigma>0) mean_free_path = 1./Sigma;
133  lastCS=Sigma;
134 
135  /*G4cout<<"Sigma "<<Sigma<<G4endl;
136  G4cout<<"mean_free_path [mm] "<<mean_free_path/mm<<G4endl;
137  */
138 
139 
140  return mean_free_path;
141 }
G4double condition(const G4ErrorSymMatrix &m)
G4ParticleDefinition * GetDefinition() const
G4double GetCrossSectionCorrection(G4ParticleDefinition *aPartDef, G4double PreStepEkin, const G4MaterialCutsCouple *aCouple, G4bool &fwd_is_used, G4double &fwd_TotCS)
void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4VAdjointReverseReaction(G4String process_name, G4bool whichScatCase)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
virtual void SampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
G4double GetKineticEnergy() const
void PreparePhysicsTable(const G4ParticleDefinition &)
bool G4bool
Definition: G4Types.hh:79
G4AdjointCSManager * theAdjointCSManager
Definition: G4Step.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual void Initialize(const G4Track &)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
static G4AdjointCSManager * GetAdjointCSManager()