Geant4 11.2.2
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;
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)
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
269 // initialize the G4V3DNucleus from G4Nucleus
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=0;
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 != 0)
306 { // free memory from previous loop event
307 ClearAndDestroy(products);
308 delete products;
309 products=0;
310 }
311
312 G4int massNumber=aNucleus.GetA_asInt();
313 the3DNucleus->Init(massNumber, aNucleus.GetZ_asInt());
314 thePropagator->Init(the3DNucleus);
315 G4KineticTrack * kt;
316 collisionLoopMaxCount = 200;
317 do // sample impact parameter until collisions are found
318 {
319 theCurrentTime=0;
320 G4double radius = the3DNucleus->GetOuterRadius()+3*fermi;
321 initialPosition=GetSpherePoint(1.1*radius, initial4Momentum); // get random position
322 kt = new G4KineticTrack(definition, 0., initialPosition, initial4Momentum);
324 // secondaries has been cleared by Propagate() in the previous loop event
325 secondaries= new G4KineticTrackVector;
326 secondaries->push_back(kt);
327 if(massNumber > 1) // 1H1 is special case
328 {
329 products = Propagate(secondaries, the3DNucleus);
330 } else {
331 products = Propagate1H1(secondaries,the3DNucleus);
332 }
333 // until we FIND a collision ... or give up
334 } while(! products && --collisionLoopMaxCount>0); /* Loop checking, 31.08.2015, G.Folger */
335
336 if(++interactionCounter>99) break;
337 // ...until we find an ALLOWED collision ... or give up
338 } while(products && products->size() == 0); /* Loop checking, 31.08.2015, G.Folger */
339
340 if(products && products->size()>0)
341 {
342 // G4cout << "BIC Applyyourself: number of products " << products->size() << G4endl;
343
344 // Fill the G4ParticleChange * with products
346 G4ReactionProductVector::iterator iter;
347
348 for(iter = products->begin(); iter != products->end(); ++iter)
349 {
350 G4DynamicParticle * aNewDP =
351 new G4DynamicParticle((*iter)->GetDefinition(),
352 (*iter)->GetTotalEnergy(),
353 (*iter)->GetMomentum());
354 G4HadSecondary aNew = G4HadSecondary(aNewDP);
355 G4double time=(*iter)->GetFormationTime();
356 if(time < 0.0) { time = 0.0; }
357 aNew.SetTime(timePrimary + time);
358 aNew.SetCreatorModelID((*iter)->GetCreatorModelID());
359 aNew.SetParentResonanceDef((*iter)->GetParentResonanceDef());
360 aNew.SetParentResonanceID((*iter)->GetParentResonanceID());
362 }
363
364 //DebugFinalEpConservation(aTrack, products);
365
366
367 } else { // no interaction, return primary
368 if(fBCDEBUG) G4cerr << " ######### Binary Cascade Reaction void, return initial state ######### "<< G4endl;
372 }
373
374 if ( products )
375 {
376 ClearAndDestroy(products);
377 delete products;
378 }
379
380 delete the3DNucleus;
381 the3DNucleus = nullptr;
382
383 if(fBCDEBUG) G4cerr << " ######### Binary Cascade Reaction ends ######### "<< G4endl;
384
385 return &theParticleChange;
386}
@ isAlive
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
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 *)
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 SetParentResonanceDef(const G4ParticleDefinition *parentDef)
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
void SetParentResonanceID(const G4int parentID)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
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
virtual G4double GetOuterRadius()=0
virtual void Init(G4int theA, G4int theZ, G4int numberOfLambdas=0)=0
virtual void Init(G4V3DNucleus *theNucleus)=0

◆ 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 ";
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}
void ModelDescription(std::ostream &outFile) const
const G4String & GetModelName() const
virtual void DeExciteModelDescription(std::ostream &outFile) const =0

◆ Propagate()

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

Implements G4VIntraNuclearTransportModel.

Definition at line 388 of file G4BinaryCascade.cc.

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