source: trunk/source/processes/electromagnetic/adjoint/src/G4VAdjointReverseReaction.cc@ 1331

Last change on this file since 1331 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 5.3 KB
Line 
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,v 1.2 2009/11/20 10:31:20 ldesorgh Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29#include "G4VAdjointReverseReaction.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
44G4VAdjointReverseReaction::
45 G4VAdjointReverseReaction(G4String process_name, G4bool whichScatCase):
46 G4VDiscreteProcess(process_name)
47{theAdjointCSManager = G4AdjointCSManager::GetAdjointCSManager();
48 IsScatProjToProjCase=whichScatCase;
49 fParticleChange=new G4ParticleChange();
50
51}
52//////////////////////////////////////////////////////////////////////////////
53//
54G4VAdjointReverseReaction::
55 ~G4VAdjointReverseReaction()
56{;
57}
58//////////////////////////////////////////////////////////////////////////////
59//
60void G4VAdjointReverseReaction::PreparePhysicsTable(const G4ParticleDefinition&)
61{;
62}
63//////////////////////////////////////////////////////////////////////////////
64//
65void G4VAdjointReverseReaction::BuildPhysicsTable(const G4ParticleDefinition&)
66{
67
68 theAdjointCSManager->BuildCrossSectionMatrices(); //do not worry it will be done just once
69 theAdjointCSManager->BuildTotalSigmaTables();
70
71}
72//////////////////////////////////////////////////////////////////////////////
73//
74G4VParticleChange* G4VAdjointReverseReaction::PostStepDoIt(const G4Track& track, const G4Step& )
75{
76
77
78
79 fParticleChange->Initialize(track);
80
81 /* if (IsFwdCSUsed && IsIntegralModeUsed){ //INtegral mode still unstable
82 G4double Tkin = step.GetPostStepPoint()->GetKineticEnergy();
83 G4double fwdCS = theAdjointCSManager->GetTotalForwardCS(track.GetDefinition(), Tkin, track.GetMaterialCutsCouple());
84 //G4cout<<"lastCS "<<lastCS<<G4endl;
85 if (fwdCS<lastCS*G4UniformRand()) { // the reaction does not take place, same integral method as the one used for forward ionisation in G4
86 ClearNumberOfInteractionLengthLeft();
87 return fParticleChange;
88 }
89
90 }
91 */
92
93 theAdjointEMModel->SampleSecondaries(track,
94 IsScatProjToProjCase,
95 fParticleChange);
96
97 ClearNumberOfInteractionLengthLeft();
98 return fParticleChange;
99
100
101
102}
103//////////////////////////////////////////////////////////////////////////////
104//
105G4double G4VAdjointReverseReaction::GetMeanFreePath(const G4Track& track,
106 G4double ,
107 G4ForceCondition* condition)
108{ *condition = NotForced;
109 G4double preStepKinEnergy = track.GetKineticEnergy();
110
111 G4double Sigma =
112 theAdjointEMModel->AdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase);
113
114 G4double fwd_TotCS;
115 Sigma *= theAdjointCSManager->GetCrossSectionCorrection(track.GetDefinition(),preStepKinEnergy,track.GetMaterialCutsCouple(),IsFwdCSUsed, fwd_TotCS);
116 //G4cout<<fwd_TotCS<<G4endl;
117 /*if (IsFwdCSUsed && IsIntegralModeUsed){ //take the maximum cross section only for charged particle
118 G4double e_sigma_max, sigma_max;
119 theAdjointCSManager->GetMaxFwdTotalCS(track.GetDefinition(),
120 track.GetMaterialCutsCouple(), e_sigma_max, sigma_max);
121 if (e_sigma_max > preStepKinEnergy){
122 Sigma*=sigma_max/fwd_TotCS;
123 }
124 }
125 */
126
127 G4double mean_free_path = 1.e60 *mm;
128 if (Sigma>0) mean_free_path = 1./Sigma;
129 lastCS=Sigma;
130
131 /*G4cout<<"Sigma "<<Sigma<<G4endl;
132 G4cout<<"mean_free_path [mm] "<<mean_free_path/mm<<G4endl;
133 */
134
135
136 return mean_free_path;
137}
138
Note: See TracBrowser for help on using the repository browser.