Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPInelasticCompFS.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 // 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
31 // 070606 bug fix and migrate to enable to Partial cases by T. Koi
32 // 080603 bug fix for Hadron Hyper News #932 by T. Koi
33 // 080612 bug fix contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #4,6
34 // 080717 bug fix of calculation of residual momentum by T. Koi
35 // 080801 protect negative avalable energy by T. Koi
36 // introduce theNDLDataA,Z which has A and Z of NDL data by T. Koi
37 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
38 // 090514 Fix bug in IC electron emission case
39 // Contribution from Chao Zhang (Chao.Zhang@usd.edu) and Dongming Mei(Dongming.Mei@usd.edu)
40 // 100406 "nothingWasKnownOnHadron=1" then sample mu isotropic in CM
41 // add two_body_reaction
42 // 100909 add safty
43 // 101111 add safty for _nat_ data case in Binary reaction, but break conservation
44 // 110430 add Reaction Q value and break up flag (MF3::QI and LR)
45 //
47 #include "G4PhysicalConstants.hh"
48 #include "G4SystemOfUnits.hh"
49 #include "G4Nucleus.hh"
50 #include "G4NucleiProperties.hh"
51 #include "G4He3.hh"
52 #include "G4Alpha.hh"
53 #include "G4Electron.hh"
54 #include "G4NeutronHPDataUsed.hh"
55 #include "G4ParticleTable.hh"
56 
58 {
59  // char the[100] = {""};
60  // std::ostrstream ost(the, 100, std::ios::out);
61  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
62  // G4String * aName = new G4String(the);
63  // std::ifstream from(*aName, std::ios::in);
64 
65  std::ostringstream ost;
66  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
67  G4String aName = ost.str();
68  std::ifstream from(aName, std::ios::in);
69 
70  if(!from) return; // no data found for this isotope
71  // std::ifstream theGammaData(*aName, std::ios::in);
72  std::ifstream theGammaData(aName, std::ios::in);
73 
74  theGammas.Init(theGammaData);
75  // delete aName;
76 }
77 
79 {
80  gammaPath = "/Inelastic/Gammas/";
81  if(!getenv("G4NEUTRONHPDATA"))
82  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
83  G4String tBase = getenv("G4NEUTRONHPDATA");
84  gammaPath = tBase+gammaPath;
85  G4String tString = dirName;
86  G4bool dbool;
87  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, tString, aFSType, dbool);
88  G4String filename = aFile.GetName();
89  SetAZMs( A, Z, M, aFile );
90  //theBaseA = aFile.GetA();
91  //theBaseZ = aFile.GetZ();
92  //theNDLDataA = (int)aFile.GetA();
93  //theNDLDataZ = aFile.GetZ();
94  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
95  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
96  {
97  if(getenv("NeutronHPNamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
98  hasAnyData = false;
99  hasFSData = false;
100  hasXsec = false;
101  return;
102  }
103  //theBaseA = A;
104  //theBaseZ = G4int(Z+.5);
105  std::ifstream theData(filename, std::ios::in);
106  if(!theData)
107  {
108  hasAnyData = false;
109  hasFSData = false;
110  hasXsec = false;
111  theData.close();
112  return;
113  }
114  // here we go
115  G4int infoType, dataType, dummy;
116  G4int sfType, it;
117  hasFSData = false;
118  while (theData >> infoType)
119  {
120  hasFSData = true;
121  theData >> dataType;
122  theData >> sfType >> dummy;
123  it = 50;
124  if(sfType>=600||(sfType<100&&sfType>=50)) it = sfType%50;
125  if(dataType==3)
126  {
127  //theData >> dummy >> dummy;
128  //TK110430
129  // QI and LR introudced since G4NDL3.15
130  G4double dqi;
131  G4int ilr;
132  theData >> dqi >> ilr;
133 
134  QI[ it ] = dqi*eV;
135  LR[ it ] = ilr;
136  theXsection[it] = new G4NeutronHPVector;
137  G4int total;
138  theData >> total;
139  theXsection[it]->Init(theData, total, eV);
140  //std::cout << theXsection[it]->GetXsec(1*MeV) << std::endl;
141  }
142  else if(dataType==4)
143  {
145  theAngularDistribution[it]->Init(theData);
146  }
147  else if(dataType==5)
148  {
150  theEnergyDistribution[it]->Init(theData);
151  }
152  else if(dataType==6)
153  {
155  theEnergyAngData[it]->Init(theData);
156  }
157  else if(dataType==12)
158  {
160  theFinalStatePhotons[it]->InitMean(theData);
161  }
162  else if(dataType==13)
163  {
165  theFinalStatePhotons[it]->InitPartials(theData);
166  }
167  else if(dataType==14)
168  {
169  theFinalStatePhotons[it]->InitAngular(theData);
170  }
171  else if(dataType==15)
172  {
173  theFinalStatePhotons[it]->InitEnergies(theData);
174  }
175  else
176  {
177  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4NeutronHPInelasticCompFS");
178  }
179  }
180  theData.close();
181 }
182 
184 {
185 
186 // it = 0 has without Photon
187  G4double running[50];
188  running[0] = 0;
189  unsigned int i;
190  for(i=0; i<50; i++)
191  {
192  if(i!=0) running[i]=running[i-1];
193  if(theXsection[i] != 0)
194  {
195  running[i] += std::max(0., theXsection[i]->GetXsec(eKinetic));
196  }
197  }
198  G4double random = G4UniformRand();
199  G4double sum = running[49];
200  G4int it = 50;
201  if(0!=sum)
202  {
203  G4int i0;
204  for(i0=0; i0<50; i0++)
205  {
206  it = i0;
207  if(random < running[i0]/sum) break;
208  }
209  }
210 //debug: it = 1;
211  return it;
212 }
213 
214 
215  //n,p,d,t,he3,a
217 {
218 
219 // prepare neutron
220  theResult.Clear();
221  G4double eKinetic = theTrack.GetKineticEnergy();
222  const G4HadProjectile *incidentParticle = &theTrack;
223  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
224  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
225  theNeutron.SetKineticEnergy( eKinetic );
226 
227 // prepare target
228  G4int i;
229  for(i=0; i<50; i++)
230  { if(theXsection[i] != 0) { break; } }
231 
232  G4double targetMass=0;
233  G4double eps = 0.0001;
234  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
236 // if(theEnergyAngData[i]!=0)
237 // targetMass = theEnergyAngData[i]->GetTargetMass();
238 // else if(theAngularDistribution[i]!=0)
239 // targetMass = theAngularDistribution[i]->GetTargetMass();
240 // else if(theFinalStatePhotons[50]!=0)
241 // targetMass = theFinalStatePhotons[50]->GetTargetMass();
242  G4Nucleus aNucleus;
244  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
245  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
246 
247 // prepare the residual mass
248  G4double residualMass=0;
249  G4double residualZ = theBaseZ - aDefinition->GetPDGCharge();
250  G4double residualA = theBaseA - aDefinition->GetBaryonNumber()+1;
251  residualMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(residualA+eps), static_cast<G4int>(residualZ+eps)) ) /
253 
254 // prepare energy in target rest frame
255  G4ReactionProduct boosted;
256  boosted.Lorentz(theNeutron, theTarget);
257  eKinetic = boosted.GetKineticEnergy();
258 // G4double momentumInCMS = boosted.GetTotalMomentum();
259 
260 // select exit channel for composite FS class.
261  G4int it = SelectExitChannel( eKinetic );
262 
263 // set target and neutron in the relevant exit channel
264  InitDistributionInitialState(theNeutron, theTarget, it);
265 
266  G4ReactionProductVector * thePhotons = 0;
267  G4ReactionProductVector * theParticles = 0;
268  G4ReactionProduct aHadron;
269  aHadron.SetDefinition(aDefinition); // what if only cross-sections exist ==> Na 23 11 @@@@
270  G4double availableEnergy = theNeutron.GetKineticEnergy() + theNeutron.GetMass() - aHadron.GetMass() +
271  (targetMass - residualMass)*G4Neutron::Neutron()->GetPDGMass();
272 //080730c
273  if ( availableEnergy < 0 )
274  {
275  //G4cout << "080730c Adjust availavleEnergy " << G4endl;
276  availableEnergy = 0;
277  }
278  G4int nothingWasKnownOnHadron = 0;
279  G4int dummy;
280  G4double eGamm = 0;
281  G4int iLevel=it-1;
282 
283 // TK without photon has it = 0
284  if( 50 == it )
285  {
286 
287 // TK Excitation level is not determined
288  iLevel=-1;
289  aHadron.SetKineticEnergy(availableEnergy*residualMass*G4Neutron::Neutron()->GetPDGMass()/
290  (aHadron.GetMass()+residualMass*G4Neutron::Neutron()->GetPDGMass()));
291 
292  //aHadron.SetMomentum(theNeutron.GetMomentum()*(1./theNeutron.GetTotalMomentum())*
293  // std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
294  // aHadron.GetMass()*aHadron.GetMass()));
295 
296  //TK add safty 100909
297  G4double p2 = ( aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy() - aHadron.GetMass()*aHadron.GetMass() );
298  G4double p = 0.0;
299  if ( p2 > 0.0 ) p = std::sqrt( p );
300 
301  aHadron.SetMomentum(theNeutron.GetMomentum()*(1./theNeutron.GetTotalMomentum())*p );
302 
303  }
304  else
305  {
306  while( iLevel!=-1 && theGammas.GetLevel(iLevel) == 0 ) { iLevel--; }
307  }
308 
309 
310  if ( theAngularDistribution[it] != 0 ) // MF4
311  {
312  if(theEnergyDistribution[it]!=0) // MF5
313  {
314  //************************************************************
315  /*
316  aHadron.SetKineticEnergy(theEnergyDistribution[it]->Sample(eKinetic, dummy));
317  G4double eSecN = aHadron.GetKineticEnergy();
318  */
319  //************************************************************
320  //EMendoza --> maximum allowable energy should be taken into account.
321  G4double dqi = 0.0;
322  if ( QI[it] < 0 || 849 < QI[it] ) dqi = QI[it]; //For backword compatibility QI introduced since G4NDL3.15
323  G4double MaxEne=eKinetic+dqi;
324  G4double eSecN;
325  do{
326  eSecN=theEnergyDistribution[it]->Sample(eKinetic, dummy);
327  }while(eSecN>MaxEne);
328  aHadron.SetKineticEnergy(eSecN);
329  //************************************************************
330  eGamm = eKinetic-eSecN;
331  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
332  {
333  if(theGammas.GetLevelEnergy(iLevel)<eGamm) break;
334  }
335  G4double random = 2*G4UniformRand();
336  iLevel+=G4int(random);
337  if(iLevel>theGammas.GetNumberOfLevels()-1)iLevel = theGammas.GetNumberOfLevels()-1;
338  }
339  else
340  {
341  G4double eExcitation = 0;
342  if(iLevel>=0) eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();
343  while (eKinetic-eExcitation < 0 && iLevel>0)
344  {
345  iLevel--;
346  eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();
347  }
348  //110610TK BEGIN
349  //Use QI value for calculating excitation energy of residual.
350  G4bool useQI=false;
351  G4double dqi = QI[it];
352  if ( dqi < 0 || 849 < dqi ) useQI = true; //Former libraies does not have values of this range
353 
354  if ( useQI )
355  {
356  // QI introudced since G4NDL3.15
357  eExcitation = -QI[it];
358  //Re-evluate iLevel based on this eExcitation
359  iLevel = 0;
360  G4bool find = false;
361  G4int imaxEx = 0;
362  while( !theGammas.GetLevel(iLevel+1) == 0 )
363  {
364  G4double maxEx = 0.0;
365  if ( maxEx < theGammas.GetLevel(iLevel)->GetLevelEnergy() )
366  {
367  maxEx = theGammas.GetLevel(iLevel)->GetLevelEnergy();
368  imaxEx = iLevel;
369  }
370  if ( eExcitation < theGammas.GetLevel(iLevel)->GetLevelEnergy() )
371  {
372  find = true;
373  iLevel--;
374  // very small eExcitation, iLevel becomes -1, this is protected below.
375  if ( iLevel == -1 ) iLevel = 0; // But cause energy trouble.
376  break;
377  }
378  iLevel++;
379  }
380  // In case, cannot find proper level, then use the maximum level.
381  if ( !find ) iLevel = imaxEx;
382  }
383  //110610TK END
384 
385  if(getenv("InelasticCompFSLogging") && eKinetic-eExcitation < 0)
386  {
387  throw G4HadronicException(__FILE__, __LINE__, "SEVERE: InelasticCompFS: Consistency of data not good enough, please file report");
388  }
389  if(eKinetic-eExcitation < 0) eExcitation = 0;
390  if(iLevel!= -1) aHadron.SetKineticEnergy(eKinetic - eExcitation);
391 
392  }
394 
395  if( theFinalStatePhotons[it] == 0 )
396  {
397  //G4cout << "110610 USE Gamma Level" << G4endl;
398 // TK comment Most n,n* eneter to this
399  thePhotons = theGammas.GetDecayGammas(iLevel);
400  eGamm -= theGammas.GetLevelEnergy(iLevel);
401  if(eGamm>0) // @ ok for now, but really needs an efficient way of correllated sampling @
402  {
403  G4ReactionProduct * theRestEnergy = new G4ReactionProduct;
404  theRestEnergy->SetDefinition(G4Gamma::Gamma());
405  theRestEnergy->SetKineticEnergy(eGamm);
406  G4double costh = 2.*G4UniformRand()-1.;
407  G4double phi = twopi*G4UniformRand();
408  theRestEnergy->SetMomentum(eGamm*std::sin(std::acos(costh))*std::cos(phi),
409  eGamm*std::sin(std::acos(costh))*std::sin(phi),
410  eGamm*costh);
411  if(thePhotons == 0) { thePhotons = new G4ReactionProductVector; }
412  thePhotons->push_back(theRestEnergy);
413  }
414  }
415  }
416  else if(theEnergyAngData[it] != 0) // MF6
417  {
418  theParticles = theEnergyAngData[it]->Sample(eKinetic);
419  }
420  else
421  {
422  // @@@ what to do, if we have photon data, but no info on the hadron itself
423  nothingWasKnownOnHadron = 1;
424  }
425 
426  //G4cout << "theFinalStatePhotons it " << it << G4endl;
427  //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
428  //G4cout << "theFinalStatePhotons it " << it << G4endl;
429  //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
430  //G4cout << "thePhotons " << thePhotons << G4endl;
431 
432  if ( theFinalStatePhotons[it] != 0 )
433  {
434  // the photon distributions are in the Nucleus rest frame.
435  // TK residual rest frame
436  G4ReactionProduct boosted_tmp;
437  boosted_tmp.Lorentz(theNeutron, theTarget);
438  G4double anEnergy = boosted_tmp.GetKineticEnergy();
439  thePhotons = theFinalStatePhotons[it]->GetPhotons(anEnergy);
440  G4double aBaseEnergy = theFinalStatePhotons[it]->GetLevelEnergy();
441  G4double testEnergy = 0;
442  if(thePhotons!=0 && thePhotons->size()!=0)
443  { aBaseEnergy-=thePhotons->operator[](0)->GetTotalEnergy(); }
444  if(theFinalStatePhotons[it]->NeedsCascade())
445  {
446  while(aBaseEnergy>0.01*keV)
447  {
448  // cascade down the levels
449  G4bool foundMatchingLevel = false;
450  G4int closest = 2;
451  G4double deltaEold = -1;
452  for(G4int j=1; j<it; j++)
453  {
454  if(theFinalStatePhotons[j]!=0)
455  {
456  testEnergy = theFinalStatePhotons[j]->GetLevelEnergy();
457  }
458  else
459  {
460  testEnergy = 0;
461  }
462  G4double deltaE = std::abs(testEnergy-aBaseEnergy);
463  if(deltaE<0.1*keV)
464  {
465  G4ReactionProductVector * theNext =
466  theFinalStatePhotons[j]->GetPhotons(anEnergy);
467  thePhotons->push_back(theNext->operator[](0));
468  aBaseEnergy = testEnergy-theNext->operator[](0)->GetTotalEnergy();
469  delete theNext;
470  foundMatchingLevel = true;
471  break; // ===>
472  }
473  if(theFinalStatePhotons[j]!=0 && ( deltaE<deltaEold||deltaEold<0.) )
474  {
475  closest = j;
476  deltaEold = deltaE;
477  }
478  } // <=== the break goes here.
479  if(!foundMatchingLevel)
480  {
481  G4ReactionProductVector * theNext =
482  theFinalStatePhotons[closest]->GetPhotons(anEnergy);
483  thePhotons->push_back(theNext->operator[](0));
484  aBaseEnergy = aBaseEnergy-theNext->operator[](0)->GetTotalEnergy();
485  delete theNext;
486  }
487  }
488  }
489  }
490  unsigned int i0;
491  if(thePhotons!=0)
492  {
493  for(i0=0; i0<thePhotons->size(); i0++)
494  {
495  // back to lab
496  thePhotons->operator[](i0)->Lorentz(*(thePhotons->operator[](i0)), -1.*theTarget);
497  }
498  }
499  //G4cout << "nothingWasKnownOnHadron " << nothingWasKnownOnHadron << G4endl;
500  if(nothingWasKnownOnHadron)
501  {
502 // TKDB 100405
503 // In this case, hadron should be isotropic in CM
504 // mu and p should be correlated
505 //
506  G4double totalPhotonEnergy = 0.0;
507  if ( thePhotons != 0 )
508  {
509  unsigned int nPhotons = thePhotons->size();
510  unsigned int ii0;
511  for ( ii0=0; ii0<nPhotons; ii0++)
512  {
513  //thePhotons has energies at LAB system
514  totalPhotonEnergy += thePhotons->operator[](ii0)->GetTotalEnergy();
515  }
516  }
517 
518  //isotropic distribution in CM
519  G4double mu = 1.0 - 2 * G4UniformRand();
520 
521  // need momentums in target rest frame;
522  G4LorentzVector target_in_LAB ( theTarget.GetMomentum() , theTarget.GetTotalEnergy() );
523  G4ThreeVector boostToTargetRest = -target_in_LAB.boostVector();
524  G4LorentzVector proj_in_LAB = incidentParticle->Get4Momentum();
525 
526  G4DynamicParticle* proj = new G4DynamicParticle( G4Neutron::Neutron() , proj_in_LAB.boost( boostToTargetRest ) );
528  G4DynamicParticle* hadron = new G4DynamicParticle( aHadron.GetDefinition() , G4ThreeVector(0) ); // will be fill momentum
529 
530  two_body_reaction ( proj , targ , hadron , mu );
531 
532  G4LorentzVector hadron_in_trag_rest = hadron->Get4Momentum();
533  G4LorentzVector hadron_in_LAB = hadron_in_trag_rest.boost ( -boostToTargetRest );
534  aHadron.SetMomentum( hadron_in_LAB.v() );
535  aHadron.SetKineticEnergy ( hadron_in_LAB.e() - hadron_in_LAB.m() );
536 
537  delete proj;
538  delete targ;
539  delete hadron;
540 
541 //TKDB 100405
542 /*
543  G4double totalPhotonEnergy = 0;
544  if(thePhotons!=0)
545  {
546  unsigned int nPhotons = thePhotons->size();
547  unsigned int i0;
548  for(i0=0; i0<nPhotons; i0++)
549  {
550  totalPhotonEnergy += thePhotons->operator[](i0)->GetTotalEnergy();
551  }
552  }
553  availableEnergy -= totalPhotonEnergy;
554  residualMass += totalPhotonEnergy/G4Neutron::Neutron()->GetPDGMass();
555  aHadron.SetKineticEnergy(availableEnergy*residualMass*G4Neutron::Neutron()->GetPDGMass()/
556  (aHadron.GetMass()+residualMass*G4Neutron::Neutron()->GetPDGMass()));
557  G4double CosTheta = 1.0 - 2.0*G4UniformRand();
558  G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
559  G4double Phi = twopi*G4UniformRand();
560  G4ThreeVector Vector(std::cos(Phi)*SinTheta, std::sin(Phi)*SinTheta, CosTheta);
561  //aHadron.SetMomentum(Vector* std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
562  // aHadron.GetMass()*aHadron.GetMass()));
563  G4double p2 = aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()- aHadron.GetMass()*aHadron.GetMass();
564 
565  G4double p = 0.0;
566  if ( p2 > 0.0 )
567  p = std::sqrt ( p2 );
568 
569  aHadron.SetMomentum( Vector*p );
570 */
571 
572  }
573 
574 // fill the result
575 // Beware - the recoil is not necessarily in the particles...
576 // Can be calculated from momentum conservation?
577 // The idea is that the particles ar emitted forst, and the gammas only once the
578 // recoil is on the residual; assumption is that gammas do not contribute to
579 // the recoil.
580 // This needs more design @@@
581 
582  G4int nSecondaries = 2; // the hadron and the recoil
583  G4bool needsSeparateRecoil = false;
584  G4int totalBaryonNumber = 0;
585  G4int totalCharge = 0;
586  G4ThreeVector totalMomentum(0);
587  if(theParticles != 0)
588  {
589  nSecondaries = theParticles->size();
590  G4ParticleDefinition * aDef;
591  unsigned int ii0;
592  for(ii0=0; ii0<theParticles->size(); ii0++)
593  {
594  aDef = theParticles->operator[](ii0)->GetDefinition();
595  totalBaryonNumber+=aDef->GetBaryonNumber();
596  totalCharge+=G4int(aDef->GetPDGCharge()+eps);
597  totalMomentum += theParticles->operator[](ii0)->GetMomentum();
598  }
599  if(totalBaryonNumber!=G4int(theBaseA+eps+incidentParticle->GetDefinition()->GetBaryonNumber()))
600  {
601  needsSeparateRecoil = true;
602  nSecondaries++;
603  residualA = G4int(theBaseA+eps+incidentParticle->GetDefinition()->GetBaryonNumber()
604  -totalBaryonNumber);
605  residualZ = G4int(theBaseZ+eps+incidentParticle->GetDefinition()->GetPDGCharge()
606  -totalCharge);
607  }
608  }
609 
610  G4int nPhotons = 0;
611  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
612  nSecondaries += nPhotons;
613 
614  G4DynamicParticle * theSec;
615 
616  if( theParticles==0 )
617  {
618  theSec = new G4DynamicParticle;
619  theSec->SetDefinition(aHadron.GetDefinition());
620  theSec->SetMomentum(aHadron.GetMomentum());
621  theResult.AddSecondary(theSec);
622 
623  aHadron.Lorentz(aHadron, theTarget);
624  G4ReactionProduct theResidual;
626  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0));
627  theResidual.SetKineticEnergy(aHadron.GetKineticEnergy()*aHadron.GetMass()/theResidual.GetMass());
628 
629  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #6
630  //theResidual.SetMomentum(-1.*aHadron.GetMomentum());
631  G4ThreeVector incidentNeutronMomentum = theNeutron.GetMomentum();
632  theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
633 
634  theResidual.Lorentz(theResidual, -1.*theTarget);
635  G4ThreeVector totalPhotonMomentum(0,0,0);
636  if(thePhotons!=0)
637  {
638  for(i=0; i<nPhotons; i++)
639  {
640  totalPhotonMomentum += thePhotons->operator[](i)->GetMomentum();
641  }
642  }
643  theSec = new G4DynamicParticle;
644  theSec->SetDefinition(theResidual.GetDefinition());
645  theSec->SetMomentum(theResidual.GetMomentum()-totalPhotonMomentum);
646  theResult.AddSecondary(theSec);
647  }
648  else
649  {
650  for(i0=0; i0<theParticles->size(); i0++)
651  {
652  theSec = new G4DynamicParticle;
653  theSec->SetDefinition(theParticles->operator[](i0)->GetDefinition());
654  theSec->SetMomentum(theParticles->operator[](i0)->GetMomentum());
655  theResult.AddSecondary(theSec);
656  delete theParticles->operator[](i0);
657  }
658  delete theParticles;
659  if(needsSeparateRecoil && residualZ!=0)
660  {
661  G4ReactionProduct theResidual;
663  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0));
664  G4double resiualKineticEnergy = theResidual.GetMass()*theResidual.GetMass();
665  resiualKineticEnergy += totalMomentum*totalMomentum;
666  resiualKineticEnergy = std::sqrt(resiualKineticEnergy) - theResidual.GetMass();
667 // cout << "Kinetic energy of the residual = "<<resiualKineticEnergy<<endl;
668  theResidual.SetKineticEnergy(resiualKineticEnergy);
669 
670  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #4
671  //theResidual.SetMomentum(-1.*totalMomentum);
672  //G4ThreeVector incidentNeutronMomentum = theNeutron.GetMomentum();
673  //theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
674 //080717 TK Comment still do NOT include photon's mometum which produce by thePhotons
675  theResidual.SetMomentum( theNeutron.GetMomentum() + theTarget.GetMomentum() - totalMomentum );
676 
677  theSec = new G4DynamicParticle;
678  theSec->SetDefinition(theResidual.GetDefinition());
679  theSec->SetMomentum(theResidual.GetMomentum());
680  theResult.AddSecondary(theSec);
681  }
682  }
683  if(thePhotons!=0)
684  {
685  for(i=0; i<nPhotons; i++)
686  {
687  theSec = new G4DynamicParticle;
688  //Bug reported Chao Zhang (Chao.Zhang@usd.edu), Dongming Mei(Dongming.Mei@usd.edu) Feb. 25, 2009
689  //theSec->SetDefinition(G4Gamma::Gamma());
690  theSec->SetDefinition( thePhotons->operator[](i)->GetDefinition() );
691  //But never cause real effect at least with G4NDL3.13 TK
692  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
693  theResult.AddSecondary(theSec);
694  delete thePhotons->operator[](i);
695  }
696 // some garbage collection
697  delete thePhotons;
698  }
699 
700 //080721
702  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
703  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
704  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
705  adjust_final_state ( init_4p_lab );
706 
707 // clean up the primary neutron
709 }
710 
711 
712 
713 #include "G4RotationMatrix.hh"
714 void G4NeutronHPInelasticCompFS::two_body_reaction ( G4DynamicParticle* proj, G4DynamicParticle* targ, G4DynamicParticle* hadron, G4double mu )
715 {
716 
717 // Target rest flame
718 // 4vector in targ rest frame;
719 // targ could have excitation energy (photon energy will be emiited) tricky but,,,
720 
721  G4LorentzVector before = proj->Get4Momentum() + targ->Get4Momentum();
722 
723  G4ThreeVector p3_proj = proj->GetMomentum();
724  G4ThreeVector d = p3_proj.unit();
725  G4RotationMatrix rot;
726  G4RotationMatrix rot1;
727  rot1.setPhi( pi/2 + d.phi() );
728  G4RotationMatrix rot2;
729  rot2.setTheta( d.theta() );
730  rot=rot2*rot1;
731  proj->SetMomentum( rot*p3_proj );
732 
733 // Now proj only has pz component;
734 
735 // mu in CM system
736 
737  //Valid only for neutron incidence
739 
740  G4double Q = proj->GetDefinition()->GetPDGMass() + targ->GetDefinition()->GetPDGMass()
741  - ( hadron->GetDefinition()->GetPDGMass() + residual->GetDefinition()->GetPDGMass() );
742 
743  // Non Relativistic Case
744  G4double A = targ->GetDefinition()->GetPDGMass() / proj->GetDefinition()->GetPDGMass();
745  G4double AA = hadron->GetDefinition()->GetPDGMass() / proj->GetDefinition()->GetPDGMass();
746  G4double E1 = proj->GetKineticEnergy();
747 
748 // 101111
749 // In _nat_ data (Q+E1) could become negative value, following line is safty for this case.
750  //if ( (Q+E1) < 0 )
751  if ( ( 1 + (1+A)/A*Q/E1 ) < 0 )
752  {
753 // 1.0e-6 eV is additional safty for numeric precision
754  Q = -( A/(1+A)*E1 ) + 1.0e-6*eV;
755  }
756 
757  G4double beta = std::sqrt ( A*(A+1-AA)/AA*( 1 + (1+A)/A*Q/E1 ) );
758  G4double gamma = AA/(A+1-AA)*beta;
759  G4double E3 = AA/std::pow((1+A),2)*(beta*beta+1+2*beta*mu)*E1;
760  G4double omega3 = (1+beta*mu)/std::sqrt(beta*beta+1+2*beta*mu);
761  if ( omega3 > 1.0 ) omega3 = 1.0;
762 
763  G4double E4 = (A+1-AA)/std::pow((1+A),2)*(gamma*gamma+1-2*gamma*mu)*E1;
764  G4double omega4 = (1-gamma*mu)/std::sqrt(gamma*gamma+1-2*gamma*mu);
765  if ( omega4 > 1.0 ) omega4 = 1.0;
766 
767  hadron->SetKineticEnergy ( E3 );
768 
769  G4double M = hadron->GetDefinition()->GetPDGMass();
770  G4double pmag = std::sqrt ((E3+M)*(E3+M)-M*M) ;
771  G4ThreeVector p ( 0 , pmag*std::sqrt(1-omega3*omega3), pmag*omega3 );
772 
773  G4double M4 = residual->GetDefinition()->GetPDGMass();
774  G4double pmag4 = std::sqrt ((E4+M4)*(E4+M4)-M4*M4) ;
775  G4ThreeVector p4 ( 0 , -pmag4*std::sqrt(1-omega4*omega4), pmag4*omega4 );
776 
777 // Rotate to orginal target rest flame.
778  p *= rot.inverse();
779  hadron->SetMomentum( p );
780 // Now hadron had 4 momentum in target rest flame
781 
782 // TypeA
783  p4 *= rot.inverse();
784  residual->SetMomentum ( p4 );
785 
786 //TypeB1
787  //residual->Set4Momentum ( p4_residual );
788 //TypeB2
789  //residual->SetMomentum ( p4_residual.v() );
790 
791 // Type A make difference in Momenutum
792 // Type B1 make difference in Mass of residual
793 // Type B2 make difference in total energy.
794 
795  delete residual;
796 
797 }