Geant4 11.1.1
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 119 of file G4BinaryCascade.cc.

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

◆ ~G4BinaryCascade()

G4BinaryCascade::~G4BinaryCascade ( )
virtual

Definition at line 170 of file G4BinaryCascade.cc.

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

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 251 of file G4BinaryCascade.cc.

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

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

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

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

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