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