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 | // this code implementation is the intellectual property of |
---|
27 | // neutron_hp -- source file |
---|
28 | // J.P. Wellisch, Nov-1996 |
---|
29 | // A prototype of the low energy neutron transport model. |
---|
30 | // |
---|
31 | // By copying, distributing or modifying the Program (or any work |
---|
32 | // based on the Program) you indicate your acceptance of this statement, |
---|
33 | // and all its terms. |
---|
34 | // |
---|
35 | // $Id: G4NeutronHPInelastic.cc,v 1.23 2007/06/22 09:23:48 gcosmo Exp $ |
---|
36 | // GEANT4 tag $Name: $ |
---|
37 | // |
---|
38 | // 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi) |
---|
39 | // |
---|
40 | #include "G4NeutronHPInelastic.hh" |
---|
41 | |
---|
42 | G4NeutronHPInelastic::G4NeutronHPInelastic() |
---|
43 | :G4HadronicInteraction("NeutronHPInelastic") |
---|
44 | { |
---|
45 | SetMinEnergy( 0.0 ); |
---|
46 | SetMaxEnergy( 20.*MeV ); |
---|
47 | system("echo $G4NEUTRONHPDATA"); |
---|
48 | // G4cout << " entering G4NeutronHPInelastic constructor"<<G4endl; |
---|
49 | if(!getenv("G4NEUTRONHPDATA")) |
---|
50 | throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files."); |
---|
51 | dirName = getenv("G4NEUTRONHPDATA"); |
---|
52 | G4String tString = "/Inelastic/"; |
---|
53 | dirName = dirName + tString; |
---|
54 | numEle = G4Element::GetNumberOfElements(); |
---|
55 | theInelastic = new G4NeutronHPChannelList[numEle]; |
---|
56 | |
---|
57 | for (G4int i=0; i<numEle; i++) |
---|
58 | { |
---|
59 | theInelastic[i].Init((*(G4Element::GetElementTable()))[i], dirName); |
---|
60 | do |
---|
61 | { |
---|
62 | theInelastic[i].Register(&theNFS, "F01"); // has |
---|
63 | theInelastic[i].Register(&theNXFS, "F02"); |
---|
64 | theInelastic[i].Register(&the2NDFS, "F03"); |
---|
65 | theInelastic[i].Register(&the2NFS, "F04"); // has, E Done |
---|
66 | theInelastic[i].Register(&the3NFS, "F05"); // has, E Done |
---|
67 | theInelastic[i].Register(&theNAFS, "F06"); |
---|
68 | theInelastic[i].Register(&theN3AFS, "F07"); |
---|
69 | theInelastic[i].Register(&the2NAFS, "F08"); |
---|
70 | theInelastic[i].Register(&the3NAFS, "F09"); |
---|
71 | theInelastic[i].Register(&theNPFS, "F10"); |
---|
72 | theInelastic[i].Register(&theN2AFS, "F11"); |
---|
73 | theInelastic[i].Register(&the2N2AFS, "F12"); |
---|
74 | theInelastic[i].Register(&theNDFS, "F13"); |
---|
75 | theInelastic[i].Register(&theNTFS, "F14"); |
---|
76 | theInelastic[i].Register(&theNHe3FS, "F15"); |
---|
77 | theInelastic[i].Register(&theND2AFS, "F16"); |
---|
78 | theInelastic[i].Register(&theNT2AFS, "F17"); |
---|
79 | theInelastic[i].Register(&the4NFS, "F18"); // has, E Done |
---|
80 | theInelastic[i].Register(&the2NPFS, "F19"); |
---|
81 | theInelastic[i].Register(&the3NPFS, "F20"); |
---|
82 | theInelastic[i].Register(&theN2PFS, "F21"); |
---|
83 | theInelastic[i].Register(&theNPAFS, "F22"); |
---|
84 | theInelastic[i].Register(&thePFS, "F23"); |
---|
85 | theInelastic[i].Register(&theDFS, "F24"); |
---|
86 | theInelastic[i].Register(&theTFS, "F25"); |
---|
87 | theInelastic[i].Register(&theHe3FS, "F26"); |
---|
88 | theInelastic[i].Register(&theAFS, "F27"); |
---|
89 | theInelastic[i].Register(&the2AFS, "F28"); |
---|
90 | theInelastic[i].Register(&the3AFS, "F29"); |
---|
91 | theInelastic[i].Register(&the2PFS, "F30"); |
---|
92 | theInelastic[i].Register(&thePAFS, "F31"); |
---|
93 | theInelastic[i].Register(&theD2AFS, "F32"); |
---|
94 | theInelastic[i].Register(&theT2AFS, "F33"); |
---|
95 | theInelastic[i].Register(&thePDFS, "F34"); |
---|
96 | theInelastic[i].Register(&thePTFS, "F35"); |
---|
97 | theInelastic[i].Register(&theDAFS, "F36"); |
---|
98 | theInelastic[i].RestartRegistration(); |
---|
99 | } |
---|
100 | while(!theInelastic[i].HasDataInAnyFinalState()); |
---|
101 | } |
---|
102 | } |
---|
103 | |
---|
104 | G4NeutronHPInelastic::~G4NeutronHPInelastic() |
---|
105 | { |
---|
106 | delete [] theInelastic; |
---|
107 | } |
---|
108 | |
---|
109 | #include "G4NeutronHPThermalBoost.hh" |
---|
110 | |
---|
111 | G4HadFinalState * G4NeutronHPInelastic::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& ) |
---|
112 | { |
---|
113 | const G4Material * theMaterial = aTrack.GetMaterial(); |
---|
114 | G4int n = theMaterial->GetNumberOfElements(); |
---|
115 | G4int index = theMaterial->GetElement(0)->GetIndex(); |
---|
116 | G4int it=0; |
---|
117 | if(n!=1) |
---|
118 | { |
---|
119 | xSec = new G4double[n]; |
---|
120 | G4double sum=0; |
---|
121 | G4int i; |
---|
122 | const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume(); |
---|
123 | G4double rWeight; |
---|
124 | G4NeutronHPThermalBoost aThermalE; |
---|
125 | for (i=0; i<n; i++) |
---|
126 | { |
---|
127 | index = theMaterial->GetElement(i)->GetIndex(); |
---|
128 | rWeight = NumAtomsPerVolume[i]; |
---|
129 | xSec[i] = theInelastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack, |
---|
130 | theMaterial->GetElement(i), |
---|
131 | theMaterial->GetTemperature())); |
---|
132 | xSec[i] *= rWeight; |
---|
133 | sum+=xSec[i]; |
---|
134 | } |
---|
135 | G4double random = G4UniformRand(); |
---|
136 | G4double running = 0; |
---|
137 | for (i=0; i<n; i++) |
---|
138 | { |
---|
139 | running += xSec[i]; |
---|
140 | index = theMaterial->GetElement(i)->GetIndex(); |
---|
141 | it = i; |
---|
142 | //if(random<=running/sum) break; |
---|
143 | if( sum == 0 || random<=running/sum) break; |
---|
144 | } |
---|
145 | delete [] xSec; |
---|
146 | } |
---|
147 | return theInelastic[index].ApplyYourself(theMaterial->GetElement(it), aTrack); |
---|
148 | } |
---|