Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
adjoint
src
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 69844 2013-05-16 09:19:33Z gcosmo $
27
//
28
#include "
G4VAdjointReverseReaction.hh
"
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
44
G4VAdjointReverseReaction::
45
G4VAdjointReverseReaction
(
G4String
process_name,
G4bool
whichScatCase):
46
G4VDiscreteProcess
(process_name)
47
{
theAdjointCSManager
=
G4AdjointCSManager::GetAdjointCSManager
();
48
IsScatProjToProjCase
=whichScatCase;
49
fParticleChange
=
new
G4ParticleChange
();
50
IsFwdCSUsed=
false
;
51
IsIntegralModeUsed=
false
;
52
lastCS=0.;
53
}
55
//
56
G4VAdjointReverseReaction::
57
~G4VAdjointReverseReaction
()
58
{
if
(
fParticleChange
)
delete
fParticleChange
;
59
}
61
//
62
void
G4VAdjointReverseReaction::PreparePhysicsTable
(
const
G4ParticleDefinition
&)
63
{;
64
}
66
//
67
void
G4VAdjointReverseReaction::BuildPhysicsTable
(
const
G4ParticleDefinition
&)
68
{
69
70
theAdjointCSManager
->
BuildCrossSectionMatrices
();
//do not worry it will be done just once
71
theAdjointCSManager
->
BuildTotalSigmaTables
();
72
73
}
75
//
76
G4VParticleChange
*
G4VAdjointReverseReaction::PostStepDoIt
(
const
G4Track
& track,
const
G4Step
& )
77
{
78
79
80
81
fParticleChange
->
Initialize
(track);
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
95
theAdjointEMModel
->
SampleSecondaries
(track,
96
IsScatProjToProjCase
,
97
fParticleChange
);
98
99
ClearNumberOfInteractionLengthLeft
();
100
return
fParticleChange
;
101
102
103
104
}
106
//
107
G4double
G4VAdjointReverseReaction::GetMeanFreePath
(
const
G4Track
& track,
108
G4double
,
109
G4ForceCondition
*
condition
)
110
{ *condition =
NotForced
;
111
G4double
preStepKinEnergy = track.
GetKineticEnergy
();
112
113
/*G4double Sigma =
114
theAdjointEMModel->AdjointCrossSection(track.GetMaterialCutsCouple(),preStepKinEnergy,IsScatProjToProjCase);*/
115
116
G4double
Sigma =
117
theAdjointEMModel
->
GetAdjointCrossSection
(track.
GetMaterialCutsCouple
(),preStepKinEnergy,
IsScatProjToProjCase
);
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
}
Generated on Sat May 25 2013 14:33:29 for Geant4 by
1.8.4