Geant4  10.01.p03
G4ParticleHPInelasticData.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 // particle_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 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
38 //
40 #include "G4ParticleHPManager.hh"
41 #include "G4Neutron.hh"
42 #include "G4ElementTable.hh"
43 #include "G4ParticleHPData.hh"
44 
47 {
48  const char* dataDirVariable;
49  if( projectile == G4Neutron::Neutron() ) {
50  dataDirVariable = "G4NEUTRONHPDATA";
51  }else if( projectile == G4Proton::Proton() ) {
52  dataDirVariable = "G4PROTONHPDATA";
53  }else if( projectile == G4Deuteron::Deuteron() ) {
54  dataDirVariable = "G4DEUTERONHPDATA";
55  }else if( projectile == G4Triton::Triton() ) {
56  dataDirVariable = "G4TRITONHPDATA";
57  }else if( projectile == G4He3::He3() ) {
58  dataDirVariable = "G4HE3HPDATA";
59  }else if( projectile == G4Alpha::Alpha() ) {
60  dataDirVariable = "G4ALPHAHPDATA";
61  } else {
62  G4String message("G4ParticleHPInelasticData may only be called for neutron, proton, deuteron, triton, He3 or alpha, while it is called for " + projectile->GetParticleName());
63  throw G4HadronicException(__FILE__, __LINE__,message.c_str());
64  }
65  // G4cout << this << " G4ParticleHPInelasticData::G4ParticleHPInelasticData " << projectile->GetParticleName() << " DATADIR " << dataDirVariable << G4endl;//GDEB
66  G4String dataName = projectile->GetParticleName()+"HPInelasticXS";
67  dataName.at(0) = toupper(dataName.at(0)) ;
68  SetName( dataName );
69 
70  if(!getenv(dataDirVariable)){
71  G4String message("Please setenv " + G4String(dataDirVariable) + " to point to the " + projectile->GetParticleName() + " cross-section files.");
72  throw G4HadronicException(__FILE__, __LINE__,message.c_str());
73  }
74 
75  G4cout << "@@@ G4ParticleHPInelasticData instantiated for particle " << projectile->GetParticleName() << " data directory variable is " << dataDirVariable << " pointing to " << getenv(dataDirVariable) << G4endl;
76 
79 
80  ke_cache = 0.0;
81  xs_cache = 0.0;
82  element_cache = NULL;
83  material_cache = NULL;
84 
85  onFlightDB = true;
86  theCrossSections = 0;
87  theProjectile=projectile;
88 
89 
90  BuildPhysicsTableHP(projectile, dataDirVariable);
91 
92 }
93 
95 {
96  if ( theCrossSections != NULL ) {
98  delete theCrossSections;
99  theCrossSections = NULL;
100  }
101 }
102 
104  G4int /*Z*/ , G4int /*A*/ ,
105  const G4Element* /*elm*/ ,
106  const G4Material* /*mat*/ )
107 {
108  G4double eKin = dp->GetKineticEnergy();
109  if ( eKin > GetMaxKinEnergy()
110  || eKin < GetMinKinEnergy()
111  || dp->GetDefinition() != theProjectile ) return false;
112 
113  return true;
114 }
115 
117  G4int /*Z*/ , G4int /*A*/ ,
118  const G4Isotope* /*iso*/ ,
119  const G4Element* element ,
120  const G4Material* material )
121 {
122  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
123 
124  ke_cache = dp->GetKineticEnergy();
125  element_cache = element;
126  material_cache = material;
127  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
128  xs_cache = xs;
129  return xs;
130  //return GetCrossSection( dp , element , material->GetTemperature() );
131 }
132 
133 /*
134 G4bool G4ParticleHPInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
135 {
136  G4bool result = true;
137  G4double eKin = aP->GetKineticEnergy();
138  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
139  return result;
140 }
141 */
142 
143 void G4ParticleHPInelasticData::BuildPhysicsTableHP(G4ParticleDefinition* projectile,const char* /* dataDirVariable */)
144 {
145  // if(&projectile!=G4Neutron::Neutron())
146  // throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
147 
148 //080428
149  if ( getenv( "G4PHP_NEGLECT_DOPPLER" ) )
150  {
151  G4cout << "Find environment variable of \"G4PHP_NEGLECT_DOPPLER\"." << G4endl;
152  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of inelastic scattering of neutrons (<20MeV)." << G4endl;
153  onFlightDB = false;
154  }
155 
156  size_t numberOfElements = G4Element::GetNumberOfElements();
157 // theCrossSections = new G4PhysicsTable( numberOfElements );
158 // TKDB
159  //if ( theCrossSections == 0 )
160  //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
161  if ( theCrossSections == NULL )
162  theCrossSections = new G4PhysicsTable( numberOfElements );
163  else
165 
166  // make a PhysicsVector for each element
167 
168  G4ParticleHPData* hpData = new G4ParticleHPData(projectile); //NEW
169  static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
170  if (!theElementTable) theElementTable= G4Element::GetElementTable();
171  for( size_t i=0; i<numberOfElements; ++i )
172  {
173  //NEW G4PhysicsVector* physVec = G4ParticleHPData::
174  //NEW Instance(projectile, dataDirVariable)->MakePhysicsVector((*theElementTable)[i], this);
175  G4PhysicsVector* physVec = hpData->MakePhysicsVector((*theElementTable)[i], this);
176  theCrossSections->push_back(physVec);
177  }
178 }
179 
181 {
182  if(&projectile!=theProjectile)
183  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use ParticleHP data for a wrong projectile!!!");
184 
185 //
186 // Dump element based cross section
187 // range 10e-5 eV to 20 MeV
188 // 10 point per decade
189 // in barn
190 //
191 
192  G4cout << G4endl;
193  G4cout << G4endl;
194  G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
195  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
196  G4cout << G4endl;
197  G4cout << "Name of Element" << G4endl;
198  G4cout << "Energy[eV] XS[barn]" << G4endl;
199  G4cout << G4endl;
200 
201  size_t numberOfElements = G4Element::GetNumberOfElements();
202  static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
203  if (!theElementTable) theElementTable= G4Element::GetElementTable();
204 
205  for ( size_t i = 0 ; i < numberOfElements ; ++i )
206  {
207 
208  G4cout << (*theElementTable)[i]->GetName() << G4endl;
209 
210  G4int ie = 0;
211 
212  for ( ie = 0 ; ie < 130 ; ie++ )
213  {
214  G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *CLHEP::eV;
215  G4bool outOfRange = false;
216 
217  if ( eKinetic < 20*CLHEP::MeV )
218  {
219  G4cout << eKinetic/CLHEP::eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/CLHEP::barn << G4endl;
220  }
221 
222  }
223 
224  G4cout << G4endl;
225  }
226 
227  //G4cout << "G4ParticleHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
228 }
229 
230 #include "G4NucleiProperties.hh"
231 
233 GetCrossSection(const G4DynamicParticle* projectile, const G4Element*anE, G4double aT)
234 {
235  G4double result = 0;
236  G4bool outOfRange;
237  G4int index = anE->GetIndex();
238 
239  // prepare neutron
240  G4double eKinetic = projectile->GetKineticEnergy();
241 
242  // T. K.
243 //if ( getenv( "G4PHP_NEGLECT_DOPPLER" ) )
244 //080428
245  if ( !onFlightDB )
246  {
247  G4double factor = 1.0;
248  if ( eKinetic < aT * CLHEP::k_Boltzmann )
249  {
250  // below 0.1 eV neutrons
251  // Have to do some, but now just igonre.
252  // Will take care after performance check.
253  // factor = factor * targetV;
254  }
255  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
256 
257  }
258 
259  G4ReactionProduct theNeutron( projectile->GetDefinition() );
260  theNeutron.SetMomentum( projectile->GetMomentum() );
261  theNeutron.SetKineticEnergy( eKinetic );
262 
263  // prepare thermal nucleus
264  G4Nucleus aNuc;
265  G4double eps = 0.0001;
266  G4double theA = anE->GetN();
267  G4double theZ = anE->GetZ();
268  G4double eleMass;
269  eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
270  ) / theProjectile->GetPDGMass();
271 
272  G4ReactionProduct boosted;
273  G4double aXsection;
274 
275  // MC integration loop
276  G4int counter = 0;
277  G4int failCount = 0;
278  G4double buffer = 0;
279  G4int size = G4int(std::max(10., aT/60*CLHEP::kelvin));
280  G4ThreeVector neutronVelocity = 1./theProjectile->GetPDGMass()*theNeutron.GetMomentum();
281  G4double neutronVMag = neutronVelocity.mag();
282 
283  // G4cout << " G4ParticleHPInelasticData 2 " << size << G4endl;//GDEB
284 #ifndef G4PHP_DOPPLER_LOOP_ONCE
285  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
286  {
287  if(counter) buffer = result/counter;
288  while (counter<size)
289  {
290  counter ++;
291 #endif
292  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
293  boosted.Lorentz(theNeutron, aThermalNuc);
294  G4double theEkin = boosted.GetKineticEnergy();
295  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
296  // G4cout << " G4ParticleHPInelasticData aXsection " << aXsection << " index " << index << " theEkin " << theEkin << " outOfRange " << outOfRange <<G4endl;//GDEB
297  if(aXsection <0)
298  {
299  if(failCount<1000)
300  {
301  failCount++;
302 #ifndef G4PHP_DOPPLER_LOOP_ONCE
303  counter--;
304  continue;
305 #endif
306  }
307  else
308  {
309  aXsection = 0;
310  }
311  }
312  // velocity correction.
313  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
314  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
315  result += aXsection;
316  }
317 #ifndef G4PHP_DOPPLER_LOOP_ONCE
318  size += size;
319  }
320  result /= counter;
321 #endif
322 
323 /*
324  // Checking impact of G4PHP_NEGLECT_DOPPLER
325  G4cout << " result " << result << " "
326  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
327  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
328 */
329 // G4cout << this << " G4ParticleHPInelasticData result " << result << G4endl; //GDEB
330 
331  return result;
332 }
333 
335 {
337 }
339 {
341 }
static G4ParticleHPManager * GetInstance()
static const double MeV
Definition: G4SIunits.hh:193
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetN() const
Definition: G4Element.hh:134
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4ParticleDefinition * theProjectile
void push_back(G4PhysicsVector *)
G4double GetZ() const
Definition: G4Element.hh:131
void DumpPhysicsTable(const G4ParticleDefinition &)
#define buffer
Definition: xmlparse.cc:611
static const G4double eps
G4ParticleDefinition * GetDefinition() const
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
#define G4ThreadLocal
Definition: tls.hh:89
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4ParticleHPInelasticData(G4ParticleDefinition *projectile=G4Neutron::Neutron())
void BuildPhysicsTableHP(G4ParticleDefinition *projectile, const char *dataDirVariable)
void SetName(const G4String &)
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
void SetMinKinEnergy(G4double value)
bool G4bool
Definition: G4Types.hh:79
static G4Triton * Triton()
Definition: G4Triton.cc:95
size_t GetIndex() const
Definition: G4Element.hh:181
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static const double kelvin
Definition: G4SIunits.hh:260
G4double GetKineticEnergy() const
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
static const double eV
Definition: G4SIunits.hh:194
static const G4double factor
G4double GetPDGMass() const
void SetMaxKinEnergy(G4double value)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static const double barn
Definition: G4SIunits.hh:95
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
static G4He3 * He3()
Definition: G4He3.cc:94
G4double GetMass() const
void clearAndDestroy()
G4ThreeVector GetMomentum() const