source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPInelasticData.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.3 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 "G4NeutronHPInelasticData.hh"
35#include "G4Neutron.hh"
36#include "G4ElementTable.hh"
37#include "G4NeutronHPData.hh"
38
39G4bool G4NeutronHPInelasticData::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
47G4NeutronHPInelasticData::G4NeutronHPInelasticData()
48{
49// TKDB
50   theCrossSections = 0;
51  BuildPhysicsTable(*G4Neutron::Neutron());
52}
53   
54G4NeutronHPInelasticData::~G4NeutronHPInelasticData()
55{
56// TKDB
57  if ( theCrossSections != 0 )
58  {  theCrossSections->clearAndDestroy(); }
59  delete theCrossSections;
60}
61   
62void G4NeutronHPInelasticData::BuildPhysicsTable(const G4ParticleDefinition& aP)
63{
64  if(&aP!=G4Neutron::Neutron()) 
65     throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); 
66  size_t numberOfElements = G4Element::GetNumberOfElements();
67//  theCrossSections = new G4PhysicsTable( numberOfElements );
68// TKDB
69   if ( theCrossSections == 0 )
70   { 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 G4NeutronHPInelasticData::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 << "Inelastic 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 << "G4NeutronHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
130}
131
132#include "G4NucleiPropertiesTable.hh"
133
134G4double G4NeutronHPInelasticData::
135GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
136{
137  G4double result = 0;
138  G4bool outOfRange;
139  G4int index = anE->GetIndex();
140
141  // prepare neutron
142  G4double eKinetic = aP->GetKineticEnergy();
143
144  // T. K.
145  if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
146  {
147     G4double factor = 1.0;
148     if ( eKinetic < aT * k_Boltzmann ) 
149     {
150        // below 0.1 eV neutrons
151        // Have to do some, but now just igonre.   
152        // Will take care after performance check. 
153        // factor = factor * targetV;
154     }
155     return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor; 
156  }   
157
158  G4ReactionProduct theNeutron( aP->GetDefinition() );
159  theNeutron.SetMomentum( aP->GetMomentum() );
160  theNeutron.SetKineticEnergy( eKinetic );
161
162  // prepare thermal nucleus
163  G4Nucleus aNuc;
164  G4double eps = 0.0001;
165  G4double theA = anE->GetN();
166  G4double theZ = anE->GetZ();
167  G4double eleMass; 
168  eleMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theZ+eps), static_cast<G4int>(theA+eps))
169             ) / G4Neutron::Neutron()->GetPDGMass();
170 
171  G4ReactionProduct boosted;
172  G4double aXsection;
173 
174  // MC integration loop
175  G4int counter = 0;
176  G4int failCount = 0;
177  G4double buffer = 0;
178  G4int size = G4int(std::max(10., aT/60*kelvin));
179  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
180  G4double neutronVMag = neutronVelocity.mag();
181
182  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
183  {
184    if(counter) buffer = result/counter;
185    while (counter<size)
186    {
187      counter ++;
188      G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
189      boosted.Lorentz(theNeutron, aThermalNuc);
190      G4double theEkin = boosted.GetKineticEnergy();
191      aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
192      if(aXsection <0) 
193      {
194        if(failCount<1000)
195        {
196          failCount++;
197          counter--;
198          continue;
199        }
200        else
201        {
202          aXsection = 0;
203        }
204      }
205      // velocity correction.
206      G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
207      aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
208      result += aXsection;
209    }
210    size += size;
211  }
212  result /= counter;
213/*
214  // Checking impact of  G4NEUTRONHP_NEGLECT_DOPPLER
215  G4cout << " result " << result << " "
216         << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
217         << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
218*/
219  return result;
220}
Note: See TracBrowser for help on using the repository browser.