source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPElasticData.cc @ 1340

Last change on this file since 1340 was 962, checked in by garnier, 15 years ago

update processes

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