Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LENDCrossSection.cc
Go to the documentation of this file.
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 // Class Description
27 // Cross Section for LEND (Low Energy Nuclear Data)
28 // LEND is Geant4 interface for GIDI (General Interaction Data Interface)
29 // which gives a discription of nuclear and atomic reactions, such as
30 // Binary collision cross sections
31 // Particle number multiplicity distributions of reaction products
32 // Energy and angular distributions of reaction products
33 // Derived calculational constants
34 // GIDI is developped at Lawrence Livermore National Laboratory
35 // Class Description - End
36 
37 // 071025 First implementation done by T. Koi (SLAC/SCCS)
38 // 101118 Name modifications for release T. Koi (SLAC/PPA)
39 
40 #include "G4LENDCrossSection.hh"
41 #include "G4SystemOfUnits.hh"
42 #include "G4ElementTable.hh"
43 #include "G4HadronicException.hh"
44 
45 //TK110811
47  const G4Element* /*element*/ , const G4Material* /*material*/ )
48 {
49  G4double eKin = dp->GetKineticEnergy();
50  if ( dp->GetDefinition() != proj ) return false;
51  if ( eKin > GetMaxKinEnergy() || eKin < GetMinKinEnergy() ) return false;
52 
53  return true;
54 }
55 
57  const G4Isotope* /*isotope*/ , const G4Element* /*elment*/ , const G4Material* material )
58 {
59 
60  G4double xs = 0.0;
61  G4double ke = dp->GetKineticEnergy();
62  G4double temp = material->GetTemperature();
63 
64  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
65 
66  xs = getLENDCrossSection ( aTarget , ke , temp );
67 
68  return xs;
69 }
70 
71 
72 /*
73 G4bool G4LENDCrossSection::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
74 {
75  G4bool result = true;
76  G4double eKin = aP->GetKineticEnergy();
77  if( eKin > GetMaxKinEnergy() || aP->GetDefinition() != proj ) result = false;
78  return result;
79 }
80 */
81 
84 {
85 
86  default_evaluation = "endl99";
87  allow_nat = false;
88  allow_any = false;
89 
90  SetMinKinEnergy( 0*MeV );
91  SetMaxKinEnergy( 20*MeV );
92 
93  lend_manager = G4LENDManager::GetInstance();
94 
95 }
96 
98 {
99 
100  for ( std::map< G4int , G4LENDUsedTarget* >::iterator
101  it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
102  {
103  delete it->second;
104  }
105 
106 }
107 
109 {
111 }
112 
114 {
115 
116  if ( &aP != proj )
117  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use LEND data for particles other than neutrons!!!");
118 
119  G4cout << G4endl;
120  G4cout << "Dump Cross Sections of " << GetName() << G4endl;
121  G4cout << "(Pointwise cross-section at 300 Kelvin.)" << G4endl;
122  G4cout << G4endl;
123 
124  G4cout << "Target informaiton " << G4endl;
125 
126  for ( std::map< G4int , G4LENDUsedTarget* >::iterator
127  it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
128  {
129  G4cout
130  << "Wanted " << it->second->GetWantedEvaluation()
131  << ", Z= " << it->second->GetWantedZ()
132  << ", A= " << it->second->GetWantedA()
133  << "; Actual " << it->second->GetActualEvaluation()
134  << ", Z= " << it->second->GetActualZ()
135  << ", A= " << it->second->GetActualA()
136  << ", " << it->second->GetTarget()
137  << G4endl;
138 
139  G4int ie = 0;
140 
141  G4GIDI_target* aTarget = it->second->GetTarget();
142  G4double aT = 300;
143  for ( ie = 0 ; ie < 130 ; ie++ )
144  {
145  G4double ke = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
146 
147  if ( ke < 20*MeV )
148  {
149  G4cout << " " << GetName() << ", cross section at " << ke/eV << " [eV] = " << getLENDCrossSection ( aTarget , ke , aT )/barn << " [barn] " << G4endl;
150  }
151  }
152  G4cout << G4endl;
153 
154  }
155 
156 }
157 
158 
159 /*
160 //110810
161 //G4double G4LENDCrossSection::GetCrossSection(const G4DynamicParticle* aP , const G4Element* anElement , G4double aT)
162 G4double G4LENDCrossSection::GetCrossSection(const G4DynamicParticle* aP , int iZ , const G4Material* aMat)
163 {
164 
165 //110810
166  G4double aT = aMat->GetTemperature();
167  G4Element* anElement = lend_manager->GetNistElementBuilder()->FindOrBuildElement( iZ );
168 
169  G4double ke = aP->GetKineticEnergy();
170  G4double XS = 0.0;
171 
172  G4int numberOfIsotope = anElement->GetNumberOfIsotopes();
173 
174  if ( numberOfIsotope > 0 )
175  {
176  // User Defined Abundances
177  for ( G4int i_iso = 0 ; i_iso < numberOfIsotope ; i_iso++ )
178  {
179 
180  G4int iZ = anElement->GetIsotope( i_iso )->GetZ();
181  G4int iA = anElement->GetIsotope( i_iso )->GetN();
182  G4double ratio = anElement->GetRelativeAbundanceVector()[i_iso];
183 
184  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
185  XS += ratio*getLENDCrossSection ( aTarget , ke , aT );
186 
187  }
188  }
189  else
190  {
191  // Natural Abundances
192  G4NistElementBuilder* nistElementBuild = lend_manager->GetNistElementBuilder();
193  G4int iZ = int ( anElement->GetZ() );
194  G4int numberOfNistIso = nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) );
195 
196  G4int Nfirst = nistElementBuild->GetNistFirstIsotopeN( iZ );
197  for ( G4int i = 0 ; i < numberOfNistIso ; i++ )
198  {
199  G4int iA = Nfirst + i;
200  G4double ratio = nistElementBuild->GetIsotopeAbundance( iZ , iA );
201  if ( ratio > 0.0 )
202  {
203  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
204  XS += ratio*getLENDCrossSection ( aTarget , ke , aT );
205  //G4cout << ke/eV << " " << iZ << " " << iMass << " " << aTarget << " " << getLENDCrossSection ( aTarget , ke , aT ) << G4endl;
206  }
207  }
208  }
209 
210  //G4cout << "XS= " << XS << G4endl;
211  return XS;
212 }
213 
214 
215 
216 //110810
217 //G4double G4LENDCrossSection::GetIsoCrossSection(const G4DynamicParticle* dp, const G4Isotope* isotope, G4double aT )
218 G4double G4LENDCrossSection::GetIsoCrossSection(const G4DynamicParticle* dp, const G4Isotope* isotope, const G4Material* aMat)
219 {
220 
221 //110810
222  G4double aT = aMat->GetTemperature();
223 
224  G4double ke = dp->GetKineticEnergy();
225 
226  G4int iZ = isotope->GetZ();
227  G4int iA = isotope->GetN();
228 
229  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
230 
231  return getLENDCrossSection ( aTarget , ke , aT );
232 
233 }
234 
235 
236 
237 //110810
238 //G4double G4LENDCrossSection::GetZandACrossSection(const G4DynamicParticle* dp, G4int iZ, G4int iA, G4double aT)
239 G4double G4LENDCrossSection::GetZandACrossSection(const G4DynamicParticle* dp, G4int iZ, G4int iA, const G4Material* aMat)
240 {
241 
242 //110810
243  G4double aT = aMat->GetTemperature();
244 
245  G4double ke = dp->GetKineticEnergy();
246 
247  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
248 
249  return getLENDCrossSection ( aTarget , ke , aT );
250 
251 }
252 */
253 
254 
255 
256 void G4LENDCrossSection::recreate_used_target_map()
257 {
258  for ( std::map< G4int , G4LENDUsedTarget* >::iterator
259  it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
260  {
261  delete it->second;
262  }
263  usedTarget_map.clear();
264 
266 }
267 
268 
269 
271 {
272 
273  lend_manager->RequestChangeOfVerboseLevel( verboseLevel );
274 
275  size_t numberOfElements = G4Element::GetNumberOfElements();
276  static const G4ElementTable* theElementTable = G4Element::GetElementTable();
277 
278  for ( size_t i = 0 ; i < numberOfElements ; ++i )
279  {
280 
281  const G4Element* anElement = (*theElementTable)[i];
282  G4int numberOfIsotope = anElement->GetNumberOfIsotopes();
283 
284  if ( numberOfIsotope > 0 )
285  {
286  // User Defined Abundances
287  for ( G4int i_iso = 0 ; i_iso < numberOfIsotope ; i_iso++ )
288  {
289  G4int iZ = anElement->GetIsotope( i_iso )->GetZ();
290  G4int iA = anElement->GetIsotope( i_iso )->GetN();
291 
292  //G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( G4Neutron::Neutron() , default_evaluation , iZ , iA );
293  G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iA );
294  if ( allow_nat == true ) aTarget->AllowNat();
295  if ( allow_any == true ) aTarget->AllowAny();
296  usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iA ) , aTarget ) );
297  }
298  }
299  else
300  {
301  // Natural Abundances
302  G4NistElementBuilder* nistElementBuild = lend_manager->GetNistElementBuilder();
303  G4int iZ = int ( anElement->GetZ() );
304  //G4cout << nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) ) << G4endl;
305  G4int numberOfNistIso = nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) );
306 
307  for ( G4int ii = 0 ; ii < numberOfNistIso ; ii++ )
308  {
309  //G4cout << nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + i ) << G4endl;
310  if ( nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii ) > 0 )
311  {
312  G4int iMass = nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii;
313  //G4cout << iZ << " " << nistElementBuild->GetNistFirstIsotopeN( iZ ) + i << " " << nistElementBuild->GetIsotopeAbundance ( iZ , iMass ) << G4endl;
314 
315  G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iMass );
316  if ( allow_nat == true ) aTarget->AllowNat();
317  if ( allow_any == true ) aTarget->AllowAny();
318  usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iMass ) , aTarget ) );
319 
320  }
321 
322  }
323  }
324  }
325 
326  G4cout << "Dump UsedTarget for " << GetName() << G4endl;
327  G4cout << "Requested Evaluation, Z , A -> Actual Evaluation, Z , A(0=Nat) , Pointer of Target" << G4endl;
328  for ( std::map< G4int , G4LENDUsedTarget* >::iterator
329  it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
330  {
331  G4cout
332  << " " << it->second->GetWantedEvaluation()
333  << ", " << it->second->GetWantedZ()
334  << ", " << it->second->GetWantedA()
335  << " -> " << it->second->GetActualEvaluation()
336  << ", " << it->second->GetActualZ()
337  << ", " << it->second->GetActualA()
338  << ", " << it->second->GetTarget()
339  << G4endl;
340  }
341 
342 }
343 
344  // elow ehigh xs_elow xs_ehigh ke (ke < elow)
346 {
347  //XS propotinal to 1/v at low energy -> 1/root(E)
348  //XS = a * 1/root(E) + b
349  G4double a = ( y2 - y1 ) / ( 1/std::sqrt(x2) - 1/std::sqrt(x1) );
350  G4double b = y1 - a * 1/std::sqrt(x1);
351  G4double result = a * 1/std::sqrt(ke) + b;
352  return result;
353 }