Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 "G4Pow.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4ElementTable.hh"
45
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 //G4cout << "G4LENDCrossSection::GetIsoIsIsoApplicable this->GetName() = " << this->GetName() << ", iZ = " << iZ << ", iA = " << iA << ", allow_nat = " << allow_nat << G4endl;
54 //Check existence of target data
55 if ( element != NULL ) {
56 if ( element->GetNumberOfIsotopes() != 0 ) {
57 std::vector< const G4Isotope*> vIsotope;
58 for ( G4int i = 0 ; i != (G4int)element->GetNumberOfIsotopes() ; ++i ) {
59 if ( element->GetIsotope( i )->GetN() == iA ) vIsotope.push_back( element->GetIsotope( i ) );
60 }
61 for ( std::size_t i = 0 ; i != vIsotope.size() ; ++i ) {
62 G4int iM = vIsotope[i]->Getm();
63 if ( get_target_from_map( lend_manager->GetNucleusEncoding( iZ , iA , iM ) ) != NULL ) return true;
64 }
65 //No isomer has data
66 //Check natural aboundance data for the element
67 if ( get_target_from_map( lend_manager->GetNucleusEncoding( iZ , 0 , 0 ) ) != NULL ) return true;
68 } else {
69 //Check for iZ and iA under assuming iM = 0
70 if ( get_target_from_map( lend_manager->GetNucleusEncoding( iZ , iA , 0 ) ) != NULL ) return true;
71 //Check natural aboundance data for the element
72 if ( get_target_from_map( lend_manager->GetNucleusEncoding( iZ , 0 , 0 ) ) != NULL ) return true;
73 }
74 } else {
75 //Check for iZ and iA under assuming iM = 0
76 if ( get_target_from_map( lend_manager->GetNucleusEncoding( iZ , iA , 0 ) ) != NULL ) return true;
77 //Check natural aboundance data for iZ
78 if ( get_target_from_map( lend_manager->GetNucleusEncoding( iZ , 0 , 0 ) ) != NULL ) return true;
79 }
80 return false;
81}
82
84 const G4Isotope* isotope , const G4Element* /*elment*/ , const G4Material* material )
85{
86
87 G4double xs = 0.0;
88 G4double ke = dp->GetKineticEnergy();
89 G4double temp = material->GetTemperature();
90 G4int iM = 0;
91 if ( isotope != NULL ) iM = isotope->Getm();
92
93 G4GIDI_target* aTarget = get_target_from_map( lend_manager->GetNucleusEncoding( iZ , iA , iM ) );
94 if ( aTarget != NULL ) {
95 xs = getLENDCrossSection ( aTarget , ke , temp );
96 }
97 else {
98 ;
99 /*
100 G4String message;
101 message = this->GetName();
102 message += " is unexpectedly called.";
103 //G4Exception( "G4LEND::GetIsoCrossSection(,)" , "LENDCrossSection-01" , JustWarning ,
104 G4Exception( "G4LEND::GetIsoCrossSection(,)" , "LENDCrossSection-01" , FatalException ,
105 message );
106 */
107 }
108
109 return xs;
110}
111
112
113/*
114G4bool G4LENDCrossSection::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
115{
116 G4bool result = true;
117 G4double eKin = aP->GetKineticEnergy();
118 if( eKin > GetMaxKinEnergy() || aP->GetDefinition() != proj ) result = false;
119 return result;
120}
121*/
122
125{
126
127 proj = NULL; //will be set in an inherited class
128 //default_evaluation = "endl99";
129 //default_evaluation = "ENDF.B-VII.0";
130 //default_evaluation = "ENDF/BVII.1";
131 //default_evaluation = "ENDF/B-8.0";
132 //default_evaluation = "ENDF/B-7.1";
133 default_evaluation = "";
134
135 allow_nat = false;
136 allow_any = false;
137
138 SetMinKinEnergy( 0*MeV );
139 SetMaxKinEnergy( 20*MeV );
140
141 lend_manager = G4LENDManager::GetInstance();
142
143}
144
146{
147
148 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
149 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
150 {
151 delete it->second;
152 }
153
154}
155
160
162{
163
164 if ( &aP != proj )
165 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use LEND data for particles other than neutrons!!!");
166
167 G4cout << G4endl;
168 G4cout << "Dump Cross Sections of " << GetName() << G4endl;
169 G4cout << "(Pointwise cross-section at 300 Kelvin.)" << G4endl;
170 G4cout << G4endl;
171
172 G4cout << "Target informaiton " << G4endl;
173
174 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
175 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
176 {
177 G4cout
178 << "Wanted " << it->second->GetWantedEvaluation()
179 << ", Z= " << it->second->GetWantedZ()
180 << ", A= " << it->second->GetWantedA()
181 << "; Actual " << it->second->GetActualEvaluation()
182 << ", Z= " << it->second->GetActualZ()
183 << ", A= " << it->second->GetActualA()
184 << ", " << it->second->GetTarget()
185 << G4endl;
186
187 G4int ie = 0;
188
189 G4GIDI_target* aTarget = it->second->GetTarget();
190 G4double aT = 300;
191 for ( ie = 0 ; ie < 130 ; ie++ )
192 {
193 G4double ke = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
194
195 if ( ke < 20*MeV )
196 {
197 G4cout << " " << GetName() << ", cross section at " << ke/eV << " [eV] = " << getLENDCrossSection ( aTarget , ke , aT )/barn << " [barn] " << G4endl;
198 }
199 }
200 G4cout << G4endl;
201
202 }
203
204}
205
206
207/*
208//110810
209//G4double G4LENDCrossSection::GetCrossSection(const G4DynamicParticle* aP , const G4Element* anElement , G4double aT)
210G4double G4LENDCrossSection::GetCrossSection(const G4DynamicParticle* aP , int iZ , const G4Material* aMat)
211{
212
213//110810
214 G4double aT = aMat->GetTemperature();
215 G4Element* anElement = lend_manager->GetNistElementBuilder()->FindOrBuildElement( iZ );
216
217 G4double ke = aP->GetKineticEnergy();
218 G4double XS = 0.0;
219
220 G4int numberOfIsotope = anElement->GetNumberOfIsotopes();
221
222 if ( numberOfIsotope > 0 )
223 {
224 // User Defined Abundances
225 for ( G4int i_iso = 0 ; i_iso < numberOfIsotope ; i_iso++ )
226 {
227
228 G4int iZ = anElement->GetIsotope( i_iso )->GetZ();
229 G4int iA = anElement->GetIsotope( i_iso )->GetN();
230 G4double ratio = anElement->GetRelativeAbundanceVector()[i_iso];
231
232 G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
233 XS += ratio*getLENDCrossSection ( aTarget , ke , aT );
234
235 }
236 }
237 else
238 {
239 // Natural Abundances
240 G4NistElementBuilder* nistElementBuild = lend_manager->GetNistElementBuilder();
241 G4int iZ = int ( anElement->GetZ() );
242 G4int numberOfNistIso = nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) );
243
244 G4int Nfirst = nistElementBuild->GetNistFirstIsotopeN( iZ );
245 for ( G4int i = 0 ; i < numberOfNistIso ; i++ )
246 {
247 G4int iA = Nfirst + i;
248 G4double ratio = nistElementBuild->GetIsotopeAbundance( iZ , iA );
249 if ( ratio > 0.0 )
250 {
251 G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
252 XS += ratio*getLENDCrossSection ( aTarget , ke , aT );
253 //G4cout << ke/eV << " " << iZ << " " << iMass << " " << aTarget << " " << getLENDCrossSection ( aTarget , ke , aT ) << G4endl;
254 }
255 }
256 }
257
258 //G4cout << "XS= " << XS << G4endl;
259 return XS;
260}
261
262
263
264//110810
265//G4double G4LENDCrossSection::GetIsoCrossSection(const G4DynamicParticle* dp, const G4Isotope* isotope, G4double aT )
266G4double G4LENDCrossSection::GetIsoCrossSection(const G4DynamicParticle* dp, const G4Isotope* isotope, const G4Material* aMat)
267{
268
269//110810
270 G4double aT = aMat->GetTemperature();
271
272 G4double ke = dp->GetKineticEnergy();
273
274 G4int iZ = isotope->GetZ();
275 G4int iA = isotope->GetN();
276
277 G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
278
279 return getLENDCrossSection ( aTarget , ke , aT );
280
281}
282
283
284
285//110810
286//G4double G4LENDCrossSection::GetZandACrossSection(const G4DynamicParticle* dp, G4int iZ, G4int iA, G4double aT)
287G4double G4LENDCrossSection::GetZandACrossSection(const G4DynamicParticle* dp, G4int iZ, G4int iA, const G4Material* aMat)
288{
289
290//110810
291 G4double aT = aMat->GetTemperature();
292
293 G4double ke = dp->GetKineticEnergy();
294
295 G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
296
297 return getLENDCrossSection ( aTarget , ke , aT );
298
299}
300*/
301
302
303
304void G4LENDCrossSection::recreate_used_target_map()
305{
306 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
307 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
308 {
309 delete it->second;
310 }
311 usedTarget_map.clear();
312
314}
315
316
317
319{
320
321 lend_manager->RequestChangeOfVerboseLevel( verboseLevel );
322
323 std::size_t numberOfElements = G4Element::GetNumberOfElements();
324 static const G4ElementTable* theElementTable = G4Element::GetElementTable();
325
326 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
327 {
328
329 const G4Element* anElement = (*theElementTable)[i];
330 G4int numberOfIsotope = (G4int)anElement->GetNumberOfIsotopes();
331
332 if ( numberOfIsotope > 0 )
333 {
334 // User Defined Abundances
335 for ( G4int i_iso = 0 ; i_iso < numberOfIsotope ; i_iso++ )
336 {
337 G4int iZ = anElement->GetIsotope( i_iso )->GetZ();
338 G4int iA = anElement->GetIsotope( i_iso )->GetN();
339 G4int iIsomer = anElement->GetIsotope( i_iso )->Getm();
340
341 //G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( G4Neutron::Neutron() , default_evaluation , iZ , iA );
342 G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iA , iIsomer );
343 if ( allow_nat == true ) aTarget->AllowNat();
344 if ( allow_any == true ) aTarget->AllowAny();
345 usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iA , iIsomer ) , aTarget ) );
346 }
347 }
348 else
349 {
350 // Natural Abundances
351 G4NistElementBuilder* nistElementBuild = lend_manager->GetNistElementBuilder();
352 G4int iZ = int ( anElement->GetZ() );
353 //G4cout << nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) ) << G4endl;
354 G4int numberOfNistIso = nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) );
355
356 for ( G4int ii = 0 ; ii < numberOfNistIso ; ii++ )
357 {
358 //G4cout << nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + i ) << G4endl;
359 if ( nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii ) > 0 )
360 {
361 G4int iMass = nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii;
362 //G4cout << iZ << " " << nistElementBuild->GetNistFirstIsotopeN( iZ ) + i << " " << nistElementBuild->GetIsotopeAbundance ( iZ , iMass ) << G4endl;
363 G4int iIsomer = 0;
364
365 G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iMass );
366 if ( allow_nat == true ) aTarget->AllowNat();
367 if ( allow_any == true ) aTarget->AllowAny();
368 usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iMass , iIsomer ) , aTarget ) );
369
370 }
371
372 }
373 }
374 }
376}
377
378 // elow ehigh xs_elow xs_ehigh ke (ke < elow)
380{
381 //XS propotinal to 1/v at low energy -> 1/root(E)
382 //XS = a * 1/root(E) + b
383 G4double a = ( y2 - y1 ) / ( 1/std::sqrt(x2) - 1/std::sqrt(x1) );
384 G4double b = y1 - a * 1/std::sqrt(x1);
385 G4double result = a * 1/std::sqrt(ke) + b;
386 return result;
387}
388
390 G4GIDI_target* target = NULL;
391 if ( usedTarget_map.find( nuclear_code ) != usedTarget_map.end() ) {
392 target = usedTarget_map.find( nuclear_code )->second->GetTarget();
393 }
394 return target;
395}
396
398
399 if ( lend_manager->GetVerboseLevel() >= 1 || force ) {
400 if ( usedTarget_map.size() == 0 ) create_used_target_map();
401 G4cout << "Dumping UsedTarget of " << GetName() << " for " << proj->GetParticleName() << G4endl;
402 G4cout << "Requested Evaluation, Z , A -> Actual Evaluation, Z , A(0=Nat) " << G4endl;
403 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
404 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ ) {
405 G4cout
406 << " " << it->second->GetWantedEvaluation()
407 << ", " << it->second->GetWantedZ()
408 << ", " << it->second->GetWantedA()
409 << " -> " << it->second->GetActualEvaluation()
410 << ", " << it->second->GetActualZ()
411 << ", " << it->second->GetActualA()
412 << G4endl;
413 }
414 }
415}
416
std::vector< G4Element * > G4ElementTable
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition G4Element.hh:119
static std::size_t GetNumberOfElements()
Definition G4Element.cc:408
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
G4int GetZ() const
Definition G4Isotope.hh:80
G4int Getm() const
Definition G4Isotope.hh:89
G4int GetN() const
Definition G4Isotope.hh:83
G4LENDCrossSection(const G4String name="")
G4GIDI_target * get_target_from_map(G4int nuclear_code)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4ParticleDefinition * proj
G4double GetUltraLowEnergyExtrapolatedXS(G4double, G4double, G4double, G4double, G4double)
virtual G4double getLENDCrossSection(G4GIDI_target *, G4double, G4double)
void DumpPhysicsTable(const G4ParticleDefinition &)
void DumpLENDTargetInfo(G4bool force=false)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
static G4LENDManager * GetInstance()
G4double GetTemperature() const
G4int GetNumberOfNistIsotopes(G4int Z) const
G4double GetIsotopeAbundance(G4int Z, G4int N) const
G4int GetNistFirstIsotopeN(G4int Z) const
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
const G4String & GetName() const