Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPorLFissionData.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 //
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 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
32 //
33 
35 #include "G4SystemOfUnits.hh"
36 #include "G4Neutron.hh"
37 #include "G4ElementTable.hh"
38 #include "G4NeutronHPData.hh"
39 
40 #include "G4PhysicsVector.hh"
41 
42 
44 {
45  SetMinKinEnergy( 0*MeV );
46  SetMaxKinEnergy( 20*MeV );
47 
48  ke_cache = 0.0;
49  xs_cache = 0.0;
50  element_cache = NULL;
51  material_cache = NULL;
52 // BuildPhysicsTable(*G4Neutron::Neutron());
53 }
54 
56 {
57 // delete theCrossSections;
58 }
59 
61  G4int /*Z*/ , G4int /*A*/ ,
62  const G4Element* element ,
63  const G4Material* /*mat*/ )
64 {
65  G4double eKin = dp->GetKineticEnergy();
66  if ( eKin > GetMaxKinEnergy()
67  || eKin < GetMinKinEnergy()
68  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
69  if ( unavailable_elements->find( element->GetName() ) != unavailable_elements->end() ) return false;
70 
71  return true;
72 }
73 
75  G4int /*Z*/ , G4int /*A*/ ,
76  const G4Isotope* /*iso*/ ,
77  const G4Element* element ,
78  const G4Material* material )
79 {
80  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
81 
82  ke_cache = dp->GetKineticEnergy();
83  element_cache = element;
84  material_cache = material;
85  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
86  xs_cache = xs;
87  return xs;
88  //return GetCrossSection( dp , element , material->GetTemperature() );
89 }
90 
92 :G4VCrossSectionDataSet("NeutronHPorLFissionXS")
93 {
94  theFissionChannel = pChannel;
95  unavailable_elements = pSet;
96 
97  SetMinKinEnergy( 0*MeV );
98  SetMaxKinEnergy( 20*MeV );
99 
100  ke_cache = 0.0;
101  xs_cache = 0.0;
102  element_cache = NULL;
103  material_cache = NULL;
104 }
105 
106 /*
107 G4bool G4NeutronHPorLFissionData::IsApplicable(const G4DynamicParticle*aP, const G4Element* anElement)
108 {
109  G4bool result = true;
110  G4double eKin = aP->GetKineticEnergy();
111  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
112  if ( unavailable_elements->find( anElement->GetName() ) != unavailable_elements->end() ) result = false;
113  return result;
114 }
115 */
116 
117 
118 
120 {
121  if( &aP!=G4Neutron::Neutron() )
122  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
123 }
124 
125 
126 
128 {
129  if(&aP!=G4Neutron::Neutron())
130  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
131 // G4cout << "G4NeutronHPorLFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
132 }
133 
134 
135 
136 #include "G4Nucleus.hh"
137 #include "G4NucleiProperties.hh"
138 #include "G4Neutron.hh"
139 #include "G4Electron.hh"
140 
143 {
144  G4double result = 0;
145  if ( anE->GetZ() < 90 ) return result;
146 
147 // G4bool outOfRange;
148  G4int index = anE->GetIndex();
149 
150  // prepare neutron
151  G4double eKinetic = aP->GetKineticEnergy();
152  G4ReactionProduct theNeutron( aP->GetDefinition() );
153  theNeutron.SetMomentum( aP->GetMomentum() );
154  theNeutron.SetKineticEnergy( eKinetic );
155 
156  // prepare thermal nucleus
157  G4Nucleus aNuc;
158  G4double eps = 0.0001;
159  G4double theA = anE->GetN();
160  G4double theZ = anE->GetZ();
161  G4double eleMass;
162  eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
164 
165  G4ReactionProduct boosted;
166  G4double aXsection;
167 
168  // MC integration loop
169  G4int counter = 0;
170  G4double buffer = 0;
171  G4int size = G4int(std::max(10., aT/60*kelvin));
172  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
173  G4double neutronVMag = neutronVelocity.mag();
174  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
175  {
176  if(counter) buffer = result/counter;
177  while (counter<size)
178  {
179  counter ++;
180  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
181  boosted.Lorentz(theNeutron, aThermalNuc);
182  G4double theEkin = boosted.GetKineticEnergy();
183  //aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
184  aXsection = theFissionChannel[index].GetXsec( theEkin );
185  // velocity correction.
186  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
187  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
188  result += aXsection;
189  }
190  size += size;
191  }
192  result /= counter;
193  return result;
194 }