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