Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPCaptureFS.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 // neutron_hp -- source file
27 // J.P. Wellisch, Nov-1996
28 // A prototype of the low energy neutron transport model.
29 //
30 // 12-April-06 Enable IC electron emissions T. Koi
31 // 26-January-07 Add G4NEUTRONHP_USE_ONLY_PHOTONEVAPORATION flag
32 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
33 // 101203 Bugzilla/Geant4 Problem 1155 Lack of residual in some case
34 // 110430 Temporary solution in the case of being MF6 final state in Capture reaction (MT102)
35 //
36 #include "G4NeutronHPCaptureFS.hh"
37 #include "G4PhysicalConstants.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4Gamma.hh"
40 #include "G4ReactionProduct.hh"
41 #include "G4Nucleus.hh"
42 #include "G4PhotonEvaporation.hh"
43 #include "G4Fragment.hh"
44 #include "G4ParticleTable.hh"
45 #include "G4NeutronHPDataUsed.hh"
46 
48  {
49 
50  G4int i;
51  theResult.Clear();
52 // prepare neutron
53  G4double eKinetic = theTrack.GetKineticEnergy();
54  const G4HadProjectile *incidentParticle = &theTrack;
55  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
56  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
57  theNeutron.SetKineticEnergy( eKinetic );
58 
59 // prepare target
61  G4Nucleus aNucleus;
62  G4double eps = 0.0001;
63  if(targetMass<500*MeV)
64  targetMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theBaseA+eps) , static_cast<G4int>(theBaseZ+eps) )) /
66  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
67  G4double temperature = theTrack.GetMaterial()->GetTemperature();
68  theTarget = aNucleus.GetBiasedThermalNucleus(targetMass, neutronVelocity, temperature);
69 
70 // go to nucleus rest system
71  theNeutron.Lorentz(theNeutron, -1*theTarget);
72  eKinetic = theNeutron.GetKineticEnergy();
73 
74 // dice the photons
75 
76  G4ReactionProductVector * thePhotons = 0;
77  if ( HasFSData() && !getenv ( "G4NEUTRONHP_USE_ONLY_PHOTONEVAPORATION" ) )
78  {
79  //TK110430
80  if ( hasExactMF6 )
81  {
82  theMF6FinalState.SetTarget(theTarget);
83  theMF6FinalState.SetNeutron(theNeutron);
84  thePhotons = theMF6FinalState.Sample( eKinetic );
85  }
86  else
87  thePhotons = theFinalStatePhotons.GetPhotons(eKinetic);
88  }
89  else
90  {
91  G4ThreeVector aCMSMomentum = theNeutron.GetMomentum()+theTarget.GetMomentum();
92  G4LorentzVector p4(aCMSMomentum, theTarget.GetTotalEnergy() + theNeutron.GetTotalEnergy());
93  G4Fragment nucleus(static_cast<G4int>(theBaseA+1), static_cast<G4int>(theBaseZ) ,p4);
94  G4PhotonEvaporation photonEvaporation;
95  // T. K. add
96  photonEvaporation.SetICM( TRUE );
97  G4FragmentVector* products = photonEvaporation.BreakItUp(nucleus);
98  G4FragmentVector::iterator it;
99  thePhotons = new G4ReactionProductVector;
100  for(it=products->begin(); it!=products->end(); it++)
101  {
102  G4ReactionProduct * theOne = new G4ReactionProduct;
103  // T. K. add
104  if ( (*it)->GetParticleDefinition() != 0 )
105  theOne->SetDefinition( (*it)->GetParticleDefinition() );
106  else
107  theOne->SetDefinition( G4Gamma::Gamma() ); // this definiion will be over writen
108 
109  // T. K. comment out below line
110  //theOne->SetDefinition( G4Gamma::Gamma() );
112  if( (*it)->GetMomentum().mag() > 10*MeV ) theOne->SetDefinition( theTable->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ)) );
113 
114  //if ( (*i)->GetExcitationEnergy() > 0 )
115  if ( (*it)->GetExcitationEnergy() > 1.0e-2*eV )
116  {
117  G4double ex = (*it)->GetExcitationEnergy();
118  G4ReactionProduct* aPhoton = new G4ReactionProduct;
119  aPhoton->SetDefinition( G4Gamma::Gamma() );
120  aPhoton->SetMomentum( (*it)->GetMomentum().vect().unit() * ex );
121  //aPhoton->SetTotalEnergy( ex ); //will be calculated from momentum
122  thePhotons->push_back(aPhoton);
123  }
124 
125  theOne->SetMomentum( (*it)->GetMomentum().vect() * ( (*it)->GetMomentum().t() - (*it)->GetExcitationEnergy() ) / (*it)->GetMomentum().t() ) ;
126  //theOne->SetTotalEnergy( (*i)->GetMomentum().t() - (*i)->GetExcitationEnergy() ); //will be calculated from momentum
127  thePhotons->push_back(theOne);
128  delete *it;
129  }
130  delete products;
131  }
132 
133 
134 
135 // add them to the final state
136 
137  G4int nPhotons = 0;
138  if(thePhotons!=0) nPhotons=thePhotons->size();
139 
140 
141  //110527TKDB Unused codes, Detected by gcc4.6 compiler
142  //G4int nParticles = nPhotons;
143  //if(1==nPhotons) nParticles = 2;
144 
145 //Make at least one photon
146 //101203 TK
147  if ( nPhotons == 0 )
148  {
149  G4ReactionProduct * theOne = new G4ReactionProduct;
150  theOne->SetDefinition( G4Gamma::Gamma() );
151  G4double theta = pi*G4UniformRand();
152  G4double phi = twopi*G4UniformRand();
153  G4double sinth = std::sin(theta);
154  G4ThreeVector direction( sinth*std::cos(phi), sinth*std::sin(phi), std::cos(theta) );
155  theOne->SetMomentum( direction ) ;
156  thePhotons->push_back(theOne);
157  nPhotons++; // 0 -> 1
158  }
159 //One photon case: energy set to Q-value
160 //101203 TK
161  //if ( nPhotons == 1 )
162  if ( nPhotons == 1 && thePhotons->operator[](0)->GetDefinition()->GetBaryonNumber() == 0 )
163  {
164  G4ThreeVector direction = thePhotons->operator[](0)->GetMomentum().unit();
165  G4double Q = G4ParticleTable::GetParticleTable()->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA), 0, static_cast<G4int>(theBaseZ))->GetPDGMass() + G4Neutron::Neutron()->GetPDGMass()
166  - G4ParticleTable::GetParticleTable()->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ))->GetPDGMass();
167  thePhotons->operator[](0)->SetMomentum( Q*direction );
168  }
169 //
170 
171  // back to lab system
172  for(i=0; i<nPhotons; i++)
173  {
174  thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), theTarget);
175  }
176 
177  // Recoil, if only one gamma
178  //if (1==nPhotons)
179  if ( nPhotons == 1 && thePhotons->operator[](0)->GetDefinition()->GetBaryonNumber() == 0 )
180  {
181  G4DynamicParticle * theOne = new G4DynamicParticle;
183  ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ));
184  theOne->SetDefinition(aRecoil);
185  // Now energy;
186  // Can be done slightly better @
187  G4ThreeVector aMomentum = theTrack.Get4Momentum().vect()
188  +theTarget.GetMomentum()
189  -thePhotons->operator[](0)->GetMomentum();
190 
191  G4ThreeVector theMomUnit = aMomentum.unit();
192  G4double aKinEnergy = theTrack.GetKineticEnergy()
193  +theTarget.GetKineticEnergy(); // gammas come from Q-value
194  G4double theResMass = aRecoil->GetPDGMass();
195  G4double theResE = aRecoil->GetPDGMass()+aKinEnergy;
196  G4double theAbsMom = std::sqrt(theResE*theResE - theResMass*theResMass);
197  G4ThreeVector theMomentum = theAbsMom*theMomUnit;
198  theOne->SetMomentum(theMomentum);
199  theResult.AddSecondary(theOne);
200  }
201 
202  // Now fill in the gammas.
203  for(i=0; i<nPhotons; i++)
204  {
205  // back to lab system
206  G4DynamicParticle * theOne = new G4DynamicParticle;
207  theOne->SetDefinition(thePhotons->operator[](i)->GetDefinition());
208  theOne->SetMomentum(thePhotons->operator[](i)->GetMomentum());
209  theResult.AddSecondary(theOne);
210  delete thePhotons->operator[](i);
211  }
212  delete thePhotons;
213 
214 //101203TK
215  G4bool residual = false;
217  ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ));
218  for ( G4int j = 0 ; j != theResult.GetNumberOfSecondaries() ; j++ )
219  {
220  if ( theResult.GetSecondary(j)->GetParticle()->GetDefinition() == aRecoil ) residual = true;
221  }
222 
223  if ( residual == false )
224  {
225  //G4ParticleDefinition * aRecoil = G4ParticleTable::GetParticleTable()
226  // ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ));
227  G4int nNonZero = 0;
228  G4LorentzVector p_photons(0,0,0,0);
229  for ( G4int j = 0 ; j != theResult.GetNumberOfSecondaries() ; j++ )
230  {
231  p_photons += theResult.GetSecondary(j)->GetParticle()->Get4Momentum();
232  // To many 0 momentum photons -> Check PhotonDist
233  if ( theResult.GetSecondary(j)->GetParticle()->Get4Momentum().e() > 0 ) nNonZero++;
234  }
235 
236  // Can we include kinetic energy here?
237  G4double deltaE = ( theTrack.Get4Momentum().e() + theTarget.GetTotalEnergy() )
238  - ( p_photons.e() + aRecoil->GetPDGMass() );
239 
240 //Add photons
241  if ( nPhotons - nNonZero > 0 )
242  {
243  //G4cout << "TKDB G4NeutronHPCaptureFS::ApplyYourself we will create additional " << nPhotons - nNonZero << " photons" << G4endl;
244  std::vector<G4double> vRand;
245  vRand.push_back( 0.0 );
246  for ( G4int j = 0 ; j != nPhotons - nNonZero - 1 ; j++ )
247  {
248  vRand.push_back( G4UniformRand() );
249  }
250  vRand.push_back( 1.0 );
251  std::sort( vRand.begin(), vRand.end() );
252 
253  std::vector<G4double> vEPhoton;
254  for ( G4int j = 0 ; j < (G4int)vRand.size() - 1 ; j++ )
255  {
256  vEPhoton.push_back( deltaE * ( vRand[j+1] - vRand[j] ) );
257  }
258  std::sort( vEPhoton.begin(), vEPhoton.end() );
259 
260  for ( G4int j = 0 ; j < nPhotons - nNonZero - 1 ; j++ )
261  {
262  //Isotopic in LAB OK?
263  G4double theta = pi*G4UniformRand();
264  G4double phi = twopi*G4UniformRand();
265  G4double sinth = std::sin(theta);
266  G4double en = vEPhoton[j];
267  G4ThreeVector tempVector(en*sinth*std::cos(phi), en*sinth*std::sin(phi), en*std::cos(theta) );
268 
269  p_photons += G4LorentzVector ( tempVector, tempVector.mag() );
270  G4DynamicParticle * theOne = new G4DynamicParticle;
271  theOne->SetDefinition( G4Gamma::Gamma() );
272  theOne->SetMomentum( tempVector );
273  theResult.AddSecondary(theOne);
274  }
275 
276 // Add last photon
277  G4DynamicParticle * theOne = new G4DynamicParticle;
278  theOne->SetDefinition( G4Gamma::Gamma() );
279 // For better momentum conservation
280  G4ThreeVector lastPhoton = -p_photons.vect().unit()*vEPhoton.back();
281  p_photons += G4LorentzVector( lastPhoton , lastPhoton.mag() );
282  theOne->SetMomentum( lastPhoton );
283  theResult.AddSecondary(theOne);
284  }
285 
286 //Add residual
287  G4DynamicParticle * theOne = new G4DynamicParticle;
288  G4ThreeVector aMomentum = theTrack.Get4Momentum().vect() + theTarget.GetMomentum()
289  - p_photons.vect();
290  theOne->SetDefinition(aRecoil);
291  theOne->SetMomentum( aMomentum );
292  theResult.AddSecondary(theOne);
293 
294  }
295 //101203TK END
296 
297 // clean up the primary neutron
299  return &theResult;
300  }
301 
302 #include <sstream>
304  {
305 
306  //TK110430 BEGIN
307  std::stringstream ss;
308  ss << static_cast<G4int>(Z);
309  G4String sZ;
310  ss >> sZ;
311  ss.clear();
312  ss << static_cast<G4int>(A);
313  G4String sA;
314  ss >> sA;
315 
316  ss.clear();
317  G4String sM;
318  if ( M > 0 )
319  {
320  ss << "m";
321  ss << M;
322  ss >> sM;
323  ss.clear();
324  }
325 
326  G4String element_name = theNames.GetName( static_cast<G4int>(Z)-1 );
327  G4String filenameMF6 = dirName+"/FSMF6/"+sZ+"_"+sA+sM+"_"+element_name;
328  std::ifstream dummyIFS(filenameMF6, std::ios::in);
329  if ( dummyIFS.good() == true ) hasExactMF6=true;
330 
331  //TK110430 Just for checking
332  //ENDF-VII.0 no case (check done at 110430
333  /*
334  if ( hasExactMF6 == true )
335  {
336  G4String filename = dirName+"FS/"+sZ+"_"+sA+"_"+element_name;
337  std::ifstream dummyIFS(filename, std::ios::in);
338  if ( dummyIFS.good() == true )
339  {
340  G4cout << "TKDB Capture Both FS and FSMF6 are exist for Z = " << sZ << ", A = " << sA << G4endl;;
341  }
342  }
343  */
344 
345  //TK110430 Only use MF6MT102 which has exactly same A and Z
346  //Even _nat_ do not select and there is no _nat_ case in ENDF-VII.0
347  if ( hasExactMF6 == true )
348  {
349  std::ifstream theData(filenameMF6, std::ios::in);
350  theMF6FinalState.Init(theData);
351  theData.close();
352  return;
353  }
354  //TK110430 END
355 
356 
357  G4String tString = "/FS";
358  G4bool dbool;
359  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
360 
361  G4String filename = aFile.GetName();
362  SetAZMs( A, Z, M, aFile );
363  //theBaseA = A;
364  //theBaseZ = G4int(Z+.5);
365  if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
366  {
367  hasAnyData = false;
368  hasFSData = false;
369  hasXsec = false;
370  return;
371  }
372  std::ifstream theData(filename, std::ios::in);
373 
374  hasFSData = theFinalStatePhotons.InitMean(theData);
375  if(hasFSData)
376  {
377  targetMass = theFinalStatePhotons.GetTargetMass();
378  theFinalStatePhotons.InitAngular(theData);
379  theFinalStatePhotons.InitEnergies(theData);
380  }
381  theData.close();
382  }