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

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

update processes

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