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