Geant4 10.7.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 G4ReactionProductVectorPropagate (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0
 
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
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void PropagateModelDescription (std::ostream &outFile) 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 G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 ModelDescription (std::ostream &outFile) const
 
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 if ( theBIC_ID == -1 ) {
169#ifdef G4MULTITHREADED
170 G4MUTEXLOCK(&G4BinaryCascade::BICMutex);
171 if ( theBIC_ID == -1 ) {
172#endif
173 theBIC_ID = G4PhysicsModelCatalog::Register("Binary Cascade");
174#ifdef G4MULTITHREADED
175 }
176 G4MUTEXUNLOCK(&G4BinaryCascade::BICMutex);
177#endif
178 }
179
180}
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
void SetMaxEnergy(const G4double anEnergy)
static G4int Register(const G4String &)
G4VPreCompoundModel * GetDeExcitation() const
void SetDeExcitation(G4VPreCompoundModel *ptr)
G4ExcitationHandler * GetExcitationHandler() const

◆ ~G4BinaryCascade()

G4BinaryCascade::~G4BinaryCascade ( )
virtual

Definition at line 182 of file G4BinaryCascade.cc.

183{
184 ClearAndDestroy(&theTargetList);
185 ClearAndDestroy(&theSecondaryList);
186 ClearAndDestroy(&theCapturedList);
187 delete thePropagator;
188 delete theCollisionMgr;
189 for(auto & ptr : theImR) { delete ptr; }
190 theImR.clear();
191 delete theLateParticle;
192 delete theH1Scatterer;
193}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 263 of file G4BinaryCascade.cc.

266{
267 if(std::getenv("BCDEBUG") ) G4cerr << " ######### Binary Cascade Reaction starts ######### "<< G4endl;
268
269 G4LorentzVector initial4Momentum = aTrack.Get4Momentum();
270 const G4ParticleDefinition * definition = aTrack.GetDefinition();
271
272 if(initial4Momentum.e()-initial4Momentum.m()<theBCminP &&
273 ( definition==G4Neutron::NeutronDefinition() || definition==G4Proton::ProtonDefinition() ) )
274 {
275 return theDeExcitation->ApplyYourself(aTrack, aNucleus);
276 }
277
279 // initialize the G4V3DNucleus from G4Nucleus
281
282 // Build a KineticTrackVector with the G4Track
283 G4KineticTrackVector * secondaries;// = new G4KineticTrackVector;
284 G4ThreeVector initialPosition(0., 0., 0.); // will be set later
285
286 if(!std::getenv("I_Am_G4BinaryCascade_Developer") )
287 {
288 if(definition!=G4Neutron::NeutronDefinition() &&
289 definition!=G4Proton::ProtonDefinition() &&
290 definition!=G4PionPlus::PionPlusDefinition() &&
292 {
293 G4cerr << "You are trying to use G4BinaryCascade with " <<definition->GetParticleName()<<" as projectile."<<G4endl;
294 G4cerr << "G4BinaryCascade should not be used for projectiles other than nucleons or pions."<<G4endl;
295 G4cerr << "If you want to continue, please switch on the developer environment: "<<G4endl;
296 G4cerr << "setenv I_Am_G4BinaryCascade_Developer 1 "<<G4endl<<G4endl;
297 throw G4HadronicException(__FILE__, __LINE__, "G4BinaryCascade - used for unvalid particle type - Fatal");
298 }
299 }
300
301 // keep primary
302 thePrimaryType = definition;
303 thePrimaryEscape = false;
304
305 G4double timePrimary=aTrack.GetGlobalTime();
306
307 // try until an interaction will happen
308 G4ReactionProductVector * products=0;
309 G4int interactionCounter = 0,collisionLoopMaxCount;
310 do
311 {
312 // reset status that could be changed in previous loop event
313 theCollisionMgr->ClearAndDestroy();
314
315 if(products != 0)
316 { // free memory from previous loop event
317 ClearAndDestroy(products);
318 delete products;
319 products=0;
320 }
321
322 G4int massNumber=aNucleus.GetA_asInt();
323 the3DNucleus->Init(massNumber, aNucleus.GetZ_asInt());
324 thePropagator->Init(the3DNucleus);
325 G4KineticTrack * kt;
326 collisionLoopMaxCount = 200;
327 do // sample impact parameter until collisions are found
328 {
329 theCurrentTime=0;
330 G4double radius = the3DNucleus->GetOuterRadius()+3*fermi;
331 initialPosition=GetSpherePoint(1.1*radius, initial4Momentum); // get random position
332 kt = new G4KineticTrack(definition, 0., initialPosition, initial4Momentum);
334 // secondaries has been cleared by Propagate() in the previous loop event
335 secondaries= new G4KineticTrackVector;
336 secondaries->push_back(kt);
337 if(massNumber > 1) // 1H1 is special case
338 {
339 products = Propagate(secondaries, the3DNucleus);
340 } else {
341 products = Propagate1H1(secondaries,the3DNucleus);
342 }
343 // until we FIND a collision ... or give up
344 } while(! products && --collisionLoopMaxCount>0); /* Loop checking, 31.08.2015, G.Folger */
345
346 if(++interactionCounter>99) break;
347 // ...until we find an ALLOWED collision ... or give up
348 } while(products && products->size() == 0); /* Loop checking, 31.08.2015, G.Folger */
349
350 if(products && products->size()>0)
351 {
352 // G4cout << "BIC Applyyourself: number of products " << products->size() << G4endl;
353
354 // Fill the G4ParticleChange * with products
356 G4ReactionProductVector::iterator iter;
357
358 for(iter = products->begin(); iter != products->end(); ++iter)
359 {
360 G4DynamicParticle * aNewDP =
361 new G4DynamicParticle((*iter)->GetDefinition(),
362 (*iter)->GetTotalEnergy(),
363 (*iter)->GetMomentum());
364 G4HadSecondary aNew = G4HadSecondary(aNewDP);
365 G4double time=(*iter)->GetFormationTime();
366 if(time < 0.0) { time = 0.0; }
367 aNew.SetTime(timePrimary + time);
368 aNew.SetCreatorModelType((*iter)->GetCreatorModel());
370 }
371
372 //DebugFinalEpConservation(aTrack, products);
373
374
375 } else { // no interaction, return primary
376 if(std::getenv("BCDEBUG") ) G4cerr << " ######### Binary Cascade Reaction void, return initial state ######### "<< G4endl;
380 }
381
382 if ( products )
383 {
384 ClearAndDestroy(products);
385 delete products;
386 }
387
388 delete the3DNucleus;
389 the3DNucleus = NULL;
390
391 if(std::getenv("BCDEBUG") ) G4cerr << " ######### Binary Cascade Reaction ends ######### "<< G4endl;
392
393 return &theParticleChange;
394}
@ 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:57
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 SetTime(G4double aT)
void SetCreatorModelType(G4int idx)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
CascadeState SetState(const CascadeState new_state)
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
const G4String & GetParticleName() const
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:92
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:92
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
virtual G4double GetOuterRadius()=0
virtual void Init(G4int theA, G4int theZ)=0
virtual void Init(G4V3DNucleus *theNucleus)=0

◆ ModelDescription()

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

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 195 of file G4BinaryCascade.cc.

196{
197 outFile << "G4BinaryCascade is an intra-nuclear cascade model in which\n"
198 << "an incident hadron collides with a nucleon, forming two\n"
199 << "final-state particles, one or both of which may be resonances.\n"
200 << "The resonances then decay hadronically and the decay products\n"
201 << "are then propagated through the nuclear potential along curved\n"
202 << "trajectories until they re-interact or leave the nucleus.\n"
203 << "This model is valid for incident pions up to 1.5 GeV and\n"
204 << "nucleons up to 10 GeV.\n"
205 << "The remaining excited nucleus is handed on to ";
206 if (theDeExcitation) // pre-compound
207 {
208 outFile << theDeExcitation->GetModelName() << " : \n ";
210 }
211 else if (theExcitationHandler) // de-excitation
212 {
213 outFile << "G4ExcitationHandler"; //theExcitationHandler->GetModelName();
214 theExcitationHandler->ModelDescription(outFile);
215 }
216 else
217 {
218 outFile << "void.\n";
219 }
220 outFile<< " \n";
221}
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 396 of file G4BinaryCascade.cc.

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

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

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