Geant4_10
G4BinaryLightIonReaction.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 #include <algorithm>
27 #include <vector>
28 #include <cmath>
29 #include <numeric>
30 
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "G4LorentzVector.hh"
35 #include "G4LorentzRotation.hh"
37 #include "G4ping.hh"
38 #include "G4Delete.hh"
39 #include "G4Neutron.hh"
40 #include "G4VNuclearDensity.hh"
41 #include "G4FermiMomentum.hh"
42 #include "G4HadTmpUtil.hh"
43 #include "G4PreCompoundModel.hh"
45 
46 //#define debug_G4BinaryLightIonReaction
47 //#define debug_BLIR_finalstate
48 
50 : G4HadronicInteraction("Binary Light Ion Cascade"),
51  theProjectileFragmentation(ptr),
52  pA(0),pZ(0), tA(0),tZ(0),spectatorA(0),spectatorZ(0),
53  projectile3dNucleus(0),target3dNucleus(0)
54 {
55  if(!ptr) {
58  G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
59  if(!pre) { pre = new G4PreCompoundModel(); }
60  theProjectileFragmentation = pre;
61  }
62  theModel = new G4BinaryCascade(theProjectileFragmentation);
63  theHandler = theProjectileFragmentation->GetExcitationHandler();
64 
65  debug_G4BinaryLightIonReactionResults=getenv("debug_G4BinaryLightIonReactionResults")!=0;
66 }
67 
69 {}
70 
72 {
73  outFile << "G4Binary Light Ion Cascade is an intra-nuclear cascade model\n"
74  << "using G4BinaryCasacde to model the interaction of a light\n"
75  << "nucleus with a nucleus.\n"
76  << "The lighter of the two nuclei is treated like a set of projectiles\n"
77  << "which are transported simultanously through the heavier nucleus.\n";
78 }
79 
80 //--------------------------------------------------------------------------------
82 {
84 };
85 
87 ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus & targetNucleus )
88 {
89  if(getenv("BLICDEBUG") ) G4cerr << " ######### Binary Light Ion Reaction starts ######### " << G4endl;
90  G4ping debug("debug_G4BinaryLightIonReaction");
91  pA=aTrack.GetDefinition()->GetBaryonNumber();
92  pZ=G4lrint(aTrack.GetDefinition()->GetPDGCharge()/eplus);
93  tA=targetNucleus.GetA_asInt();
94  tZ=targetNucleus.GetZ_asInt();
95 
96  G4LorentzVector mom(aTrack.Get4Momentum());
97  //G4cout << "proj mom : " << mom << G4endl;
98  G4LorentzRotation toBreit(mom.boostVector());
99 
100  G4bool swapped=SetLighterAsProjectile(mom, toBreit);
101  //G4cout << "after swap, swapped? / mom " << swapped << " / " << mom <<G4endl;
103  G4ReactionProductVector * cascaders=0; //new G4ReactionProductVector;
104 // G4double m_nucl(0); // to check energy balance
105 
106 
107  // G4double m1=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(pZ,pA);
108  // G4cout << "Entering the decision point "
109  // << (mom.t()-mom.mag())/pA << " "
110  // << pA<<" "<< pZ<<" "
111  // << tA<<" "<< tZ<<G4endl
112  // << " "<<mom.t()-mom.mag()<<" "
113  // << mom.t()- m1<<G4endl;
114  if( (mom.t()-mom.mag())/pA < 50*MeV )
115  {
116  // G4cout << "Using pre-compound only, E= "<<mom.t()-mom.mag()<<G4endl;
117  // m_nucl = mom.mag();
118  cascaders=FuseNucleiAndPrompound(mom);
119  }
120  else
121  {
122  result=Interact(mom,toBreit);
123 
124  if(! result )
125  {
126  {
127  // abort!!
128 
129  G4cerr << "G4BinaryLightIonReaction no final state for: " << G4endl;
130  G4cerr << " Primary " << aTrack.GetDefinition()
131  << ", (A,Z)=(" << aTrack.GetDefinition()->GetBaryonNumber()
132  << "," << aTrack.GetDefinition()->GetPDGCharge()/eplus << ") "
133  << ", kinetic energy " << aTrack.GetKineticEnergy()
134  << G4endl;
135  G4cerr << " Target nucleus (A,Z)=("
136  << (swapped?pA:tA) << ","
137  << (swapped?pZ:tZ) << ")" << G4endl;
138  G4cerr << " if frequent, please submit above information as bug report"
139  << G4endl << G4endl;
140 
141  theResult.Clear();
142  theResult.SetStatusChange(isAlive);
143  theResult.SetEnergyChange(aTrack.GetKineticEnergy());
144  theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
145  return &theResult;
146 
147  }
148  }
149 
150  // Calculate excitation energy,
151  G4double theStatisticalExEnergy = GetProjectileExcitation();
152 
153 
154  pInitialState = mom;
155  //G4cout << "pInitialState from aTrack : " << pInitialState;
156  pInitialState.setT(pInitialState.getT() +
158  //G4cout << "target nucleus added : " << pInitialState << G4endl;
159 
160  delete target3dNucleus;target3dNucleus=0;
161  delete projectile3dNucleus;projectile3dNucleus=0;
162 
164 
165  cascaders = new G4ReactionProductVector;
166 
167  G4LorentzVector pspectators=SortResult(result,spectators,cascaders);
168 
169  // pFinalState=std::accumulate(cascaders->begin(),cascaders->end(),pFinalState,ReactionProduct4Mom);
170  std::vector<G4ReactionProduct *>::iterator iter;
171 
172  //G4cout << "pInitialState, pFinalState / pspectators"<< pInitialState << " / " << pFinalState << " / " << pspectators << G4endl;
173  // if ( spectA-spectatorA !=0 || spectZ-spectatorZ !=0)
174  // {
175  // G4cout << "spect Nucl != spectators: nucl a,z; spect a,z" <<
176  // spectatorA <<" "<< spectatorZ <<" ; " << spectA <<" "<< spectZ << G4endl;
177  // }
178  delete result;
179  result=0;
180  G4LorentzVector momentum(pInitialState-pFinalState);
181  G4int loopcount(0);
182  //G4cout << "momentum, pspectators : " << momentum << " / " << pspectators << G4endl;
183  while (std::abs(momentum.e()-pspectators.e()) > 10*MeV)
184  {
185  G4LorentzVector pCorrect(pInitialState-pspectators);
186  //G4cout << "BIC nonconservation? (pInitialState-pFinalState) / spectators :" << momentum << " / " << pspectators << "pCorrect "<< pCorrect<< G4endl;
187  // Correct outgoing casacde particles.... to have momentum of (initial state - spectators)
188  G4bool EnergyIsCorrect=EnergyAndMomentumCorrector(cascaders, pCorrect);
189  if ( ! EnergyIsCorrect && debug_G4BinaryLightIonReactionResults)
190  {
191  G4cout << "Warning - G4BinaryLightIonReaction E/P correction for cascaders failed" << G4endl;
192  }
193  pFinalState=G4LorentzVector(0,0,0,0);
194  unsigned int i;
195  for(i=0; i<cascaders->size(); i++)
196  {
197  pFinalState += G4LorentzVector( (*cascaders)[i]->GetMomentum(), (*cascaders)[i]->GetTotalEnergy() );
198  }
199  momentum=pInitialState-pFinalState;
200  if (++loopcount > 10 )
201  {
202  if ( momentum.vect().mag() - momentum.e()> 10*keV )
203  {
204  G4cerr << "G4BinaryLightIonReaction.cc: Cannot correct 4-momentum of cascade particles" << G4endl;
205  throw G4HadronicException(__FILE__, __LINE__, "G4BinaryCasacde::ApplyCollision()");
206  } else {
207  break;
208  }
209 
210  }
211  }
212 
213  if (spectatorA > 0 )
214  {
215  // check spectator momentum
216  if ( momentum.vect().mag() - momentum.e()> 10*keV )
217  {
218 
219  G4ReactionProductVector::iterator ispectator;
220  for (ispectator=spectators->begin();ispectator!=spectators->end();ispectator++)
221  {
222  delete *ispectator;
223  }
224  delete spectators;
225 
226  G4cout << "G4BinaryLightIonReaction.cc: mom check: " << momentum
227  << " 3.mag "<< momentum.vect().mag() << G4endl
228  << " .. pInitialState/pFinalState/spectators " << pInitialState <<" "
229  << pFinalState << " " << pspectators << G4endl
230  << " .. A,Z " << spectatorA <<" "<< spectatorZ << G4endl;
231  G4cout << "G4BinaryLightIonReaction invalid final state for: " << G4endl;
232  G4cout << " Primary " << aTrack.GetDefinition()
233  << ", (A,Z)=(" << aTrack.GetDefinition()->GetBaryonNumber()
234  << "," << aTrack.GetDefinition()->GetPDGCharge()/eplus << ") "
235  << ", kinetic energy " << aTrack.GetKineticEnergy()
236  << G4endl;
237  G4cout << " Target nucleus (A,Z)=(" << targetNucleus.GetA_asInt()
238  << "," << targetNucleus.GetZ_asInt() << ")" << G4endl;
239  G4cout << " if frequent, please submit above information as bug report"
240  << G4endl << G4endl;
241 
242  theResult.Clear();
243  theResult.SetStatusChange(isAlive);
244  theResult.SetEnergyChange(aTrack.GetKineticEnergy());
245  theResult.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
246  return &theResult;
247  }
248 
249 
250  DeExciteSpectatorNucleus(spectators, cascaders, theStatisticalExEnergy, momentum);
251  } else {
252  delete spectators;
253  }
254  }
255  // Rotate to lab
256  G4LorentzRotation toZ;
257  toZ.rotateZ(-1*mom.phi());
258  toZ.rotateY(-1*mom.theta());
259  G4LorentzRotation toLab(toZ.inverse());
260 
261  // Fill the particle change, while rotating. Boost from projectile breit-frame in case we swapped.
262  // theResult.Clear();
263  theResult.Clear();
264  theResult.SetStatusChange(stopAndKill);
265  G4double Etot(0);
266  size_t i=0;
267  for(i=0; i<cascaders->size(); i++)
268  {
269  if((*cascaders)[i]->GetNewlyAdded())
270  {
271  G4DynamicParticle * aNew =
272  new G4DynamicParticle((*cascaders)[i]->GetDefinition(),
273  (*cascaders)[i]->GetTotalEnergy(),
274  (*cascaders)[i]->GetMomentum() );
275  G4LorentzVector tmp = aNew->Get4Momentum();
276  if(swapped)
277  {
278  tmp*=toBreit.inverse();
279  tmp.setVect(-tmp.vect());
280  }
281  tmp *= toLab;
282  aNew->Set4Momentum(tmp);
283  //G4cout << "result[" << i << "], 4vect: " << tmp << G4endl;
284  theResult.AddSecondary(aNew);
285  Etot += tmp.e();
286  // G4cout << "LIBIC: Secondary " << aNew->GetDefinition()->GetParticleName()
287  // <<" "<< aNew->GetMomentum()
288  // <<" "<< aNew->GetTotalEnergy()
289  // << G4endl;
290  }
291  delete (*cascaders)[i];
292  }
293  delete cascaders;
294 
295 #ifdef debug_BLIR_result
296  G4cout << "Result analysis, secondaries" << theResult.GetNumberOfSecondaries() << G4endl;
297  G4cout << " Energy conservation initial/primary/nucleus/final/delta(init-final) "
298  << aTrack.GetTotalEnergy() + m_nucl << aTrack.GetTotalEnergy() << m_nucl <<Etot
299  << aTrack.GetTotalEnergy() + m_nucl - Etot;
300 #endif
301 
302  if(getenv("BLICDEBUG") ) G4cerr << " ######### Binary Light Ion Reaction number ends ######### " << G4endl;
303 
304  return &theResult;
305 }
306 
307 //--------------------------------------------------------------------------------
308 
309 //****************************************************************************
310 G4bool G4BinaryLightIonReaction::EnergyAndMomentumCorrector(
311  G4ReactionProductVector* Output, G4LorentzVector& TotalCollisionMom)
312 //****************************************************************************
313 {
314  const int nAttemptScale = 2500;
315  const double ErrLimit = 1.E-6;
316  if (Output->empty())
317  return TRUE;
318  G4LorentzVector SumMom(0,0,0,0);
319  G4double SumMass = 0;
320  G4double TotalCollisionMass = TotalCollisionMom.m();
321  size_t i = 0;
322  // Calculate sum hadron 4-momenta and summing hadron mass
323  for(i = 0; i < Output->size(); i++)
324  {
325  SumMom += G4LorentzVector((*Output)[i]->GetMomentum(),(*Output)[i]->GetTotalEnergy());
326  SumMass += (*Output)[i]->GetDefinition()->GetPDGMass();
327  }
328  // G4cout << " E/P corrector, SumMass, SumMom.m2, TotalMass "
329  // << SumMass <<" "<< SumMom.m2() <<" "<<TotalCollisionMass<< G4endl;
330  if (SumMass > TotalCollisionMass) return FALSE;
331  SumMass = SumMom.m2();
332  if (SumMass < 0) return FALSE;
333  SumMass = std::sqrt(SumMass);
334 
335  // Compute c.m.s. hadron velocity and boost KTV to hadron c.m.s.
336  G4ThreeVector Beta = -SumMom.boostVector();
337  //G4cout << " == pre boost 2 "<< SumMom.e()<< " "<< SumMom.mag()<<" "<< Beta <<G4endl;
338  //--old Output->Boost(Beta);
339  for(i = 0; i < Output->size(); i++)
340  {
341  G4LorentzVector mom = G4LorentzVector((*Output)[i]->GetMomentum(),(*Output)[i]->GetTotalEnergy());
342  mom *= Beta;
343  (*Output)[i]->SetMomentum(mom.vect());
344  (*Output)[i]->SetTotalEnergy(mom.e());
345  }
346 
347  // Scale total c.m.s. hadron energy (hadron system mass).
348  // It should be equal interaction mass
349  G4double Scale = 0,OldScale=0;
350  G4double factor = 1.;
351  G4int cAttempt = 0;
352  G4double Sum = 0;
353  G4bool success = false;
354  for(cAttempt = 0; cAttempt < nAttemptScale; cAttempt++)
355  {
356  Sum = 0;
357  for(i = 0; i < Output->size(); i++)
358  {
359  G4LorentzVector HadronMom = G4LorentzVector((*Output)[i]->GetMomentum(),(*Output)[i]->GetTotalEnergy());
360  HadronMom.setVect(HadronMom.vect()+ factor*Scale*HadronMom.vect());
361  G4double E = std::sqrt(HadronMom.vect().mag2() + sqr((*Output)[i]->GetDefinition()->GetPDGMass()));
362  HadronMom.setE(E);
363  (*Output)[i]->SetMomentum(HadronMom.vect());
364  (*Output)[i]->SetTotalEnergy(HadronMom.e());
365  Sum += E;
366  }
367  OldScale=Scale;
368  Scale = TotalCollisionMass/Sum - 1;
369  // G4cout << "E/P corr - " << cAttempt << " " << Scale << G4endl;
370  if (std::abs(Scale) <= ErrLimit
371  || OldScale == Scale) // protect 'frozen' situation and divide by 0 in calculating new factor below
372  {
373  if (debug_G4BinaryLightIonReactionResults) G4cout << "E/p corrector: " << cAttempt << G4endl;
374  success = true;
375  break;
376  }
377  if ( cAttempt > 10 )
378  {
379  // G4cout << " speed it up? " << std::abs(OldScale/(OldScale-Scale)) << G4endl;
380  factor=std::max(1.,std::log(std::abs(OldScale/(OldScale-Scale))));
381  // G4cout << " ? factor ? " << factor << G4endl;
382  }
383  }
384 
385  if( (!success) && debug_G4BinaryLightIonReactionResults)
386  {
387  G4cout << "G4G4BinaryLightIonReaction::EnergyAndMomentumCorrector - Warning"<<G4endl;
388  G4cout << " Scale not unity at end of iteration loop: "<<TotalCollisionMass<<" "<<Sum<<" "<<Scale<<G4endl;
389  G4cout << " Increase number of attempts or increase ERRLIMIT"<<G4endl;
390  }
391 
392  // Compute c.m.s. interaction velocity and KTV back boost
393  Beta = TotalCollisionMom.boostVector();
394  //--old Output->Boost(Beta);
395  for(i = 0; i < Output->size(); i++)
396  {
397  G4LorentzVector mom = G4LorentzVector((*Output)[i]->GetMomentum(),(*Output)[i]->GetTotalEnergy());
398  mom *= Beta;
399  (*Output)[i]->SetMomentum(mom.vect());
400  (*Output)[i]->SetTotalEnergy(mom.e());
401  }
402  return TRUE;
403 }
404 G4bool G4BinaryLightIonReaction::SetLighterAsProjectile(G4LorentzVector & mom,const G4LorentzRotation & toBreit)
405 {
406  G4bool swapped = false;
407  if(tA<pA)
408  {
409  swapped = true;
410  G4int tmp(0);
411  tmp = tA; tA=pA; pA=tmp;
412  tmp = tZ; tZ=pZ; pZ=tmp;
414  G4LorentzVector it(m1, G4ThreeVector(0,0,0));
415  mom = toBreit*it;
416  }
417  return swapped;
418 }
419 G4ReactionProductVector * G4BinaryLightIonReaction::FuseNucleiAndPrompound(const G4LorentzVector & mom)
420 {
421  G4Fragment aPreFrag;
422  aPreFrag.SetA(pA+tA);
423  aPreFrag.SetZ(pZ+tZ);
424  aPreFrag.SetNumberOfParticles(pA);
425  aPreFrag.SetNumberOfCharged(pZ);
426  aPreFrag.SetNumberOfHoles(0);
427  G4ThreeVector plop(0.,0., mom.vect().mag());
429  G4LorentzVector aL(mom.t()+m_nucl, plop);
430  aPreFrag.SetMomentum(aL);
431 
432 
433  // G4cout << "Fragment INFO "<< pA+tA <<" "<<pZ+tZ<<" "
434  // << aL <<" "<<preFragDef->GetParticleName()<<G4endl;
435  G4ReactionProductVector * cascaders = theProjectileFragmentation->DeExcite(aPreFrag);
436  //G4double tSum = 0;
437  for(size_t count = 0; count<cascaders->size(); count++)
438  {
439  cascaders->operator[](count)->SetNewlyAdded(true);
440  //tSum += cascaders->operator[](count)->GetKineticEnergy();
441  }
442  // G4cout << "Exiting pre-compound only, E= "<<tSum<<G4endl;
443  return cascaders;
444 }
445 G4ReactionProductVector * G4BinaryLightIonReaction::Interact(G4LorentzVector & mom, const G4LorentzRotation & toBreit)
446 {
448  G4double projectileMass(0);
449  G4LorentzVector it;
450 
451  G4int tryCount(0);
452  do
453  {
454  ++tryCount;
455  projectile3dNucleus = new G4Fancy3DNucleus;
456  projectile3dNucleus->Init(pA, pZ);
457  projectile3dNucleus->CenterNucleons();
459  projectile3dNucleus->GetCharge(),projectile3dNucleus->GetMassNumber());
460  it=toBreit * G4LorentzVector(projectileMass,G4ThreeVector(0,0,0));
461 
462  target3dNucleus = new G4Fancy3DNucleus;
463  target3dNucleus->Init(tA, tZ);
464  G4double impactMax = target3dNucleus->GetOuterRadius()+projectile3dNucleus->GetOuterRadius();
465  // G4cout << "out radius - nucleus - projectile " << target3dNucleus->GetOuterRadius()/fermi << " - " << projectile3dNucleus->GetOuterRadius()/fermi << G4endl;
466  G4double aX=(2.*G4UniformRand()-1.)*impactMax;
467  G4double aY=(2.*G4UniformRand()-1.)*impactMax;
468  G4ThreeVector pos(aX, aY, -2.*impactMax-5.*fermi);
469 
470  G4KineticTrackVector * initalState = new G4KineticTrackVector;
471  projectile3dNucleus->StartLoop();
472  G4Nucleon * aNuc;
473  G4LorentzVector tmpV(0,0,0,0);
474  G4LorentzVector nucleonMom(1./pA*mom);
475  nucleonMom.setZ(nucleonMom.vect().mag());
476  nucleonMom.setX(0);
477  nucleonMom.setY(0);
478  theFermi.Init(pA,pZ);
479  while( (aNuc=projectile3dNucleus->GetNextNucleon()) )
480  {
481  G4LorentzVector p4 = aNuc->GetMomentum();
482  tmpV+=p4;
483  G4ThreeVector nucleonPosition(aNuc->GetPosition());
484  G4double density=(projectile3dNucleus->GetNuclearDensity())->GetDensity(nucleonPosition);
485  nucleonPosition += pos;
486  G4KineticTrack * it1 = new G4KineticTrack(aNuc, nucleonPosition, nucleonMom );
488  G4double pfermi= theFermi.GetFermiMomentum(density);
489  G4double mass = aNuc->GetDefinition()->GetPDGMass();
490  G4double Efermi= std::sqrt( sqr(mass) + sqr(pfermi)) - mass;
491  it1->SetProjectilePotential(-Efermi);
492  initalState->push_back(it1);
493  }
494 
495  result=theModel->Propagate(initalState, target3dNucleus);
496  if( result && result->size()==0)
497  {
498  delete result;
499  result=0;
500  }
501  if ( ! result )
502  {
503  delete target3dNucleus;
504  delete projectile3dNucleus;
505  }
506 
507  // std::for_each(initalState->begin(), initalState->end(), Delete<G4KineticTrack>());
508  // delete initalState;
509 
510  } while (! result && tryCount< 150);
511  return result;
512 }
513 G4double G4BinaryLightIonReaction::GetProjectileExcitation()
514 {
515  spectatorA=spectatorZ=0;
516 
517  G4Nucleon * aNuc;
518  // targetNucleus->StartLoop();
519  // while( (aNuc=targetNucleus->GetNextNucleon()) )
520  // {
521  // G4cout << " tgt Nucleon : " << aNuc->GetDefinition()->GetParticleName() <<" "<< aNuc->AreYouHit() <<" "<<aNuc->GetMomentum()<<G4endl;
522  // }
523  // the projectileNucleus excitation energy estimate...
524  G4double theStatisticalExEnergy = 0;
525  projectile3dNucleus->StartLoop();
526  while( (aNuc=projectile3dNucleus->GetNextNucleon()) )
527  {
528  // G4cout << " Nucleon : " << aNuc->GetDefinition()->GetParticleName() <<" "<< aNuc->AreYouHit() <<" "<<aNuc->GetMomentum()<<G4endl;
529  if(!aNuc->AreYouHit())
530  {
531  spectatorA++;
532  spectatorZ+=G4lrint(aNuc->GetDefinition()->GetPDGCharge()/eplus);
533  }
534  else
535  {
536  G4ThreeVector aPosition(aNuc->GetPosition());
537  G4double localDensity = projectile3dNucleus->GetNuclearDensity()->GetDensity(aPosition);
538  G4double localPfermi = theFermi.GetFermiMomentum(localDensity);
539  G4double nucMass = aNuc->GetDefinition()->GetPDGMass();
540  G4double localFermiEnergy = std::sqrt(nucMass*nucMass + localPfermi*localPfermi) - nucMass;
541  G4double deltaE = localFermiEnergy - (aNuc->GetMomentum().t()-aNuc->GetMomentum().mag());
542  theStatisticalExEnergy += deltaE;
543  }
544  }
545  return theStatisticalExEnergy;
546 }
547 
548 G4LorentzVector G4BinaryLightIonReaction::SortResult(G4ReactionProductVector * result, G4ReactionProductVector * spectators,G4ReactionProductVector * cascaders)
549 {
550  unsigned int i(0);
551  // G4int spectA(0),spectZ(0);
552  G4LorentzVector pspectators(0,0,0,0);
553  pFinalState=G4LorentzVector(0,0,0,0);
554  for(i=0; i<result->size(); i++)
555  {
556  if( (*result)[i]->GetNewlyAdded() )
557  {
558  pFinalState += G4LorentzVector( (*result)[i]->GetMomentum(), (*result)[i]->GetTotalEnergy() );
559  cascaders->push_back((*result)[i]);
560  }
561  else {
562  // G4cout <<" spectator ... ";
563  pspectators += G4LorentzVector( (*result)[i]->GetMomentum(), (*result)[i]->GetTotalEnergy() );
564  spectators->push_back((*result)[i]);
565  // spectA++;
566  // spectZ+= G4lrint((*result)[i]->GetDefinition()->GetPDGCharge()/eplus);
567  }
568 
569  // G4cout << (*result)[i]<< " "
570  // << (*result)[i]->GetDefinition()->GetParticleName() << " "
571  // << (*result)[i]->GetMomentum()<< " "
572  // << (*result)[i]->GetTotalEnergy() << G4endl;
573  }
574  //G4cout << "pFinalState / pspectators" << pFinalState << " / " << pspectators << G4endl;
575  return pspectators;
576 }
577 
578 void G4BinaryLightIonReaction::DeExciteSpectatorNucleus(G4ReactionProductVector * spectators, G4ReactionProductVector * cascaders,
579  G4double theStatisticalExEnergy, G4LorentzVector & pSpectators)
580 {
581  // call precompound model
582  G4ReactionProductVector * proFrag = 0;
583  G4LorentzVector pFragment(0.,0.,0.,0.);
584  // G4cout << " == pre boost 1 "<< momentum.e()<< " "<< momentum.mag()<<G4endl;
585  G4LorentzRotation boost_fragments;
586  // G4cout << " == post boost 1 "<< momentum.e()<< " "<< momentum.mag()<<G4endl;
587  // G4LorentzRotation boost_spectator_mom(-momentum.boostVector());
588  // G4cout << "- momentum " << boost_spectator_mom * momentum << G4endl;
589  G4LorentzVector pFragments(0,0,0,0);
590 
591  if(spectatorZ>0 && spectatorA>1)
592  {
593  // Make the fragment
594  G4Fragment aProRes;
595  aProRes.SetA(spectatorA);
596  aProRes.SetZ(spectatorZ);
597  aProRes.SetNumberOfParticles(0);
598  aProRes.SetNumberOfCharged(0);
599  aProRes.SetNumberOfHoles(pA-spectatorA);
600  G4double mFragment=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(spectatorZ,spectatorA);
601  pFragment=G4LorentzVector(0,0,0,mFragment+std::max(0.,theStatisticalExEnergy) );
602  aProRes.SetMomentum(pFragment);
603 
604  proFrag = theHandler->BreakItUp(aProRes);
605 
606  boost_fragments = G4LorentzRotation(pSpectators.boostVector());
607 
608  // G4cout << " Fragment a,z, Mass Fragment, mass spect-mom, exitationE "
609  // << spectatorA <<" "<< spectatorZ <<" "<< mFragment <<" "
610  // << momentum.mag() <<" "<< momentum.mag() - mFragment
611  // << " "<<theStatisticalExEnergy
612  // << " "<< boost_fragments*pFragment<< G4endl;
613  G4ReactionProductVector::iterator ispectator;
614  for (ispectator=spectators->begin();ispectator!=spectators->end();ispectator++)
615  {
616  delete *ispectator;
617  }
618  }
619  else if(spectatorA!=0)
620  {
621  G4ReactionProductVector::iterator ispectator;
622  for (ispectator=spectators->begin();ispectator!=spectators->end();ispectator++)
623  {
624  (*ispectator)->SetNewlyAdded(true);
625  cascaders->push_back(*ispectator);
626  pFragments+=G4LorentzVector((*ispectator)->GetMomentum(),(*ispectator)->GetTotalEnergy());
627  // G4cout << "from spectator "
628  // << (*ispectator)->GetDefinition()->GetParticleName() << " "
629  // << (*ispectator)->GetMomentum()<< " "
630  // << (*ispectator)->GetTotalEnergy() << G4endl;
631  }
632  }
633  // / if (spectators)
634  delete spectators;
635 
636  // collect the evaporation part and boost to spectator frame
637  G4ReactionProductVector::iterator ii;
638  if(proFrag)
639  {
640  for(ii=proFrag->begin(); ii!=proFrag->end(); ii++)
641  {
642  (*ii)->SetNewlyAdded(true);
643  G4LorentzVector tmp((*ii)->GetMomentum(),(*ii)->GetTotalEnergy());
644  tmp *= boost_fragments;
645  (*ii)->SetMomentum(tmp.vect());
646  (*ii)->SetTotalEnergy(tmp.e());
647  // result->push_back(*ii);
648  pFragments += tmp;
649  }
650  }
651 
652  // G4cout << "Fragmented p, momentum, delta " << pFragments <<" "<<momentum
653  // <<" "<< pFragments-momentum << G4endl;
654 
655  // correct p/E of Cascade secondaries
656  G4LorentzVector pCas=pInitialState - pFragments;
657 
658  //G4cout <<" Going to correct from " << pFinalState << " to " << pCas << G4endl;
659  G4bool EnergyIsCorrect=EnergyAndMomentumCorrector(cascaders, pCas);
660  if ( ! EnergyIsCorrect && debug_G4BinaryLightIonReactionResults)
661  {
662  G4cout << "G4BinaryLightIonReaction E/P correction for nucleus failed, will try to correct overall" << G4endl;
663  }
664 
665  // Add deexcitation secondaries
666  if(proFrag)
667  {
668  for(ii=proFrag->begin(); ii!=proFrag->end(); ii++)
669  {
670  cascaders->push_back(*ii);
671  }
672  delete proFrag;
673  }
674  //G4cout << "EnergyIsCorrect? " << EnergyIsCorrect << G4endl;
675  if ( ! EnergyIsCorrect )
676  {
677  //G4cout <<" ! EnergyIsCorrect " << pFinalState << " to " << pInitialState << G4endl;
678  if (! EnergyAndMomentumCorrector(cascaders,pInitialState))
679  {
680  if(debug_G4BinaryLightIonReactionResults)
681  G4cout << "G4BinaryLightIonReaction E/P corrections failed" << G4endl;
682  }
683  }
684 
685 }
686 
Definition: G4ping.hh:33
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
Hep3Vector boostVector() const
Float_t tmp
Definition: plot.C:37
G4BinaryLightIonReaction(G4VPreCompoundModel *ptr=0)
virtual void ModelDescription(std::ostream &) const
const G4LorentzVector & GetMomentum() const
Definition: G4Nucleon.hh:71
tuple a
Definition: test.py:11
CLHEP::Hep3Vector G4ThreeVector
std::ofstream outFile
Definition: GammaRayTel.cc:68
CLHEP::HepLorentzRotation G4LorentzRotation
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)=0
const char * p
Definition: xmltok.h:285
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4double G4NeutronHPJENDLHEData::G4double result
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:335
virtual const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:68
double getT() const
void SetProjectilePotential(const G4double aPotential)
int G4int
Definition: G4Types.hh:78
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState) const
HepLorentzRotation & rotateY(double delta)
void SetA(G4double value)
Definition: G4Fragment.hh:294
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4double density
Definition: TRTMaterials.hh:39
tuple b
Definition: test.py:12
G4ExcitationHandler * GetExcitationHandler() const
G4IonTable * GetIonTable() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
CascadeState SetState(const CascadeState new_state)
const G4ParticleDefinition * GetDefinition() const
double mag() const
bool G4bool
Definition: G4Types.hh:79
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1232
G4bool AreYouHit() const
Definition: G4Nucleon.hh:97
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:256
G4double GetKineticEnergy() const
#define FALSE
Definition: globals.hh:52
void Init(G4int theA, G4int theZ)
G4double GetFermiMomentum(G4double density)
#define TRUE
Definition: globals.hh:55
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:344
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *, G4V3DNucleus *)
const G4LorentzVector & Get4Momentum() const
G4LorentzVector Get4Momentum() const
HepLorentzRotation & rotateZ(double delta)
G4HadronicInteraction * FindModel(const G4String &name)
void Set4Momentum(const G4LorentzVector &momentum)
void SetEnergyChange(G4double anEnergy)
G4double GetTotalEnergy() const
void Init(G4int anA, G4int aZ)
G4double GetPDGMass() const
G4double GetDensity(const G4ThreeVector &aPosition) const
static G4ParticleTable * GetParticleTable()
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define debug
static G4HadronicInteractionRegistry * Instance()
Hep3Vector unit() const
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
simulatedPeak Scale(1/simulationNormalisationFactor)
double mag2() const
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
virtual G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
void SetZ(G4double value)
Definition: G4Fragment.hh:288
G4double GetOuterRadius()
HepLorentzRotation inverse() const
G4LorentzVector operator()(G4LorentzVector a, G4ReactionProduct *b)
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:349
T sqr(const T &x)
Definition: templates.hh:145
void setVect(const Hep3Vector &)
const G4VNuclearDensity * GetNuclearDensity() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
double mag() const
G4Nucleon * GetNextNucleon()
void SetMomentumChange(const G4ThreeVector &aV)
G4int GetNumberOfSecondaries() const
void AddSecondary(G4DynamicParticle *aP)
G4GLOB_DLL std::ostream G4cerr
G4double GetTotalEnergy() const
CLHEP::HepLorentzVector G4LorentzVector