source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPCaptureData.cc @ 829

Last change on this file since 829 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 7.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// 070523 add neglecting doppler broadening on the fly. T. Koi
31// 070613 fix memory leaking by T. Koi
32// 071002 enable cross section dump by T. Koi
33//
34#include "G4NeutronHPCaptureData.hh"
35#include "G4Neutron.hh"
36#include "G4ElementTable.hh"
37#include "G4NeutronHPData.hh"
38
39G4bool G4NeutronHPCaptureData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
40{
41  G4bool result = true;
42  G4double eKin = aP->GetKineticEnergy();
43  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
44  return result;
45}
46
47G4NeutronHPCaptureData::G4NeutronHPCaptureData()
48{
49// TKDB
50   theCrossSections = 0;
51  BuildPhysicsTable(*G4Neutron::Neutron());
52}
53   
54G4NeutronHPCaptureData::~G4NeutronHPCaptureData()
55{
56// TKDB
57   if ( theCrossSections != 0 )
58      theCrossSections->clearAndDestroy();
59
60  delete theCrossSections;
61}
62   
63void G4NeutronHPCaptureData::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  // G4cout << "CALLED G4NeutronHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
69   // TKDB
70   if ( theCrossSections == 0 ) 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     if(getenv("CaptureDataIndexDebug"))
78     {
79       G4int index_debug = ((*theElementTable)[i])->GetIndex();
80       G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
81     }
82     G4PhysicsVector* physVec = G4NeutronHPData::
83      Instance()->MakePhysicsVector((*theElementTable)[i], this);
84     theCrossSections->push_back(physVec);
85  }
86}
87
88void G4NeutronHPCaptureData::DumpPhysicsTable(const G4ParticleDefinition& aP)
89{
90  if(&aP!=G4Neutron::Neutron()) 
91     throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); 
92
93//
94// Dump element based cross section
95// range 10e-5 eV to 20 MeV
96// 10 point per decade
97// in barn
98//
99
100   G4cout << G4endl;
101   G4cout << G4endl;
102   G4cout << "Capture Cross Section of Neutron HP"<< G4endl;
103   G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
104   G4cout << G4endl;
105   G4cout << "Name of Element" << G4endl;
106   G4cout << "Energy[eV]  XS[barn]" << G4endl;
107   G4cout << G4endl;
108
109   size_t numberOfElements = G4Element::GetNumberOfElements();
110   static const G4ElementTable *theElementTable = G4Element::GetElementTable();
111
112   for ( size_t i = 0 ; i < numberOfElements ; ++i )
113   {
114
115      G4cout << (*theElementTable)[i]->GetName() << G4endl;
116
117      G4int ie = 0;
118
119      for ( ie = 0 ; ie < 130 ; ie++ )
120      {
121         G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
122         G4bool outOfRange = false;
123
124         if ( eKinetic < 20*MeV )
125         {
126            G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
127         }
128
129      }
130
131      G4cout << G4endl;
132   }
133
134
135//  G4cout << "G4NeutronHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
136}
137
138#include "G4NucleiPropertiesTable.hh"
139
140G4double G4NeutronHPCaptureData::
141GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
142{
143  G4double result = 0;
144  G4bool outOfRange;
145  G4int index = anE->GetIndex();
146
147  // prepare neutron
148  G4double eKinetic = aP->GetKineticEnergy();
149
150  if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
151  {
152     G4double factor = 1.0;
153     if ( eKinetic < aT * k_Boltzmann ) 
154     {
155        // below 0.1 eV neutrons
156        // Have to do some, but now just igonre.   
157        // Will take care after performance check. 
158        // factor = factor * targetV;
159     }
160     return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor; 
161  }
162
163  G4ReactionProduct theNeutron( aP->GetDefinition() );
164  theNeutron.SetMomentum( aP->GetMomentum() );
165  theNeutron.SetKineticEnergy( eKinetic );
166
167  // prepare thermal nucleus
168  G4Nucleus aNuc;
169  G4double eps = 0.0001;
170  G4double theA = anE->GetN();
171  G4double theZ = anE->GetZ();
172  G4double eleMass; 
173  eleMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theZ+eps), static_cast<G4int>(theA+eps))
174             ) / G4Neutron::Neutron()->GetPDGMass();
175 
176  G4ReactionProduct boosted;
177  G4double aXsection;
178 
179  // MC integration loop
180  G4int counter = 0;
181  G4double buffer = 0;
182  G4int size = G4int(std::max(10., aT/60*kelvin));
183  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
184  G4double neutronVMag = neutronVelocity.mag();
185
186  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
187  {
188    if(counter) buffer = result/counter;
189    while (counter<size)
190    {
191      counter ++;
192      G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
193      boosted.Lorentz(theNeutron, aThermalNuc);
194      G4double theEkin = boosted.GetKineticEnergy();
195      aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
196      // velocity correction, or luminosity factor...
197      G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
198      aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
199      result += aXsection;
200    }
201    size += size;
202  }
203  result /= counter;
204/*
205  // Checking impact of  G4NEUTRONHP_NEGLECT_DOPPLER
206  G4cout << " result " << result << " "
207         << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
208         << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
209*/
210  return result;
211}
Note: See TracBrowser for help on using the repository browser.