source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPFissionData.cc@ 893

Last change on this file since 893 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.4 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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 070618 fix memory leaking by T. Koi
31// 071002 enable cross section dump by T. Koi
32
33#include "G4NeutronHPFissionData.hh"
34#include "G4Neutron.hh"
35#include "G4ElementTable.hh"
36#include "G4NeutronHPData.hh"
37
38G4bool G4NeutronHPFissionData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
39{
40 G4bool result = true;
41 G4double eKin = aP->GetKineticEnergy();
42 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
43 return result;
44}
45
46G4NeutronHPFissionData::G4NeutronHPFissionData()
47{
48 theCrossSections = 0;
49 BuildPhysicsTable(*G4Neutron::Neutron());
50}
51
52G4NeutronHPFissionData::~G4NeutronHPFissionData()
53{
54
55// TKDB
56 if ( theCrossSections != NULL )
57 {
58 theCrossSections->clearAndDestroy();
59 delete theCrossSections;
60 }
61}
62
63void G4NeutronHPFissionData::BuildPhysicsTable(const G4ParticleDefinition& aP)
64{
65 if(&aP!=G4Neutron::Neutron())
66 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
67 size_t numberOfElements = G4Element::GetNumberOfElements();
68 //theCrossSections = new G4PhysicsTable( numberOfElements );
69 // TKDB
70 if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
71
72 // make a PhysicsVector for each element
73
74 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
75 for( size_t i=0; i<numberOfElements; ++i )
76 {
77 G4PhysicsVector* physVec = G4NeutronHPData::
78 Instance()->MakePhysicsVector((*theElementTable)[i], this);
79 theCrossSections->push_back(physVec);
80 }
81}
82
83void G4NeutronHPFissionData::DumpPhysicsTable(const G4ParticleDefinition& aP)
84{
85 if(&aP!=G4Neutron::Neutron())
86 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
87
88//
89// Dump element based cross section
90// range 10e-5 eV to 20 MeV
91// 10 point per decade
92// in barn
93//
94
95 G4cout << G4endl;
96 G4cout << G4endl;
97 G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
98 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
99 G4cout << G4endl;
100 G4cout << "Name of Element" << G4endl;
101 G4cout << "Energy[eV] XS[barn]" << G4endl;
102 G4cout << G4endl;
103
104 size_t numberOfElements = G4Element::GetNumberOfElements();
105 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
106
107 for ( size_t i = 0 ; i < numberOfElements ; ++i )
108 {
109
110 G4cout << (*theElementTable)[i]->GetName() << G4endl;
111
112 G4int ie = 0;
113
114 for ( ie = 0 ; ie < 130 ; ie++ )
115 {
116 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
117 G4bool outOfRange = false;
118
119 if ( eKinetic < 20*MeV )
120 {
121 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
122 }
123
124 }
125
126 G4cout << G4endl;
127 }
128
129 //G4cout << "G4NeutronHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
130}
131
132#include "G4NucleiPropertiesTable.hh"
133
134G4double G4NeutronHPFissionData::
135GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
136{
137 G4double result = 0;
138 if(anE->GetZ()<90) return result;
139 G4bool outOfRange;
140 G4int index = anE->GetIndex();
141
142 // prepare neutron
143 G4double eKinetic = aP->GetKineticEnergy();
144 G4ReactionProduct theNeutron( aP->GetDefinition() );
145 theNeutron.SetMomentum( aP->GetMomentum() );
146 theNeutron.SetKineticEnergy( eKinetic );
147
148 // prepare thermal nucleus
149 G4Nucleus aNuc;
150 G4double eps = 0.0001;
151 G4double theA = anE->GetN();
152 G4double theZ = anE->GetZ();
153 G4double eleMass;
154 eleMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theZ+eps), static_cast<G4int>(theA+eps))
155 ) / G4Neutron::Neutron()->GetPDGMass();
156
157 G4ReactionProduct boosted;
158 G4double aXsection;
159
160 // MC integration loop
161 G4int counter = 0;
162 G4double buffer = 0;
163 G4int size = G4int(std::max(10., aT/60*kelvin));
164 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
165 G4double neutronVMag = neutronVelocity.mag();
166
167 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
168 {
169 if(counter) buffer = result/counter;
170 while (counter<size)
171 {
172 counter ++;
173 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
174 boosted.Lorentz(theNeutron, aThermalNuc);
175 G4double theEkin = boosted.GetKineticEnergy();
176 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
177 // velocity correction.
178 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
179 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
180 result += aXsection;
181 }
182 size += size;
183 }
184 result /= counter;
185 return result;
186}
Note: See TracBrowser for help on using the repository browser.