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