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