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