Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ParticleHPInelasticCompFS.cc
이 파일의 문서화 페이지로 가기
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 
59 #include "G4NRESP71M03.hh" // nresp71_m03.hh and nresp71_m02.hh are alike. The only difference between m02 and m03 is in the total carbon cross section that is properly included in the latter. These data are not used in nresp71_m0*.hh.
60 
62 {
63  // char the[100] = {""};
64  // std::ostrstream ost(the, 100, std::ios::out);
65  // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
66  // G4String * aName = new G4String(the);
67  // std::ifstream from(*aName, std::ios::in);
68 
69  std::ostringstream ost;
70  ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
71  G4String aName = ost.str();
72  std::ifstream from(aName, std::ios::in);
73 
74  if(!from) return; // no data found for this isotope
75  // std::ifstream theGammaData(*aName, std::ios::in);
76  std::ifstream theGammaData(aName, std::ios::in);
77 
78  theGammas.Init(theGammaData);
79  // delete aName;
80 
81 }
82 
84 {
85  gammaPath = "/Inelastic/Gammas/"; //only in neutron data base
86  if(!getenv("G4NEUTRONHPDATA"))
87  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files where Inelastic/Gammas data is found.");
88  G4String tBase = getenv("G4NEUTRONHPDATA");
89  gammaPath = tBase+gammaPath;
90  G4String tString = dirName;
91  G4bool dbool;
92  G4ParticleHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, tString, aFSType, dbool);
93  G4String filename = aFile.GetName();
94 #ifdef G4PHPDEBUG
95  if( getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPInelasticCompFS::Init FILE " << filename << G4endl;
96 #endif
97 
98  SetAZMs( A, Z, M, aFile );
99  //theBaseA = aFile.GetA();
100  //theBaseZ = aFile.GetZ();
101  //theNDLDataA = (int)aFile.GetA();
102  //theNDLDataZ = aFile.GetZ();
103  //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
104  if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) )
105  {
106 #ifdef G4PHPDEBUG
107  if(getenv("G4ParticleHPDebug_NamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl;
108 #endif
109  hasAnyData = false;
110  hasFSData = false;
111  hasXsec = false;
112  return;
113  }
114  // theBaseA = A;
115  // theBaseZ = G4int(Z+.5);
116 //std::ifstream theData(filename, std::ios::in);
117  std::istringstream theData(std::ios::in);
118  G4ParticleHPManager::GetInstance()->GetDataStream(filename,theData);
119  if(!theData) //"!" is a operator of ios
120  {
121  hasAnyData = false;
122  hasFSData = false;
123  hasXsec = false;
124  // theData.close();
125  return;
126  }
127  // here we go
128  G4int infoType, dataType, dummy;
129  G4int sfType, it;
130  hasFSData = false;
131  while (theData >> infoType) // Loop checking, 11.05.2015, T. Koi
132  {
133  hasFSData = true;
134  theData >> dataType;
135  theData >> sfType >> dummy;
136  it = 50;
137  if(sfType>=600||(sfType<100&&sfType>=50)) it = sfType%50;
138  if(dataType==3)
139  {
140  //theData >> dummy >> dummy;
141  //TK110430
142  // QI and LR introudced since G4NDL3.15
143  G4double dqi;
144  G4int ilr;
145  theData >> dqi >> ilr;
146 
147  QI[ it ] = dqi*CLHEP::eV;
148  LR[ it ] = ilr;
150  G4int total;
151  theData >> total;
152  theXsection[it]->Init(theData, total, CLHEP::eV);
153  //std::cout << theXsection[it]->GetXsec(1*MeV) << std::endl;
154  }
155  else if(dataType==4)
156  {
158  theAngularDistribution[it]->Init(theData);
159  }
160  else if(dataType==5)
161  {
163  theEnergyDistribution[it]->Init(theData);
164  }
165  else if(dataType==6)
166  {
168  // G4cout << this << " CompFS theEnergyAngData " << it << theEnergyAngData[it] << G4endl; //GDEB
169  theEnergyAngData[it]->Init(theData);
170  }
171  else if(dataType==12)
172  {
174  theFinalStatePhotons[it]->InitMean(theData);
175  }
176  else if(dataType==13)
177  {
179  theFinalStatePhotons[it]->InitPartials(theData);
180  }
181  else if(dataType==14)
182  {
183  theFinalStatePhotons[it]->InitAngular(theData);
184  }
185  else if(dataType==15)
186  {
187  theFinalStatePhotons[it]->InitEnergies(theData);
188  }
189  else
190  {
191  throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4ParticleHPInelasticCompFS");
192  }
193  }
194  // theData.close();
195 }
196 
198 {
199 
200 // it = 0 has without Photon
201  G4double running[50];
202  running[0] = 0;
203  unsigned int i;
204  for(i=0; i<50; i++)
205  {
206  if(i!=0) running[i]=running[i-1];
207  if(theXsection[i] != 0)
208  {
209  running[i] += std::max(0., theXsection[i]->GetXsec(eKinetic));
210  }
211  }
212  G4double random = G4UniformRand();
213  G4double sum = running[49];
214  G4int it = 50;
215  if(0!=sum)
216  {
217  G4int i0;
218  for(i0=0; i0<50; i0++)
219  {
220  it = i0;
221  // G4cout << " SelectExitChannel " << it << " " << random << " " << running[i0]/sum << " " << running[i0] << G4endl; //GDEB
222  if(random < running[i0]/sum) break;
223  }
224  }
225 //debug: it = 1;
226 // G4cout << " SelectExitChannel " << it << " " << sum << G4endl; //GDEB
227  return it;
228 }
229 
230 
231  //n,p,d,t,he3,a
233 {
234 
235 // prepare neutron
236  if ( theResult.Get() == NULL ) theResult.Put( new G4HadFinalState );
237  theResult.Get()->Clear();
238  G4double eKinetic = theTrack.GetKineticEnergy();
239  const G4HadProjectile *hadProjectile = &theTrack;
240  G4ReactionProduct incidReactionProduct( const_cast<G4ParticleDefinition *>(hadProjectile->GetDefinition()) ); // incidReactionProduct
241  incidReactionProduct.SetMomentum( hadProjectile->Get4Momentum().vect() );
242  incidReactionProduct.SetKineticEnergy( eKinetic );
243 
244 // prepare target
245  G4int i;
246  for(i=0; i<50; i++)
247  { if(theXsection[i] != 0) { break; } }
248 
249  G4double targetMass=0;
250  G4double eps = 0.0001;
251  targetMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
252 #ifdef G4PHPDEBUG
253  if( getenv("G4ParticleHPDebug")) G4cout <<this <<" G4ParticleHPInelasticCompFS::CompositeApply A " <<theBaseA <<" Z " <<theBaseZ <<" incident " <<hadProjectile->GetDefinition()->GetParticleName() <<G4endl;
254 #endif
255 // if(theEnergyAngData[i]!=0)
256 // targetMass = theEnergyAngData[i]->GetTargetMass();
257 // else if(theAngularDistribution[i]!=0)
258 // targetMass = theAngularDistribution[i]->GetTargetMass();
259 // else if(theFinalStatePhotons[50]!=0)
260 // targetMass = theFinalStatePhotons[50]->GetTargetMass();
262  G4Nucleus aNucleus;
263  //G4ThreeVector neuVelo = (1./hadProjectile->GetDefinition()->GetPDGMass())*incidReactionProduct.GetMomentum();
264  //theTarget = aNucleus.GetBiasedThermalNucleus( targetMass/hadProjectile->GetDefinition()->GetPDGMass() , neuVelo, theTrack.GetMaterial()->GetTemperature());
265  //G4Nucleus::GetBiasedThermalNucleus requests normalization of mass and velocity in neutron mass
266  G4ThreeVector neuVelo = ( 1./G4Neutron::Neutron()->GetPDGMass() )*incidReactionProduct.GetMomentum();
267  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass/G4Neutron::Neutron()->GetPDGMass()
268  , neuVelo, theTrack.GetMaterial()->GetTemperature() );
269 
270  theTarget.SetDefinition( G4IonTable::GetIonTable()->GetIon( G4int(theBaseZ), G4int(theBaseA) , 0.0 ) ); //XX
271 
272 // prepare the residual mass
273  G4double residualMass=0;
274  G4double residualZ = theBaseZ + theProjectile->GetPDGCharge() - aDefinition->GetPDGCharge();
275  G4double residualA = theBaseA + theProjectile->GetBaryonNumber() - aDefinition->GetBaryonNumber();
276  residualMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(residualA+eps), static_cast<G4int>(residualZ+eps));
277 
278 // prepare energy in target rest frame
279  G4ReactionProduct boosted;
280  boosted.Lorentz(incidReactionProduct, theTarget);
281  eKinetic = boosted.GetKineticEnergy();
282 // G4double momentumInCMS = boosted.GetTotalMomentum();
283 
284 // select exit channel for composite FS class.
285  G4int it = SelectExitChannel( eKinetic );
286 
287 // set target and neutron in the relevant exit channel
288  InitDistributionInitialState(incidReactionProduct, theTarget, it);
289 
290  //---------------------------------------------------------------------//
291  //Hook for NRESP71MODEL
292  if ( G4ParticleHPManager::GetInstance()->GetUseNRESP71Model() ) {
293  if ( (G4int)(theBaseZ+0.1) == 6 ) // If the reaction is with Carbon...
294  {
296  if ( use_nresp71_model( aDefinition , it , theTarget , boosted ) ) return;
297  }
298  }
299  }
300  //---------------------------------------------------------------------//
301 
302  G4ReactionProductVector * thePhotons = 0;
303  G4ReactionProductVector * theParticles = 0;
304  G4ReactionProduct aHadron;
305  aHadron.SetDefinition(aDefinition); // what if only cross-sections exist ==> Na 23 11 @@@@
306  G4double availableEnergy = incidReactionProduct.GetKineticEnergy() + incidReactionProduct.GetMass() - aHadron.GetMass() +
307  (targetMass - residualMass);
308 //080730c
309  if ( availableEnergy < 0 )
310  {
311  //G4cout << "080730c Adjust availavleEnergy " << G4endl;
312  availableEnergy = 0;
313  }
314  G4int nothingWasKnownOnHadron = 0;
315  G4int dummy;
316  G4double eGamm = 0;
317  G4int iLevel=it-1;
318 
319 // TK without photon has it = 0
320  if( 50 == it )
321  {
322 
323 // TK Excitation level is not determined
324  iLevel=-1;
325  aHadron.SetKineticEnergy(availableEnergy*residualMass/
326  (aHadron.GetMass()+residualMass));
327 
328  //aHadron.SetMomentum(incidReactionProduct.GetMomentum()*(1./incidReactionProduct.GetTotalMomentum())*
329  // std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
330  // aHadron.GetMass()*aHadron.GetMass()));
331 
332  //TK add safty 100909
333  G4double p2 = ( aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy() - aHadron.GetMass()*aHadron.GetMass() );
334  G4double p = 0.0;
335  if ( p2 > 0.0 ) p = std::sqrt( p );
336 
337  aHadron.SetMomentum(incidReactionProduct.GetMomentum()*(1./incidReactionProduct.GetTotalMomentum())*p );
338 
339  }
340  else
341  {
342  while ( iLevel!=-1 && theGammas.GetLevel(iLevel) == 0 ) { iLevel--; } // Loop checking, 11.05.2015, T. Koi
343  }
344 
345 
346  if ( theAngularDistribution[it] != 0 ) // MF4
347  {
348  if(theEnergyDistribution[it]!=0) // MF5
349  {
350  //************************************************************
351  /*
352  aHadron.SetKineticEnergy(theEnergyDistribution[it]->Sample(eKinetic, dummy));
353  G4double eSecN = aHadron.GetKineticEnergy();
354  */
355  //************************************************************
356  //EMendoza --> maximum allowable energy should be taken into account.
357  G4double dqi = 0.0;
358  if ( QI[it] < 0 || 849 < QI[it] ) dqi = QI[it]; //For backword compatibility QI introduced since G4NDL3.15
359  G4double MaxEne=eKinetic+dqi;
360  G4double eSecN;
361 
362  G4int icounter=0;
363  G4int icounter_max=1024;
364  do {
365  icounter++;
366  if ( icounter > icounter_max ) {
367  G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
368  break;
369  }
370  eSecN=theEnergyDistribution[it]->Sample(eKinetic, dummy);
371  }while(eSecN>MaxEne); // Loop checking, 11.05.2015, T. Koi
372  aHadron.SetKineticEnergy(eSecN);
373  //************************************************************
374  eGamm = eKinetic-eSecN;
375  for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
376  {
377  if(theGammas.GetLevelEnergy(iLevel)<eGamm) break;
378  }
379  G4double random = 2*G4UniformRand();
380  iLevel+=G4int(random);
381  if(iLevel>theGammas.GetNumberOfLevels()-1)iLevel = theGammas.GetNumberOfLevels()-1;
382  }
383  else
384  {
385  G4double eExcitation = 0;
386  if(iLevel>=0) eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();
387  while (eKinetic-eExcitation < 0 && iLevel>0) // Loop checking, 11.05.2015, T. Koi
388  {
389  iLevel--;
390  eExcitation = theGammas.GetLevel(iLevel)->GetLevelEnergy();
391  }
392  //110610TK BEGIN
393  //Use QI value for calculating excitation energy of residual.
394  G4bool useQI=false;
395  G4double dqi = QI[it];
396  if ( dqi < 0 || 849 < dqi ) useQI = true; //Former libraies does not have values of this range
397 
398  if ( useQI )
399  {
400  // QI introudced since G4NDL3.15
401  eExcitation = -QI[it];
402  //Re-evluate iLevel based on this eExcitation
403  iLevel = 0;
404  G4bool find = false;
405  G4int imaxEx = 0;
406  while( theGammas.GetLevel(iLevel+1) != 0 ) // Loop checking, 11.05.2015, T. Koi
407  {
408  G4double maxEx = 0.0;
409  if ( maxEx < theGammas.GetLevel(iLevel)->GetLevelEnergy() )
410  {
411  maxEx = theGammas.GetLevel(iLevel)->GetLevelEnergy();
412  imaxEx = iLevel;
413  }
414  if ( eExcitation < theGammas.GetLevel(iLevel)->GetLevelEnergy() )
415  {
416  find = true;
417  iLevel--;
418  // very small eExcitation, iLevel becomes -1, this is protected below.
419  if ( iLevel == -1 ) iLevel = 0; // But cause energy trouble.
420  break;
421  }
422  iLevel++;
423  }
424  // In case, cannot find proper level, then use the maximum level.
425  if ( !find ) iLevel = imaxEx;
426  }
427  //110610TK END
428 
429  if(getenv("G4ParticleHPDebug") && eKinetic-eExcitation < 0)
430  {
431  throw G4HadronicException(__FILE__, __LINE__, "SEVERE: InelasticCompFS: Consistency of data not good enough, please file report");
432  }
433  if(eKinetic-eExcitation < 0) eExcitation = 0;
434  if(iLevel!= -1) aHadron.SetKineticEnergy(eKinetic - eExcitation);
435 
436  }
438 
439  if( theFinalStatePhotons[it] == 0 )
440  {
441  //G4cout << "110610 USE Gamma Level" << G4endl;
442 // TK comment Most n,n* eneter to this
443  thePhotons = theGammas.GetDecayGammas(iLevel);
444  eGamm -= theGammas.GetLevelEnergy(iLevel);
445  if(eGamm>0) // @ ok for now, but really needs an efficient way of correllated sampling @
446  {
447  G4ReactionProduct * theRestEnergy = new G4ReactionProduct;
448  theRestEnergy->SetDefinition(G4Gamma::Gamma());
449  theRestEnergy->SetKineticEnergy(eGamm);
450  G4double costh = 2.*G4UniformRand()-1.;
452  theRestEnergy->SetMomentum(eGamm*std::sin(std::acos(costh))*std::cos(phi),
453  eGamm*std::sin(std::acos(costh))*std::sin(phi),
454  eGamm*costh);
455  if(thePhotons == 0) { thePhotons = new G4ReactionProductVector; }
456  thePhotons->push_back(theRestEnergy);
457  }
458  }
459  }
460  else if(theEnergyAngData[it] != 0) // MF6
461  {
462 
463  theParticles = theEnergyAngData[it]->Sample(eKinetic);
464 
465  //141017 Fix BEGIN
466  //Adjust A and Z in the case of miss much between selected data and target nucleus
467  if ( theParticles != NULL ) {
468  G4int sumA = 0;
469  G4int sumZ = 0;
470  G4int maxA = 0;
471  G4int jAtMaxA = 0;
472  for ( G4int j = 0 ; j != (G4int)theParticles->size() ; j++ ) {
473  if ( theParticles->at(j)->GetDefinition()->GetBaryonNumber() > maxA ) {
474  maxA = theParticles->at(j)->GetDefinition()->GetBaryonNumber();
475  jAtMaxA = j;
476  }
477  sumA += theParticles->at(j)->GetDefinition()->GetBaryonNumber();
478  sumZ += G4int( theParticles->at(j)->GetDefinition()->GetPDGCharge() + eps );
479  }
480  G4int dA = (G4int)theBaseA + hadProjectile->GetDefinition()->GetBaryonNumber() - sumA;
481  G4int dZ = (G4int)theBaseZ + G4int( hadProjectile->GetDefinition()->GetPDGCharge() + eps ) - sumZ;
482  if ( dA < 0 || dZ < 0 ) {
483  G4int newA = theParticles->at(jAtMaxA)->GetDefinition()->GetBaryonNumber() + dA ;
484  G4int newZ = G4int( theParticles->at(jAtMaxA)->GetDefinition()->GetPDGCharge() + eps ) + dZ;
485  G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon ( newZ , newA );
486  theParticles->at( jAtMaxA )->SetDefinition( pd );
487  }
488  }
489  //141017 Fix END
490 
491  }
492  else
493  {
494  // @@@ what to do, if we have photon data, but no info on the hadron itself
495  nothingWasKnownOnHadron = 1;
496  }
497 
498  //G4cout << "theFinalStatePhotons it " << it << G4endl;
499  //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
500  //G4cout << "theFinalStatePhotons it " << it << G4endl;
501  //G4cout << "theFinalStatePhotons[it] " << theFinalStatePhotons[it] << G4endl;
502  //G4cout << "thePhotons " << thePhotons << G4endl;
503 
504  if ( theFinalStatePhotons[it] != 0 )
505  {
506  // the photon distributions are in the Nucleus rest frame.
507  // TK residual rest frame
508  G4ReactionProduct boosted_tmp;
509  boosted_tmp.Lorentz(incidReactionProduct, theTarget);
510  G4double anEnergy = boosted_tmp.GetKineticEnergy();
511  thePhotons = theFinalStatePhotons[it]->GetPhotons(anEnergy);
512  G4double aBaseEnergy = theFinalStatePhotons[it]->GetLevelEnergy();
513  G4double testEnergy = 0;
514  if(thePhotons!=0 && thePhotons->size()!=0)
515  { aBaseEnergy-=thePhotons->operator[](0)->GetTotalEnergy(); }
516  if(theFinalStatePhotons[it]->NeedsCascade())
517  {
518  while(aBaseEnergy>0.01*CLHEP::keV) // Loop checking, 11.05.2015, T. Koi
519  {
520  // cascade down the levels
521  G4bool foundMatchingLevel = false;
522  G4int closest = 2;
523  G4double deltaEold = -1;
524  for(G4int j=1; j<it; j++)
525  {
526  if(theFinalStatePhotons[j]!=0)
527  {
528  testEnergy = theFinalStatePhotons[j]->GetLevelEnergy();
529  }
530  else
531  {
532  testEnergy = 0;
533  }
534  G4double deltaE = std::abs(testEnergy-aBaseEnergy);
535  if(deltaE<0.1*CLHEP::keV)
536  {
537  G4ReactionProductVector * theNext =
538  theFinalStatePhotons[j]->GetPhotons(anEnergy);
539  if ( thePhotons != NULL ) thePhotons->push_back(theNext->operator[](0));
540  aBaseEnergy = testEnergy-theNext->operator[](0)->GetTotalEnergy();
541  delete theNext;
542  foundMatchingLevel = true;
543  break; // ===>
544  }
545  if(theFinalStatePhotons[j]!=0 && ( deltaE<deltaEold||deltaEold<0.) )
546  {
547  closest = j;
548  deltaEold = deltaE;
549  }
550  } // <=== the break goes here.
551  if(!foundMatchingLevel)
552  {
553  G4ReactionProductVector * theNext =
554  theFinalStatePhotons[closest]->GetPhotons(anEnergy);
555  if ( thePhotons != NULL ) thePhotons->push_back(theNext->operator[](0));
556  aBaseEnergy = aBaseEnergy-theNext->operator[](0)->GetTotalEnergy();
557  delete theNext;
558  }
559  }
560  }
561  }
562  unsigned int i0;
563  if(thePhotons!=0)
564  {
565  for(i0=0; i0<thePhotons->size(); i0++)
566  {
567  // back to lab
568  thePhotons->operator[](i0)->Lorentz(*(thePhotons->operator[](i0)), -1.*theTarget);
569  }
570  }
571  //G4cout << "nothingWasKnownOnHadron " << nothingWasKnownOnHadron << G4endl;
572  if(nothingWasKnownOnHadron)
573  {
574 // TKDB 100405
575 // In this case, hadron should be isotropic in CM
576 // mu and p should be correlated
577 //
578  G4double totalPhotonEnergy = 0.0;
579  if ( thePhotons != 0 )
580  {
581  unsigned int nPhotons = thePhotons->size();
582  unsigned int ii0;
583  for ( ii0=0; ii0<nPhotons; ii0++)
584  {
585  //thePhotons has energies at LAB system
586  totalPhotonEnergy += thePhotons->operator[](ii0)->GetTotalEnergy();
587  }
588  }
589 
590  //isotropic distribution in CM
591  G4double mu = 1.0 - 2 * G4UniformRand();
592 
593  // need momentums in target rest frame;
594  G4LorentzVector target_in_LAB ( theTarget.GetMomentum() , theTarget.GetTotalEnergy() );
595  G4ThreeVector boostToTargetRest = -target_in_LAB.boostVector();
596  G4LorentzVector proj_in_LAB = hadProjectile->Get4Momentum();
597 
598  G4DynamicParticle* proj = new G4DynamicParticle( theProjectile , proj_in_LAB.boost( boostToTargetRest ) );
599  G4DynamicParticle* targ = new G4DynamicParticle( G4IonTable::GetIonTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , totalPhotonEnergy ) , G4ThreeVector(0) );
600  G4DynamicParticle* hadron = new G4DynamicParticle( aHadron.GetDefinition() , G4ThreeVector(0) ); // will be fill momentum
601 
602  two_body_reaction ( proj , targ , hadron , mu );
603 
604  G4LorentzVector hadron_in_trag_rest = hadron->Get4Momentum();
605  G4LorentzVector hadron_in_LAB = hadron_in_trag_rest.boost ( -boostToTargetRest );
606  aHadron.SetMomentum( hadron_in_LAB.v() );
607  aHadron.SetKineticEnergy ( hadron_in_LAB.e() - hadron_in_LAB.m() );
608 
609  delete proj;
610  delete targ;
611  delete hadron;
612 
613 //TKDB 100405
614 /*
615  G4double totalPhotonEnergy = 0;
616  if(thePhotons!=0)
617  {
618  unsigned int nPhotons = thePhotons->size();
619  unsigned int i0;
620  for(i0=0; i0<nPhotons; i0++)
621  {
622  totalPhotonEnergy += thePhotons->operator[](i0)->GetTotalEnergy();
623  }
624  }
625  availableEnergy -= totalPhotonEnergy;
626  residualMass += totalPhotonEnergy/theProjectile->GetPDGMass();
627  aHadron.SetKineticEnergy(availableEnergy*residualMass*theProjectile->GetPDGMass()/
628  (aHadron.GetMass()+residualMass*theProjectile->GetPDGMass()));
629  G4double CosTheta = 1.0 - 2.0*G4UniformRand();
630  G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
631  G4double Phi = twopi*G4UniformRand();
632  G4ThreeVector Vector(std::cos(Phi)*SinTheta, std::sin(Phi)*SinTheta, CosTheta);
633  //aHadron.SetMomentum(Vector* std::sqrt(aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()-
634  // aHadron.GetMass()*aHadron.GetMass()));
635  G4double p2 = aHadron.GetTotalEnergy()*aHadron.GetTotalEnergy()- aHadron.GetMass()*aHadron.GetMass();
636 
637  G4double p = 0.0;
638  if ( p2 > 0.0 )
639  p = std::sqrt ( p2 );
640 
641  aHadron.SetMomentum( Vector*p );
642 */
643 
644  }
645 
646 // fill the result
647 // Beware - the recoil is not necessarily in the particles...
648 // Can be calculated from momentum conservation?
649 // The idea is that the particles ar emitted forst, and the gammas only once the
650 // recoil is on the residual; assumption is that gammas do not contribute to
651 // the recoil.
652 // This needs more design @@@
653 
654  G4int nSecondaries = 2; // the hadron and the recoil
655  G4bool needsSeparateRecoil = false;
656  G4int totalBaryonNumber = 0;
657  G4int totalCharge = 0;
658  G4ThreeVector totalMomentum(0);
659  if(theParticles != 0)
660  {
661  nSecondaries = theParticles->size();
662  const G4ParticleDefinition * aDef;
663  unsigned int ii0;
664  for(ii0=0; ii0<theParticles->size(); ii0++)
665  {
666  aDef = theParticles->operator[](ii0)->GetDefinition();
667  totalBaryonNumber+=aDef->GetBaryonNumber();
668  totalCharge+=G4int(aDef->GetPDGCharge()+eps);
669  totalMomentum += theParticles->operator[](ii0)->GetMomentum();
670  }
671  if(totalBaryonNumber!=G4int(theBaseA+eps+hadProjectile->GetDefinition()->GetBaryonNumber()))
672  {
673  needsSeparateRecoil = true;
674  nSecondaries++;
675  residualA = G4int(theBaseA+eps+hadProjectile->GetDefinition()->GetBaryonNumber()
676  -totalBaryonNumber);
677  residualZ = G4int(theBaseZ+eps+hadProjectile->GetDefinition()->GetPDGCharge()
678  -totalCharge);
679  }
680  }
681 
682  G4int nPhotons = 0;
683  if(thePhotons!=0) { nPhotons = thePhotons->size(); }
684  nSecondaries += nPhotons;
685 
686  G4DynamicParticle * theSec;
687 
688  if( theParticles==0 )
689  {
690  theSec = new G4DynamicParticle;
691  theSec->SetDefinition(aHadron.GetDefinition());
692  theSec->SetMomentum(aHadron.GetMomentum());
693  theResult.Get()->AddSecondary(theSec);
694 #ifdef G4PHPDEBUG
695  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary1 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
696 #endif
697 
698  aHadron.Lorentz(aHadron, theTarget);
699  G4ReactionProduct theResidual;
701  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0));
702  theResidual.SetKineticEnergy(aHadron.GetKineticEnergy()*aHadron.GetMass()/theResidual.GetMass());
703 
704  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #6
705  //theResidual.SetMomentum(-1.*aHadron.GetMomentum());
706  G4ThreeVector incidentNeutronMomentum = incidReactionProduct.GetMomentum();
707  theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
708 
709  theResidual.Lorentz(theResidual, -1.*theTarget);
710  G4ThreeVector totalPhotonMomentum(0,0,0);
711  if(thePhotons!=0)
712  {
713  for(i=0; i<nPhotons; i++)
714  {
715  totalPhotonMomentum += thePhotons->operator[](i)->GetMomentum();
716  }
717  }
718  theSec = new G4DynamicParticle;
719  theSec->SetDefinition(theResidual.GetDefinition());
720  theSec->SetMomentum(theResidual.GetMomentum()-totalPhotonMomentum);
721  theResult.Get()->AddSecondary(theSec);
722 #ifdef G4PHPDEBUG
723  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary2 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
724 #endif
725  }
726  else
727  {
728  for(i0=0; i0<theParticles->size(); i0++)
729  {
730  theSec = new G4DynamicParticle;
731  theSec->SetDefinition(theParticles->operator[](i0)->GetDefinition());
732  theSec->SetMomentum(theParticles->operator[](i0)->GetMomentum());
733  theResult.Get()->AddSecondary(theSec);
734 #ifdef G4PHPDEBUG
735  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary3 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
736 #endif
737  delete theParticles->operator[](i0);
738  }
739  delete theParticles;
740  if(needsSeparateRecoil && residualZ!=0)
741  {
742  G4ReactionProduct theResidual;
744  ->GetIon(static_cast<G4int>(residualZ), static_cast<G4int>(residualA), 0));
745  G4double resiualKineticEnergy = theResidual.GetMass()*theResidual.GetMass();
746  resiualKineticEnergy += totalMomentum*totalMomentum;
747  resiualKineticEnergy = std::sqrt(resiualKineticEnergy) - theResidual.GetMass();
748 // cout << "Kinetic energy of the residual = "<<resiualKineticEnergy<<endl;
749  theResidual.SetKineticEnergy(resiualKineticEnergy);
750 
751  //080612TK contribution from Benoit Pirard and Laurent Desorgher (Univ. Bern) #4
752  //theResidual.SetMomentum(-1.*totalMomentum);
753  //G4ThreeVector incidentNeutronMomentum = incidReactionProduct.GetMomentum();
754  //theResidual.SetMomentum(incidentNeutronMomentum - aHadron.GetMomentum());
755 //080717 TK Comment still do NOT include photon's mometum which produce by thePhotons
756  theResidual.SetMomentum( incidReactionProduct.GetMomentum() + theTarget.GetMomentum() - totalMomentum );
757 
758  theSec = new G4DynamicParticle;
759  theSec->SetDefinition(theResidual.GetDefinition());
760  theSec->SetMomentum(theResidual.GetMomentum());
761  theResult.Get()->AddSecondary(theSec);
762 #ifdef G4PHPDEBUG
763  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary4 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
764 #endif
765 
766  }
767  }
768  if(thePhotons!=0)
769  {
770  for(i=0; i<nPhotons; i++)
771  {
772  theSec = new G4DynamicParticle;
773  //Bug reported Chao Zhang (Chao.Zhang@usd.edu), Dongming Mei(Dongming.Mei@usd.edu) Feb. 25, 2009
774  //theSec->SetDefinition(G4Gamma::Gamma());
775  theSec->SetDefinition( thePhotons->operator[](i)->GetDefinition() );
776  //But never cause real effect at least with G4NDL3.13 TK
777  theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
778  theResult.Get()->AddSecondary(theSec);
779 #ifdef G4PHPDEBUG
780  if( getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticCompFS::BaseApply add secondary5 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
781 #endif
782 
783  delete thePhotons->operator[](i);
784  }
785 // some garbage collection
786  delete thePhotons;
787  }
788 
789 //080721
791  G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
792  G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
793  G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
794  adjust_final_state ( init_4p_lab );
795 
796 // clean up the primary neutron
798 }
799 
800 
801 
802 #include "G4RotationMatrix.hh"
804 {
805 
806 // Target rest flame
807 // 4vector in targ rest frame;
808 // targ could have excitation energy (photon energy will be emiited) tricky but,,,
809 
810  G4LorentzVector before = proj->Get4Momentum() + targ->Get4Momentum();
811 
812  G4ThreeVector p3_proj = proj->GetMomentum();
813  G4ThreeVector d = p3_proj.unit();
814  G4RotationMatrix rot;
815  G4RotationMatrix rot1;
816  rot1.setPhi( CLHEP::pi/2 + d.phi() );
817  G4RotationMatrix rot2;
818  rot2.setTheta( d.theta() );
819  rot=rot2*rot1;
820  proj->SetMomentum( rot*p3_proj );
821 
822 // Now proj only has pz component;
823 
824 // mu in CM system
825 
826  //Valid only for neutron incidence
827  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) );
828 
829  G4double Q = proj->GetDefinition()->GetPDGMass() + targ->GetDefinition()->GetPDGMass()
830  - ( hadron->GetDefinition()->GetPDGMass() + residual->GetDefinition()->GetPDGMass() );
831 
832  // Non Relativistic Case
833  G4double A = targ->GetDefinition()->GetPDGMass() / proj->GetDefinition()->GetPDGMass();
834  G4double AA = hadron->GetDefinition()->GetPDGMass() / proj->GetDefinition()->GetPDGMass();
835  G4double E1 = proj->GetKineticEnergy();
836 
837 // 101111
838 // In _nat_ data (Q+E1) could become negative value, following line is safty for this case.
839  //if ( (Q+E1) < 0 )
840  if ( ( 1 + (1+A)/A*Q/E1 ) < 0 )
841  {
842 // 1.0e-6 eV is additional safty for numeric precision
843  Q = -( A/(1+A)*E1 ) + 1.0e-6*CLHEP::eV;
844  }
845 
846  G4double beta = std::sqrt ( A*(A+1-AA)/AA*( 1 + (1+A)/A*Q/E1 ) );
847  G4double gamma = AA/(A+1-AA)*beta;
848  G4double E3 = AA/G4Pow::GetInstance()->powN((1+A),2)*(beta*beta+1+2*beta*mu)*E1;
849  G4double omega3 = (1+beta*mu)/std::sqrt(beta*beta+1+2*beta*mu);
850  if ( omega3 > 1.0 ) omega3 = 1.0;
851 
852  G4double E4 = (A+1-AA)/G4Pow::GetInstance()->powN((1+A),2)*(gamma*gamma+1-2*gamma*mu)*E1;
853  G4double omega4 = (1-gamma*mu)/std::sqrt(gamma*gamma+1-2*gamma*mu);
854  if ( omega4 > 1.0 ) omega4 = 1.0;
855 
856  hadron->SetKineticEnergy ( E3 );
857 
858  G4double M = hadron->GetDefinition()->GetPDGMass();
859  G4double pmag = std::sqrt ((E3+M)*(E3+M)-M*M) ;
860  G4ThreeVector p ( 0 , pmag*std::sqrt(1-omega3*omega3), pmag*omega3 );
861 
862  G4double M4 = residual->GetDefinition()->GetPDGMass();
863  G4double pmag4 = std::sqrt ((E4+M4)*(E4+M4)-M4*M4) ;
864  G4ThreeVector p4 ( 0 , -pmag4*std::sqrt(1-omega4*omega4), pmag4*omega4 );
865 
866 // Rotate to orginal target rest flame.
867  p *= rot.inverse();
868  hadron->SetMomentum( p );
869 // Now hadron had 4 momentum in target rest flame
870 
871 // TypeA
872  p4 *= rot.inverse();
873  residual->SetMomentum ( p4 );
874 
875 //TypeB1
876  //residual->Set4Momentum ( p4_residual );
877 //TypeB2
878  //residual->SetMomentum ( p4_residual.v() );
879 
880 // Type A make difference in Momenutum
881 // Type B1 make difference in Mass of residual
882 // Type B2 make difference in total energy.
883 
884  delete residual;
885 
886 }
887 
888 
889 
891 {
892  if ( aDefinition == G4Neutron::Definition() ) // If the outgoing particle is a neutron...
893  {
894  // LR: flag LR in ENDF. It indicates whether there is breakup of the residual nucleus or not.
895  // it: exit channel (index of the carbon excited state)
896 
897  //if ( (G4int)(theBaseZ+0.1) == 6 ) G4cout << "LR[" << it << "] = " << LR[it] << G4endl;
898 
899  // Added by A. R. García (CIEMAT) to include the physics of C(N,N'3A) reactions from NRESP71.
900 
901  if ( LR[it] > 0 ) // If there is breakup of the residual nucleus LR(flag LR in ENDF)>0 (i.e. Z=6 MT=52-91 (it=MT-50)).
902  {
903  // Defining carbon as the target in the reference frame at rest.
904  G4ReactionProduct theCarbon(theTarget);
905 
906  theCarbon.SetMomentum(G4ThreeVector());
907  theCarbon.SetKineticEnergy(0.);
908 
909  // Creating four reaction products.
910  G4ReactionProduct theProds[4];
911 
912  // Applying C(N,N'3A) reaction mechanisms in the target rest frame.
913  if ( it == 41 )
914  {
915  // QI=QM=-7.275 MeV for C-0(N,N')C-C(3A) in ENDF/B-VII.1. This is not the value of the QI of the first step according to the model. So we don't take it. Instead, we set the one we have calculated: QI=(mn+m12C)-(ma+m9Be+Ex9Be)=-8.130 MeV.
916  nresp71_model.ApplyMechanismI_NBeA2A(boosted, theCarbon, theProds, -8.130/*QI[it]*/); // N+C --> A[0]+9BE* | 9BE* --> N[1]+8BE | 8BE --> 2*A[2,3].
917  //printf("- QI=%f\n", QI[it]);
918  }
919  else
920  {
921  nresp71_model.ApplyMechanismII_ACN2A(boosted, theCarbon, theProds, QI[it]); // N+C --> N'[0]+C* | C* --> A[1]+8BE | 8BE --> 2*A[2,3].
922  }
923 
924  //printf("it=%d qi=%f \n", it, QI[it]);
925 
926  // Returning to the reference frame where the target was in motion.
927  for ( G4int j=0; j<4; j++ )
928  {
929  theProds[j].Lorentz(theProds[j], -1.*theTarget);
930  theResult.Get()->AddSecondary(new G4DynamicParticle(theProds[j].GetDefinition(), theProds[j].GetMomentum()));
931  }
932 
933  /*G4double EN0 = theNeutron.GetKineticEnergy();
934  G4double EN1 = theProds[0].GetKineticEnergy();
935 
936  G4double EA1 = theProds[1].GetKineticEnergy();
937  G4double EA2 = theProds[2].GetKineticEnergy();
938  G4double EA3 = theProds[3].GetKineticEnergy();
939 
940  printf("Q=%f\n", EN1+EA1+EA2+EA3-EN0);*/
941 
942  // Killing the primary neutron.
944 
945  return true;
946  }
947  }
948  else if ( aDefinition == G4Alpha::Definition() ) // If the outgoing particle is an alpha, ...
949  {
950  // Added by A. R. García (CIEMAT) to include the physics of C(N,A)9BE reactions from NRESP71.
951 
952  if ( LR[it] == 0 ) // If Z=6, an alpha particle is emitted and there is no breakup of the residual nucleus LR(flag LR in ENDF)==0.
953  {
954  // Defining carbon as the target in the reference frame at rest.
955  G4ReactionProduct theCarbon(theTarget);
956 
957  theCarbon.SetMomentum(G4ThreeVector());
958  theCarbon.SetKineticEnergy(0.);
959 
960  // Creating four reaction products.
961  G4ReactionProduct theProds[2];
962 
963  // Applying C(N,A)9BE reaction mechanism.
964  nresp71_model.ApplyMechanismABE(boosted, theCarbon, theProds); // N+C --> A[0]+9BE[1].
965 
966  //G4DynamicParticle *theSec;
967  for ( G4int j=0; j<2; j++ )
968  {
969  // Returning to the system of reference where the target was in motion.
970  theProds[j].Lorentz(theProds[j], -1.*theTarget);
971  theResult.Get()->AddSecondary(new G4DynamicParticle(theProds[j].GetDefinition(), theProds[j].GetMomentum()));
972  }
973 
974  // Killing the primary neutron.
976 
977  return true;
978  }
979  else
980  {
981  G4Exception("G4ParticleHPInelasticCompFS::CompositeApply()", "G4ParticleInelasticCompFS.cc", FatalException, "Alpha production with LR!=0.");
982  }
983  }
984 
985  return false;
986 }
G4double total(Particle const *const p1, Particle const *const p2)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4ParticleHPManager * GetInstance()
CLHEP::Hep3Vector G4ThreeVector
G4double Sample(G4double anEnergy, G4int &it)
void Init(std::istream &aDataFile)
static G4Alpha * Definition()
Definition: G4Alpha.cc:49
static G4Neutron * Definition()
Definition: G4Neutron.cc:54
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
void SetKineticEnergy(G4double aEnergy)
void InitGammas(G4double AR, G4double ZR)
TString targ[ntarg]
Definition: Style.C:27
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4int ApplyMechanismABE(G4ReactionProduct &neut, G4ReactionProduct &carb, G4ReactionProduct *theProds)
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * theProjectile
G4ReactionProductVector * Sample(G4double anEnergy)
void two_body_reaction(G4DynamicParticle *, G4DynamicParticle *, G4DynamicParticle *, G4double mu)
const G4String & GetParticleName() const
G4double GetPDGCharge() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4ParticleDefinition * GetDefinition() const
virtual G4ParticleHPVector * GetXsec()
Double_t beta
G4double GetPDGMass() const
G4ParticleHPPhotonDist * theFinalStatePhotons[51]
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:177
double theta() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ReactionProductVector * GetDecayGammas(G4int aLevel)
G4bool use_nresp71_model(const G4ParticleDefinition *aDefinition, const G4int it, const G4ReactionProduct &theTarget, G4ReactionProduct &boosted)
Float_t Z
void GetDataStream(G4String, std::istringstream &iss)
G4double GetTotalEnergy() const
G4int SelectExitChannel(G4double eKinetic)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void SetKineticEnergy(const G4double en)
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
G4double GetKineticEnergy() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:57
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
G4ParticleDefinition * GetDefinition() const
const G4Material * GetMaterial() const
G4double GetLevelEnergy(G4int aLevel)
G4Cache< G4HadFinalState * > theResult
G4bool InitMean(std::istream &aDataFile)
#define G4UniformRand()
Definition: Randomize.hh:53
G4ErrorTarget * theTarget
Definition: errprop.cc:59
G4double GetTotalMomentum() const
G4LorentzVector Get4Momentum() const
double A(double temperature)
G4ParticleHPAngular * theAngularDistribution[51]
Float_t d
static const G4int maxA
G4ParticleHPEnAngCorrelation * theEnergyAngData[51]
static G4double GetNuclearMass(const G4double A, const G4double Z)
double mag2() const
value_type & Get() const
Definition: G4Cache.hh:314
static constexpr double eV
G4ThreeVector GetMomentum() const
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4ParticleHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
Hep3Vector v() const
void Put(const value_type &val) const
Definition: G4Cache.hh:318
const G4LorentzVector & Get4Momentum() const
G4int GetNumberOfSecondaries() const
Hep3Vector unit() const
G4ReactionProductVector * GetPhotons(G4double anEnergy)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static constexpr double twopi
Definition: SystemOfUnits.h:55
void Init(std::istream &aDataFile)
void adjust_final_state(G4LorentzVector)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
G4double GetMass() const
const G4ParticleDefinition * GetDefinition() const
void setPhi(double phi)
Definition: RotationE.cc:262
int G4int
Definition: G4Types.hh:78
void SampleAndUpdate(G4ReactionProduct &anIncidentParticle)
double phi() const
ifstream in
Definition: comparison.C:7
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void InitPartials(std::istream &aDataFile)
void InitDistributionInitialState(G4ReactionProduct &anIncidentPart, G4ReactionProduct &aTarget, G4int it)
G4int ApplyMechanismI_NBeA2A(G4ReactionProduct &neut, G4ReactionProduct &carb, G4ReactionProduct *theProds, const G4double QI)
void InitAngular(std::istream &aDataFile)
void setTheta(double theta)
Definition: RotationE.cc:266
void SetMomentum(const G4ThreeVector &momentum)
G4double GetKineticEnergy() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4GLOB_DLL std::ostream G4cout
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aSFType, G4ParticleDefinition *)
void Init(std::istream &aDataFile)
Hep3Vector vect() const
static double Q[]
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:504
G4ParticleHPLevel * GetLevel(G4int i)
const G4ParticleDefinition * GetParticleDefinition() const
G4int ApplyMechanismII_ACN2A(G4ReactionProduct &neut, G4ReactionProduct &carb, G4ReactionProduct *theProds, const G4double QI)
static const G4double eps
static constexpr double keV
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:113
G4ParticleHPEnergyDistribution * theEnergyDistribution[51]
void CompositeApply(const G4HadProjectile &theTrack, G4ParticleDefinition *aHadron)
HepRotation inverse() const
static constexpr double pi
Definition: SystemOfUnits.h:54
HepLorentzVector & boost(double, double, double)
void SetStatusChange(G4HadFinalStateStatus aS)
G4double GetTemperature() const
Definition: G4Material.hh:183
void InitEnergies(std::istream &aDataFile)