Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPChannelList Class Reference

#include <G4NeutronHPChannelList.hh>

Public Member Functions

 G4NeutronHPChannelList (G4int n)
 
 G4NeutronHPChannelList ()
 
void Init (G4int n)
 
 ~G4NeutronHPChannelList ()
 
G4HadFinalStateApplyYourself (const G4Element *theElement, const G4HadProjectile &aTrack)
 
void Init (G4Element *anElement, const G4String &dirName)
 
void Register (G4NeutronHPFinalState *theFS, const G4String &aName)
 
G4double GetXsec (G4double anEnergy)
 
G4int GetNumberOfChannels ()
 
G4bool HasDataInAnyFinalState ()
 
void RestartRegistration ()
 

Detailed Description

Definition at line 45 of file G4NeutronHPChannelList.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPChannelList() [1/2]

G4NeutronHPChannelList::G4NeutronHPChannelList ( G4int  n)

Definition at line 40 of file G4NeutronHPChannelList.cc.

41 {
42 nChannels = n;
43 theChannels = new G4NeutronHPChannel * [n];
44 allChannelsCreated = false;
45 theInitCount = 0;
46 }

◆ G4NeutronHPChannelList() [2/2]

G4NeutronHPChannelList::G4NeutronHPChannelList ( )

Definition at line 48 of file G4NeutronHPChannelList.cc.

49 {
50 nChannels = 0;
51 theChannels = 0;
52 allChannelsCreated = false;
53 theInitCount = 0;
54 }

◆ ~G4NeutronHPChannelList()

G4NeutronHPChannelList::~G4NeutronHPChannelList ( )

Definition at line 56 of file G4NeutronHPChannelList.cc.

57 {
58 if(theChannels!=0)
59 {
60 for(G4int i=0;i<nChannels; i++)
61 {
62 delete theChannels[i];
63 }
64 delete [] theChannels;
65 }
66 }
int G4int
Definition: G4Types.hh:66

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4NeutronHPChannelList::ApplyYourself ( const G4Element theElement,
const G4HadProjectile aTrack 
)

Definition at line 70 of file G4NeutronHPChannelList.cc.

71 {
73 G4int i, ii;
74 // decide on the isotope
75 G4int numberOfIsos(0);
76 for(ii=0; ii<nChannels; ii++)
77 {
78 numberOfIsos = theChannels[ii]->GetNiso();
79 if(numberOfIsos!=0) break;
80 }
81 G4double * running= new G4double [numberOfIsos];
82 running[0] = 0;
83 for(i=0;i<numberOfIsos; i++)
84 {
85 if(i!=0) running[i] = running[i-1];
86 for(ii=0; ii<nChannels; ii++)
87 {
88 if(theChannels[ii]->HasAnyData(i))
89 {
90 running[i] +=theChannels[ii]->GetWeightedXsec(aThermalE.GetThermalEnergy(aTrack,
91 theChannels[ii]->GetN(i),
92 theChannels[ii]->GetZ(i),
93 aTrack.GetMaterial()->GetTemperature()),
94 i);
95 }
96 }
97 }
98 G4int isotope=nChannels-1;
99 G4double random=G4UniformRand();
100 for(i=0;i<numberOfIsos; i++)
101 {
102 isotope = i;
103 //if(random<running[i]/running[numberOfIsos-1]) break;
104 if(running[numberOfIsos-1] != 0) if(random<running[i]/running[numberOfIsos-1]) break;
105 }
106 delete [] running;
107
108 // decide on the channel
109 running = new G4double[nChannels];
110 running[0]=0;
111 G4int targA=-1; // For production of unChanged
112 G4int targZ=-1;
113 for(i=0; i<nChannels; i++)
114 {
115 if(i!=0) running[i] = running[i-1];
116 if(theChannels[i]->HasAnyData(isotope))
117 {
118 running[i] += theChannels[i]->GetFSCrossSection(aThermalE.GetThermalEnergy(aTrack,
119 theChannels[i]->GetN(isotope),
120 theChannels[i]->GetZ(isotope),
121 aTrack.GetMaterial()->GetTemperature()),
122 isotope);
123 targA=(G4int)theChannels[i]->GetN(isotope); //Will be simply used the last valid value
124 targZ=(G4int)theChannels[i]->GetZ(isotope);
125 }
126 }
127
128 //TK120607
129 if ( running[nChannels-1] == 0 )
130 {
131 //This happened usually by the miss match between the cross section data and model
132 if ( targA == -1 && targZ == -1 ) {
133 throw G4HadronicException(__FILE__, __LINE__, "NeutronHP model encounter lethal discrepancy with cross section data");
134 }
135
136 //TK121106
137 G4cout << "Warning from NeutronHP: could not find proper reaction channel. This may cause by inconsistency between cross section and model. Unchanged final states are returned." << G4endl;
138 unChanged.Clear();
139
140 //For Ep Check create unchanged final state including rest target
141 G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon ( targZ , targA , 0.0 );
142 G4DynamicParticle* targ_dp = new G4DynamicParticle( targ_pd , G4ThreeVector(1,0,0), 0.0 );
143 unChanged.SetEnergyChange(aTrack.GetKineticEnergy());
144 unChanged.SetMomentumChange(aTrack.Get4Momentum().vect() );
145 unChanged.AddSecondary(targ_dp);
146 //TK121106
149 return &unChanged;
150 }
151 //TK120607
152
153
154 G4int lChan=0;
155 random=G4UniformRand();
156 for(i=0; i<nChannels; i++)
157 {
158 lChan = i;
159 if(running[nChannels-1] != 0) if(random<running[i]/running[nChannels-1]) break;
160 }
161 delete [] running;
162 return theChannels[lChan]->ApplyYourself(aTrack, isotope);
163 }
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector vect() const
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTemperature() const
Definition: G4Material.hh:181
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)

Referenced by G4NeutronHPorLEInelastic::ApplyYourself().

◆ GetNumberOfChannels()

G4int G4NeutronHPChannelList::GetNumberOfChannels ( )
inline

Definition at line 75 of file G4NeutronHPChannelList.hh.

75{ return nChannels; }

◆ GetXsec()

G4double G4NeutronHPChannelList::GetXsec ( G4double  anEnergy)
inline

Definition at line 64 of file G4NeutronHPChannelList.hh.

65 {
66 G4double result=0;
67 G4int i;
68 for(i=0; i<nChannels; i++)
69 {
70 result+=std::max(0., theChannels[i]->GetXsec(anEnergy));
71 }
72 return result;
73 }
G4double GetXsec(G4double anEnergy)

Referenced by G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLEInelasticData::GetCrossSection(), and GetXsec().

◆ HasDataInAnyFinalState()

G4bool G4NeutronHPChannelList::HasDataInAnyFinalState ( )
inline

Definition at line 77 of file G4NeutronHPChannelList.hh.

78 {
79 G4bool result = false;
80 G4int i;
81 for(i=0; i<nChannels; i++)
82 {
83 if(theChannels[i]->HasDataInAnyFinalState()) result = true;
84 }
85 return result;
86 }
bool G4bool
Definition: G4Types.hh:67

Referenced by HasDataInAnyFinalState().

◆ Init() [1/2]

void G4NeutronHPChannelList::Init ( G4Element anElement,
const G4String dirName 
)

Definition at line 165 of file G4NeutronHPChannelList.cc.

166 {
167 theDir = dirName;
168// G4cout << theDir << G4endl;
169 theElement = anElement;
170// G4cout << theElement << G4endl;
171 ;
172 }

◆ Init() [2/2]

void G4NeutronHPChannelList::Init ( G4int  n)

◆ Register()

void G4NeutronHPChannelList::Register ( G4NeutronHPFinalState theFS,
const G4String aName 
)

Definition at line 174 of file G4NeutronHPChannelList.cc.

176 {
177 if(!allChannelsCreated)
178 {
179 if(nChannels!=0)
180 {
181 G4NeutronHPChannel ** theBuffer = new G4NeutronHPChannel * [nChannels+1];
182 G4int i;
183 for(i=0; i<nChannels; i++)
184 {
185 theBuffer[i] = theChannels[i];
186 }
187 delete [] theChannels;
188 theChannels = theBuffer;
189 }
190 else
191 {
192 theChannels = new G4NeutronHPChannel * [nChannels+1];
193 }
194 G4String name;
195 name = aName+"/";
196 theChannels[nChannels] = new G4NeutronHPChannel;
197 theChannels[nChannels]->Init(theElement, theDir, name);
198 nChannels++;
199 }
200
201 //110527TKDB Unnessary codes, Detected by gcc4.6 compiler
202 //G4bool result;
203 //result = theChannels[theInitCount]->Register(theFS);
204 theChannels[theInitCount]->Register(theFS);
205 theInitCount++;
206 }
G4bool Register(G4NeutronHPFinalState *theFS)
void Init(G4Element *theElement, const G4String dirName)

Referenced by G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic().

◆ RestartRegistration()

void G4NeutronHPChannelList::RestartRegistration ( )
inline

Definition at line 87 of file G4NeutronHPChannelList.hh.

88 {
89 allChannelsCreated = true;
90 theInitCount = 0;
91 }

Referenced by G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic().


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