Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BinaryCascade Class Reference

#include <G4BinaryCascade.hh>

Inheritance diagram for G4BinaryCascade:

Public Member Functions

 G4BinaryCascade (G4VPreCompoundModel *ptr=0)
virtual ~G4BinaryCascade ()
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
virtual G4ReactionProductVectorPropagate (G4KineticTrackVector *, G4V3DNucleus *)
virtual void ModelDescription (std::ostream &) const
virtual void PropagateModelDescription (std::ostream &) const
Public Member Functions inherited from G4VIntraNuclearTransportModel
 G4VIntraNuclearTransportModel (const G4String &mName="CascadeModel", G4VPreCompoundModel *ptr=nullptr)
virtual ~G4VIntraNuclearTransportModel ()
virtual G4ReactionProductVectorPropagateNuclNucl (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
void SetDeExcitation (G4VPreCompoundModel *ptr)
void Set3DNucleus (G4V3DNucleus *const value)
void SetPrimaryProjectile (const G4HadProjectile &aPrimary)
const G4StringGetModelName () const
 G4VIntraNuclearTransportModel (const G4VIntraNuclearTransportModel &right)=delete
const G4VIntraNuclearTransportModeloperator= (const G4VIntraNuclearTransportModel &right)=delete
G4bool operator== (const G4VIntraNuclearTransportModel &right) const =delete
G4bool operator!= (const G4VIntraNuclearTransportModel &right) const =delete
Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
virtual ~G4HadronicInteraction ()
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4double GetMinEnergy () const
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
void SetMinEnergy (G4double anEnergy)
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
G4double GetMaxEnergy () const
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
void SetMaxEnergy (const G4double anEnergy)
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
G4int GetVerboseLevel () const
void SetVerboseLevel (G4int value)
const G4StringGetModelName () const
void DeActivateFor (const G4Material *aMaterial)
void ActivateFor (const G4Material *aMaterial)
void DeActivateFor (const G4Element *anElement)
void ActivateFor (const G4Element *anElement)
G4bool IsBlocked (const G4Material *aMaterial) const
G4bool IsBlocked (const G4Element *anElement) const
void SetRecoilEnergyThreshold (G4double val)
G4double GetRecoilEnergyThreshold () const
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual void InitialiseModel ()
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
G4bool operator== (const G4HadronicInteraction &right) const =delete
G4bool operator!= (const G4HadronicInteraction &right) const =delete

Additional Inherited Members

Protected Member Functions inherited from G4VIntraNuclearTransportModel
G4V3DNucleusGet3DNucleus () const
G4VPreCompoundModelGetDeExcitation () const
const G4HadProjectileGetPrimaryProjectile () const
Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
G4bool IsBlocked () const
void Block ()
Protected Attributes inherited from G4VIntraNuclearTransportModel
G4String theTransportModelName
G4V3DNucleusthe3DNucleus
G4VPreCompoundModeltheDeExcitation
const G4HadProjectilethePrimaryProjectile
Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
G4int verboseLevel
G4double theMinEnergy
G4double theMaxEnergy
G4bool isBlocked

Detailed Description

Definition at line 72 of file G4BinaryCascade.hh.

Constructor & Destructor Documentation

◆ G4BinaryCascade()

G4BinaryCascade::G4BinaryCascade ( G4VPreCompoundModel * ptr = 0)

Definition at line 120 of file G4BinaryCascade.cc.

120 :
121G4VIntraNuclearTransportModel("Binary Cascade", ptr)
122{
123 // initialise the resonance sector
124 G4ShortLivedConstructor ShortLived;
125 ShortLived.ConstructParticle();
126
127 theCollisionMgr = new G4CollisionManager;
128 theDecay=new G4BCDecay;
129 theImR.push_back(theDecay);
130 theLateParticle= new G4BCLateParticle;
131 G4MesonAbsorption * aAb=new G4MesonAbsorption;
132 theImR.push_back(aAb);
133 G4Scatterer * aSc=new G4Scatterer;
134 theH1Scatterer = new G4Scatterer;
135 theImR.push_back(aSc);
136
137 thePropagator = new G4RKPropagation;
138 theCurrentTime = 0.;
139 theBCminP = 45*MeV;
140 theCutOnP = 90*MeV;
141 theCutOnPAbsorb= 0*MeV; // No Absorption of slow Mesons, other than above G4MesonAbsorption
142
143 // reuse existing pre-compound model
144 if(!ptr) {
147 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
148 if(!pre) { pre = new G4PreCompoundModel(); }
149 SetDeExcitation(pre);
150 }
151 theExcitationHandler = GetDeExcitation()->GetExcitationHandler();
152 SetMinEnergy(0.0*GeV);
153 SetMaxEnergy(10.1*GeV);
154 //PrintWelcomeMessage();
155 thePrimaryEscape = true;
156 thePrimaryType = 0;
157
158 SetEnergyMomentumCheckLevels(1.0*perCent, 1.0*MeV);
159
160 // init data members
161 currentA=currentZ=0;
162 lateA=lateZ=0;
163 initialA=initialZ=0;
164 projectileA=projectileZ=0;
165 currentInitialEnergy=initial_nuclear_mass=0.;
166 massInNucleus=0.;
167 theOuterRadius=0.;
168 theBIC_ID = G4PhysicsModelCatalog::GetModelID("model_G4BinaryCascade");
170}
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4int GetModelID(const G4int modelIndex)
G4VIntraNuclearTransportModel(const G4String &mName="CascadeModel", G4VPreCompoundModel *ptr=nullptr)
G4VPreCompoundModel * GetDeExcitation() const
void SetDeExcitation(G4VPreCompoundModel *ptr)
G4ExcitationHandler * GetExcitationHandler() const

◆ ~G4BinaryCascade()

G4BinaryCascade::~G4BinaryCascade ( )
virtual

Definition at line 172 of file G4BinaryCascade.cc.

173{
174 ClearAndDestroy(&theTargetList);
175 ClearAndDestroy(&theSecondaryList);
176 ClearAndDestroy(&theCapturedList);
177 delete thePropagator;
178 delete theCollisionMgr;
179 for(auto & ptr : theImR) { delete ptr; }
180 theImR.clear();
181 delete theLateParticle;
182 delete theH1Scatterer;
183}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4BinaryCascade::ApplyYourself ( const G4HadProjectile & aTrack,
G4Nucleus & theNucleus )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 253 of file G4BinaryCascade.cc.

256{
257 if(fBCDEBUG) G4cerr << " ######### Binary Cascade Reaction starts ######### "<< G4endl;
258
259 G4LorentzVector initial4Momentum = aTrack.Get4Momentum();
260 const G4ParticleDefinition * definition = aTrack.GetDefinition();
261
262 if(initial4Momentum.e()-initial4Momentum.m()<theBCminP &&
263 ( definition==G4Neutron::NeutronDefinition() || definition==G4Proton::ProtonDefinition() ) )
264 {
265 return theDeExcitation->ApplyYourself(aTrack, aNucleus);
266 }
267
268 theParticleChange.Clear();
269 // initialize the G4V3DNucleus from G4Nucleus
270 the3DNucleus = new G4Fancy3DNucleus;
271
272 // Build a KineticTrackVector with the G4Track
273 G4KineticTrackVector * secondaries;// = new G4KineticTrackVector;
274 G4ThreeVector initialPosition(0., 0., 0.); // will be set later
275
276 if(!fBCDEBUG)
277 {
278 if(definition!=G4Neutron::NeutronDefinition() &&
279 definition!=G4Proton::ProtonDefinition() &&
280 definition!=G4PionPlus::PionPlusDefinition() &&
282 {
283 G4cerr << "You are trying to use G4BinaryCascade with " <<definition->GetParticleName()<<" as projectile."<<G4endl;
284 G4cerr << "G4BinaryCascade should not be used for projectiles other than nucleons or pions."<<G4endl;
285 G4cerr << "If you want to continue, please switch on the developer environment: "<<G4endl;
286 G4cerr << "setenv I_Am_G4BinaryCascade_Developer 1 "<<G4endl<<G4endl;
287 throw G4HadronicException(__FILE__, __LINE__, "G4BinaryCascade - used for unvalid particle type - Fatal");
288 }
289 }
290
291 // keep primary
292 thePrimaryType = definition;
293 thePrimaryEscape = false;
294
295 G4double timePrimary=aTrack.GetGlobalTime();
296
297 // try until an interaction will happen
298 G4ReactionProductVector * products = nullptr;
299 G4int interactionCounter = 0,collisionLoopMaxCount;
300 do
301 {
302 // reset status that could be changed in previous loop event
303 theCollisionMgr->ClearAndDestroy();
304
305 if(products != nullptr)
306 { // free memory from previous loop event
307 ClearAndDestroy(products);
308 delete products;
309 }
310
311 G4int massNumber=aNucleus.GetA_asInt();
312 the3DNucleus->Init(massNumber, aNucleus.GetZ_asInt());
313 thePropagator->Init(the3DNucleus);
314 G4KineticTrack * kt;
315 collisionLoopMaxCount = 200;
316 do // sample impact parameter until collisions are found
317 {
318 theCurrentTime=0;
319 G4double radius = the3DNucleus->GetOuterRadius()+3*fermi;
320 initialPosition=GetSpherePoint(1.1*radius, initial4Momentum); // get random position
321 kt = new G4KineticTrack(definition, 0., initialPosition, initial4Momentum);
323 // secondaries has been cleared by Propagate() in the previous loop event
324 secondaries= new G4KineticTrackVector;
325 secondaries->push_back(kt);
326 if(massNumber > 1) // 1H1 is special case
327 {
328 products = Propagate(secondaries, the3DNucleus);
329 } else {
330 products = Propagate1H1(secondaries,the3DNucleus);
331 }
332 // until we FIND a collision ... or give up
333 } while(! products && --collisionLoopMaxCount>0); /* Loop checking, 31.08.2015, G.Folger */
334
335 if(++interactionCounter>99) break;
336 // ...until we find an ALLOWED collision ... or give up
337 } while(products && products->size() == 0); /* Loop checking, 31.08.2015, G.Folger */
338
339 if(products && products->size()>0)
340 {
341 // G4cout << "BIC Applyyourself: number of products " << products->size() << G4endl;
342
343 // Fill the G4ParticleChange * with products
344 theParticleChange.SetStatusChange(stopAndKill);
345 G4ReactionProductVector::iterator iter;
346
347 for(iter = products->begin(); iter != products->end(); ++iter)
348 {
349 G4DynamicParticle * aNewDP =
350 new G4DynamicParticle((*iter)->GetDefinition(),
351 (*iter)->GetTotalEnergy(),
352 (*iter)->GetMomentum());
353 G4HadSecondary aNew = G4HadSecondary(aNewDP);
354 G4double time=(*iter)->GetFormationTime();
355 if(time < 0.0) { time = 0.0; }
356 aNew.SetTime(timePrimary + time);
357 aNew.SetCreatorModelID((*iter)->GetCreatorModelID());
358 aNew.SetParentResonanceDef((*iter)->GetParentResonanceDef());
359 aNew.SetParentResonanceID((*iter)->GetParentResonanceID());
360 theParticleChange.AddSecondary(aNew);
361 }
362
363 //DebugFinalEpConservation(aTrack, products);
364
365
366 } else { // no interaction, return primary
367 if(fBCDEBUG) G4cerr << " ######### Binary Cascade Reaction void, return initial state ######### "<< G4endl;
368 theParticleChange.SetStatusChange(isAlive);
369 theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
370 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
371 }
372
373 if ( products )
374 {
375 ClearAndDestroy(products);
376 delete products;
377 }
378
379 delete the3DNucleus;
380 the3DNucleus = nullptr;
381
382 if(fBCDEBUG) G4cerr << " ######### Binary Cascade Reaction ends ######### "<< G4endl;
383
384 return &theParticleChange;
385}
@ isAlive
@ stopAndKill
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
Hep3Vector unit() const
Hep3Vector vect() const
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *, G4V3DNucleus *)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetParentResonanceDef(const G4ParticleDefinition *parentDef)
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
void SetParentResonanceID(const G4int parentID)
CascadeState SetState(const CascadeState new_state)
static G4Neutron * NeutronDefinition()
Definition G4Neutron.cc:96
const G4String & GetParticleName() const
static G4PionMinus * PionMinusDefinition()
static G4PionPlus * PionPlusDefinition()
Definition G4PionPlus.cc:88
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85

◆ ModelDescription()

void G4BinaryCascade::ModelDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 185 of file G4BinaryCascade.cc.

186{
187 outFile << "G4BinaryCascade is an intra-nuclear cascade model in which\n"
188 << "an incident hadron collides with a nucleon, forming two\n"
189 << "final-state particles, one or both of which may be resonances.\n"
190 << "The resonances then decay hadronically and the decay products\n"
191 << "are then propagated through the nuclear potential along curved\n"
192 << "trajectories until they re-interact or leave the nucleus.\n"
193 << "This model is valid for incident pions up to 1.5 GeV and\n"
194 << "nucleons up to 10 GeV.\n"
195 << "The remaining excited nucleus is handed on to ";
196 if (theDeExcitation) // pre-compound
197 {
198 outFile << theDeExcitation->GetModelName() << " : \n ";
199 theDeExcitation->DeExciteModelDescription(outFile);
200 }
201 else if (theExcitationHandler) // de-excitation
202 {
203 outFile << "G4ExcitationHandler"; //theExcitationHandler->GetModelName();
204 theExcitationHandler->ModelDescription(outFile);
205 }
206 else
207 {
208 outFile << "void.\n";
209 }
210 outFile<< " \n";
211}

◆ Propagate()

G4ReactionProductVector * G4BinaryCascade::Propagate ( G4KineticTrackVector * secondaries,
G4V3DNucleus * aNucleus )
virtual

Implements G4VIntraNuclearTransportModel.

Definition at line 387 of file G4BinaryCascade.cc.

390{
391 G4ping debug("debug_G4BinaryCascade");
392#ifdef debug_BIC_Propagate
393 G4cout << "G4BinaryCascade Propagate starting -------------------------------------------------------" <<G4endl;
394#endif
395
396 the3DNucleus=aNucleus;
398 theOuterRadius = the3DNucleus->GetOuterRadius();
399 theCurrentTime=0;
400 theProjectile4Momentum=G4LorentzVector(0,0,0,0);
401 theMomentumTransfer=G4ThreeVector(0,0,0);
402 // build theSecondaryList, theProjectileList and theCapturedList
403 ClearAndDestroy(&theCapturedList);
404 ClearAndDestroy(&theSecondaryList);
405 theSecondaryList.clear();
406 ClearAndDestroy(&theFinalState);
407 std::vector<G4KineticTrack *>::iterator iter;
408 theCollisionMgr->ClearAndDestroy();
409
410 theCutOnP=90*MeV;
411 if(the3DNucleus->GetMass()>30) theCutOnP = 70*MeV;
412 if(the3DNucleus->GetMass()>60) theCutOnP = 50*MeV;
413 if(the3DNucleus->GetMass()>120) theCutOnP = 45*MeV;
414
415
416 BuildTargetList();
417
418#ifdef debug_BIC_GetExcitationEnergy
419 G4cout << "ExcitationEnergy0 " << GetExcitationEnergy() << G4endl;
420#endif
421
422 thePropagator->Init(the3DNucleus);
423
424 G4bool success = BuildLateParticleCollisions(secondaries);
425 if (! success ) // fails if no excitation energy left....
426 {
427 products=HighEnergyModelFSProducts(products, secondaries);
428 ClearAndDestroy(secondaries);
429 delete secondaries;
430
431#ifdef debug_G4BinaryCascade
432 G4cout << "G4BinaryCascade::Propagate: warning - high energy model failed energy conservation, returning unchanged high energy final state" << G4endl;
433#endif
434
435 return products;
436 }
437 // check baryon and charge ...
438
439 _CheckChargeAndBaryonNumber_("lateparticles");
440 _DebugEpConservation(" be4 findcollisions");
441
442 // if called stand alone find first collisions
443 FindCollisions(&theSecondaryList);
444
445
446 if(theCollisionMgr->Entries() == 0 ) //late particles ALWAYS create Entries
447 {
448 //G4cout << " no collsions -> return 0" << G4endl;
449 delete products;
450#ifdef debug_BIC_return
451 G4cout << "return @ begin2, no collisions "<< G4endl;
452#endif
453 return nullptr;
454 }
455
456 // end of initialization: do the job now
457 // loop until there are no more collisions
458
459
460 G4bool haveProducts = false;
461#ifdef debug_BIC_Propagate_Collisions
462 G4int collisionCount=0;
463#endif
464 G4int collisionLoopMaxCount=1000000;
465 while(theCollisionMgr->Entries() > 0 && currentZ && --collisionLoopMaxCount>0) /* Loop checking, 31.08.2015, G.Folger */
466 {
467 if(Absorb()) { // absorb secondaries, pions only
468 haveProducts = true;
469 }
470 if(Capture()) { // capture secondaries, nucleons only
471 haveProducts = true;
472 }
473
474 // propagate to the next collision if any (collisions could have been deleted
475 // by previous absorption or capture)
476 if(theCollisionMgr->Entries() > 0)
477 {
478 G4CollisionInitialState *
479 nextCollision = theCollisionMgr->GetNextCollision();
480#ifdef debug_BIC_Propagate_Collisions
481 G4cout << " NextCollision * , Time, curtime = " << nextCollision << " "
482 <<nextCollision->GetCollisionTime()<< " " <<
483 theCurrentTime<< G4endl;
484#endif
485 if (!DoTimeStep(nextCollision->GetCollisionTime()-theCurrentTime) )
486 {
487 // Check if nextCollision is still valid, ie. particle did not leave nucleus
488 if (theCollisionMgr->GetNextCollision() != nextCollision )
489 {
490 nextCollision = nullptr;
491 }
492 }
493 //_DebugEpConservation("Stepped");
494
495 if( nextCollision )
496 {
497 if (ApplyCollision(nextCollision))
498 {
499 //G4cerr << "ApplyCollision success " << G4endl;
500 haveProducts = true;
501#ifdef debug_BIC_Propagate_Collisions
502 collisionCount++;
503#endif
504
505 } else {
506 //G4cerr << "ApplyCollision failure " << G4endl;
507 theCollisionMgr->RemoveCollision(nextCollision);
508 }
509 }
510 }
511 }
512
513 //--------- end of on Collisions
514 //G4cout << "currentZ @ end loop " << currentZ << G4endl;
515 G4int nProtons(0);
516 for(iter = theTargetList.begin(); iter != theTargetList.end(); ++iter)
517 {
518 if ( (*iter)->GetDefinition() == G4Proton::Proton() ) ++nProtons;
519 }
520 if ( ! theTargetList.size() || ! nProtons ){
521 // nucleus completely destroyed, fill in ReactionProductVector
522 products = FillVoidNucleusProducts(products);
523#ifdef debug_BIC_return
524 G4cout << "return @ Z=0 after collision loop "<< G4endl;
525 PrintKTVector(&theSecondaryList,std::string(" theSecondaryList"));
526 G4cout << "theTargetList size: " << theTargetList.size() << G4endl;
527 PrintKTVector(&theTargetList,std::string(" theTargetList"));
528 PrintKTVector(&theCapturedList,std::string(" theCapturedList"));
529
530 G4cout << " ExcitE be4 Correct : " <<GetExcitationEnergy() << G4endl;
531 G4cout << " Mom Transfered to nucleus : " << theMomentumTransfer << " " << theMomentumTransfer.mag() << G4endl;
532 PrintKTVector(&theFinalState,std::string(" FinalState uncorrected"));
533 G4cout << "returned products: " << products->size() << G4endl;
534 _CheckChargeAndBaryonNumber_("destroyed Nucleus");
535 _DebugEpConservation("destroyed Nucleus");
536#endif
537
538 return products;
539 }
540
541 // No more collisions: absorb, capture and propagate the secondaries out of the nucleus
542 if(Absorb()) {
543 haveProducts = true;
544 // G4cout << "Absorb sucess " << G4endl;
545 }
546
547 if(Capture()) {
548 haveProducts = true;
549 // G4cout << "Capture sucess " << G4endl;
550 }
551
552 if(!haveProducts) // no collisions happened. Return an empty vector.
553 {
554#ifdef debug_BIC_return
555 G4cout << "return 3, no products "<< G4endl;
556#endif
557 return products;
558 }
559
560
561#ifdef debug_BIC_Propagate
562 G4cout << " Momentum transfer to Nucleus " << theMomentumTransfer << " " << theMomentumTransfer.mag() << G4endl;
563 G4cout << " Stepping particles out...... " << G4endl;
564#endif
565
566 StepParticlesOut();
567 _DebugEpConservation("stepped out");
568
569
570 if ( theSecondaryList.size() > 0 )
571 {
572#ifdef debug_G4BinaryCascade
573 G4cerr << "G4BinaryCascade: Warning, have active particles at end" << G4endl;
574 PrintKTVector(&theSecondaryList, "active particles @ end added to theFinalState");
575#endif
576 // add left secondaries to FinalSate
577 for ( iter =theSecondaryList.begin(); iter != theSecondaryList.end(); ++iter)
578 {
579 theFinalState.push_back(*iter);
580 }
581 theSecondaryList.clear();
582
583 }
584 while ( theCollisionMgr->Entries() > 0 ) /* Loop checking, 31.08.2015, G.Folger */
585 {
586#ifdef debug_G4BinaryCascade
587 G4cerr << " Warning: remove left over collision(s) " << G4endl;
588#endif
589 theCollisionMgr->RemoveCollision(theCollisionMgr->GetNextCollision());
590 }
591
592#ifdef debug_BIC_Propagate_Excitation
593
594 PrintKTVector(&theSecondaryList,std::string(" theSecondaryList"));
595 G4cout << "theTargetList size: " << theTargetList.size() << G4endl;
596 // PrintKTVector(&theTargetList,std::string(" theTargetList"));
597 PrintKTVector(&theCapturedList,std::string(" theCapturedList"));
598
599 G4cout << " ExcitE be4 Correct : " <<GetExcitationEnergy() << G4endl;
600 G4cout << " Mom Transfered to nucleus : " << theMomentumTransfer << " " << theMomentumTransfer.mag() << G4endl;
601 PrintKTVector(&theFinalState,std::string(" FinalState uncorrected"));
602#endif
603
604 //
605
606
607 G4double ExcitationEnergy=GetExcitationEnergy();
608
609#ifdef debug_BIC_Propagate_finals
610 PrintKTVector(&theFinalState,std::string(" FinalState be4 corr"));
611 G4cout << " Excitation Energy prefinal, #collisions:, out, captured "
612 << ExcitationEnergy << " "
613 << collisionCount << " "
614 << theFinalState.size() << " "
615 << theCapturedList.size()<<G4endl;
616#endif
617
618 if (ExcitationEnergy < 0 )
619 {
620 G4int maxtry=5, ntry=0;
621 do {
622 CorrectFinalPandE();
623 ExcitationEnergy=GetExcitationEnergy();
624 } while ( ++ntry < maxtry && ExcitationEnergy < 0 ); /* Loop checking, 31.08.2015, G.Folger */
625 }
626 _DebugEpConservation("corrected");
627
628#ifdef debug_BIC_Propagate_finals
629 PrintKTVector(&theFinalState,std::string(" FinalState corrected"));
630 G4cout << " Excitation Energy final, #collisions:, out, captured "
631 << ExcitationEnergy << " "
632 << collisionCount << " "
633 << theFinalState.size() << " "
634 << theCapturedList.size()<<G4endl;
635#endif
636
637
638 if ( ExcitationEnergy < 0. )
639 {
640 #ifdef debug_G4BinaryCascade
641 G4cerr << "G4BinaryCascade-Warning: negative excitation energy ";
642 G4cerr <<ExcitationEnergy<<G4endl;
643 PrintKTVector(&theFinalState,std::string("FinalState"));
644 PrintKTVector(&theCapturedList,std::string("captured"));
645 G4cout << "negative ExE:Final 4Momentum .mag: " << GetFinal4Momentum()
646 << " "<< GetFinal4Momentum().mag()<< G4endl
647 << "negative ExE:FinalNucleusMom .mag: " << GetFinalNucleusMomentum()
648 << " "<< GetFinalNucleusMomentum().mag()<< G4endl;
649 #endif
650 #ifdef debug_BIC_return
651 G4cout << " negative Excitation E return empty products " << products << G4endl;
652 G4cout << "return 4, excit < 0 "<< G4endl;
653 #endif
654
655 ClearAndDestroy(products);
656 return products; // return empty products- FixMe
657 }
658
659 G4ReactionProductVector * precompoundProducts=DeExcite();
660
661
662 G4DecayKineticTracks decay(&theFinalState);
663 _DebugEpConservation("decayed");
664
665 products= ProductsAddFinalState(products, theFinalState);
666
667 products= ProductsAddPrecompound(products, precompoundProducts);
668
669// products=ProductsAddFakeGamma(products);
670
671
672 thePrimaryEscape = true;
673
674 #ifdef debug_BIC_return
675 G4cout << "BIC: return @end, all ok "<< G4endl;
676 //G4cout << " return products " << products << G4endl;
677 #endif
678
679 return products;
680}
#define _CheckChargeAndBaryonNumber_(val)
#define _DebugEpConservation(val)
bool G4bool
Definition G4Types.hh:86
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition G4Proton.cc:90
ParticleList decay(Cluster *const c)
Carries out a cluster decay.

Referenced by ApplyYourself().

◆ PropagateModelDescription()

void G4BinaryCascade::PropagateModelDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 212 of file G4BinaryCascade.cc.

213{
214 outFile << "G4BinaryCascade propagtes secondaries produced by a high\n"
215 << "energy model through the wounded nucleus.\n"
216 << "Secondaries are followed after the formation time and if\n"
217 << "within the nucleus are propagated through the nuclear\n"
218 << "potential along curved trajectories until they interact\n"
219 << "with a nucleon, decay, or leave the nucleus.\n"
220 << "An interaction of a secondary with a nucleon produces two\n"
221 << "final-state particles, one or both of which may be resonances.\n"
222 << "Resonances decay hadronically and the decay products\n"
223 << "are in turn propagated through the nuclear potential along curved\n"
224 << "trajectories until they re-interact or leave the nucleus.\n"
225 << "This model is valid for pions up to 1.5 GeV and\n"
226 << "nucleons up to about 3.5 GeV.\n"
227 << "The remaining excited nucleus is handed on to ";
228 if (theDeExcitation) // pre-compound
229 {
230 outFile << theDeExcitation->GetModelName() << " : \n ";
231 theDeExcitation->DeExciteModelDescription(outFile);
232 }
233 else if (theExcitationHandler) // de-excitation
234 {
235 outFile << "G4ExcitationHandler"; //theExcitationHandler->GetModelName();
236 theExcitationHandler->ModelDescription(outFile);
237 }
238 else
239 {
240 outFile << "void.\n";
241 }
242outFile<< " \n";
243}

The documentation for this class was generated from the following files: