source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPorLEInelasticData.cc@ 830

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

import all except CVS

File size: 7.0 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//
27// 05-11-21 NeutronHP or Low Energy Parameterization Models
28// Implemented by T. Koi (SLAC/SCCS)
29// If NeutronHP data do not available for an element, then Low Energy
30// Parameterization models handle the interactions of the element.
31//
32
33#include "G4NeutronHPorLEInelasticData.hh"
34#include "G4Neutron.hh"
35#include "G4ElementTable.hh"
36#include "G4NeutronHPData.hh"
37
38#include "G4PhysicsVector.hh"
39
40
41
42G4NeutronHPorLEInelasticData::G4NeutronHPorLEInelasticData( G4NeutronHPChannelList* pChannel , std::set< G4String >* pSet )
43{
44 theInelasticChannel = pChannel;
45 unavailable_elements = pSet;
46 BuildPhysicsTable(*G4Neutron::Neutron());
47}
48
49
50
51G4bool G4NeutronHPorLEInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element* anElement)
52{
53 G4bool result = true;
54 G4double eKin = aP->GetKineticEnergy();
55 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
56 if ( unavailable_elements->find( anElement->GetName() ) != unavailable_elements->end() ) result = false;
57 return result;
58}
59
60G4NeutronHPorLEInelasticData::G4NeutronHPorLEInelasticData()
61{
62// BuildPhysicsTable(*G4Neutron::Neutron());
63}
64
65
66
67G4NeutronHPorLEInelasticData::~G4NeutronHPorLEInelasticData()
68{
69// delete theCrossSections;
70}
71
72
73#include "G4NeutronHPInelasticData.hh"
74#include "G4LPhysicsFreeVector.hh"
75//#include "G4NeutronHPElementData.hh"
76
77void G4NeutronHPorLEInelasticData::BuildPhysicsTable( const G4ParticleDefinition& aP )
78{
79 if( &aP!=G4Neutron::Neutron() )
80 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
81
82 size_t numberOfElements = G4Element::GetNumberOfElements();
83 theCrossSections = new G4PhysicsTable( numberOfElements );
84
85 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
86 for ( size_t i=0 ; i < numberOfElements; ++i )
87 {
88 G4PhysicsVector* thePhysVec = new G4LPhysicsFreeVector(0, 0, 0);
89
90 if ( unavailable_elements->find( (*theElementTable)[i]->GetName() ) == unavailable_elements->end() )
91 {
92
93 G4NeutronHPElementData* theElementData = new G4NeutronHPElementData();
94 theElementData->Init( (*theElementTable)[i] );
95
96 G4NeutronHPVector* theHPVector = theElementData->GetData( (G4NeutronHPInelasticData*)this );
97
98 G4int len = theHPVector->GetVectorLength();
99
100 if ( len!=0 )
101 {
102 G4double emin = theHPVector->GetX(0);
103 G4double emax = theHPVector->GetX(len-1);
104
105 G4LPhysicsFreeVector* aPhysVector= new G4LPhysicsFreeVector ( len , emin , emax );
106 for ( G4int i=0; i<len; i++ )
107 {
108 aPhysVector->PutValues( i , theHPVector->GetX(i) , theHPVector->GetY(i) );
109 }
110 delete thePhysVec;
111 thePhysVec = aPhysVector;
112 }
113
114 //G4PhysicsVector* physVec = G4NeutronHPData::
115 //Instance()->MakePhysicsVector((*theElementTable)[i], this);
116 //theCrossSections->push_back(physVec);
117 }
118
119 theCrossSections->push_back(thePhysVec);
120 }
121}
122
123
124
125void G4NeutronHPorLEInelasticData::DumpPhysicsTable(const G4ParticleDefinition& aP)
126{
127 if(&aP!=G4Neutron::Neutron())
128 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
129// G4cout << "G4NeutronHPorLEInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
130}
131
132
133
134#include "G4Nucleus.hh"
135#include "G4NucleiPropertiesTable.hh"
136#include "G4Neutron.hh"
137#include "G4Electron.hh"
138
139G4double G4NeutronHPorLEInelasticData::
140GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
141{
142
143 // G4cout << "Choice G4NeutronHPorLEInelasticData for element " << anE->GetName() << G4endl;
144 G4double result = 0;
145 //G4bool outOfRange;
146 G4int index = anE->GetIndex();
147
148 // prepare neutron
149 G4double eKinetic = aP->GetKineticEnergy();
150 G4ReactionProduct theNeutron( aP->GetDefinition() );
151 theNeutron.SetMomentum( aP->GetMomentum() );
152 theNeutron.SetKineticEnergy( eKinetic );
153
154 // prepare thermal nucleus
155 G4Nucleus aNuc;
156 G4double eps = 0.0001;
157 G4double theA = anE->GetN();
158 G4double theZ = anE->GetZ();
159 G4double eleMass;
160 eleMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theZ+eps), static_cast<G4int>(theA+eps))
161 ) / G4Neutron::Neutron()->GetPDGMass();
162
163 G4ReactionProduct boosted;
164 G4double aXsection;
165
166 // MC integration loop
167 G4int counter = 0;
168 G4double buffer = 0;
169 G4int size = G4int(std::max(10., aT/60*kelvin));
170 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
171 G4double neutronVMag = neutronVelocity.mag();
172 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
173 {
174 if(counter) buffer = result/counter;
175 while (counter<size)
176 {
177 counter ++;
178 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
179 boosted.Lorentz(theNeutron, aThermalNuc);
180 G4double theEkin = boosted.GetKineticEnergy();
181 //aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
182 aXsection = theInelasticChannel[index].GetXsec( theEkin );
183 // velocity correction.
184 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
185 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
186 result += aXsection;
187 }
188 size += size;
189 }
190 result /= counter;
191 //return result;
192 return result*barn;
193}
Note: See TracBrowser for help on using the repository browser.