Geant4 9.6.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 &modelName="CascadeModel", G4VPreCompoundModel *ptr=0)
 
virtual ~G4VIntraNuclearTransportModel ()
 
virtual G4ReactionProductVectorPropagate (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0
 
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
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
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)
 
const G4HadronicInteractionGetMyPointer () const
 
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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) 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
 

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 70 of file G4BinaryCascade.hh.

Constructor & Destructor Documentation

◆ G4BinaryCascade()

G4BinaryCascade::G4BinaryCascade ( G4VPreCompoundModel ptr = 0)

Definition at line 102 of file G4BinaryCascade.cc.

102 :
103G4VIntraNuclearTransportModel("Binary Cascade", ptr)
104{
105 // initialise the resonance sector
106 G4ShortLivedConstructor ShortLived;
107 ShortLived.ConstructParticle();
108
109 theCollisionMgr = new G4CollisionManager;
110 theDecay=new G4BCDecay;
111 theImR.push_back(theDecay);
112 theLateParticle= new G4BCLateParticle;
114 theImR.push_back(aAb);
115 G4Scatterer * aSc=new G4Scatterer;
116 theH1Scatterer = new G4Scatterer;
117 theImR.push_back(aSc);
118
119 thePropagator = new G4RKPropagation;
120 theCurrentTime = 0.;
121 theBCminP = 45*MeV;
122 theCutOnP = 90*MeV;
123 theCutOnPAbsorb= 0*MeV; // No Absorption of slow Mesons, other than above G4MesonAbsorption
124
125 // reuse existing pre-compound model
126 if(!ptr) {
129 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
130 if(!pre) { pre = new G4PreCompoundModel(); }
131 SetDeExcitation(pre);
132 }
133 theExcitationHandler = GetDeExcitation()->GetExcitationHandler();
134 SetMinEnergy(0.0*GeV);
135 SetMaxEnergy(10.1*GeV);
136 //PrintWelcomeMessage();
137 thePrimaryEscape = true;
138 thePrimaryType = 0;
139
140 SetEnergyMomentumCheckLevels(1*perCent, 1*MeV);
141
142 // init data members
143 currentA=currentZ=0;
144 lateA=lateZ=0;
145 initialA=initialZ=0;
146 projectileA=projectileZ=0;
147 currentInitialEnergy=initial_nuclear_mass=0.;
148 massInNucleus=0.;
149 theOuterRadius=0.;
150}
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
void SetMaxEnergy(const G4double anEnergy)
G4VPreCompoundModel * GetDeExcitation() const
void SetDeExcitation(G4VPreCompoundModel *ptr)
G4ExcitationHandler * GetExcitationHandler() const

◆ ~G4BinaryCascade()

G4BinaryCascade::~G4BinaryCascade ( )
virtual

Definition at line 159 of file G4BinaryCascade.cc.

160{
161 ClearAndDestroy(&theTargetList);
162 ClearAndDestroy(&theSecondaryList);
163 ClearAndDestroy(&theCapturedList);
164 delete thePropagator;
165 delete theCollisionMgr;
166 std::for_each(theImR.begin(), theImR.end(), Delete<G4BCAction>());
167 delete theLateParticle;
168 //delete theExcitationHandler;
169 delete theH1Scatterer;
170}

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 208 of file G4BinaryCascade.cc.

211{
212 static G4int eventcounter=0;
213
214 // if ( eventcounter == 0 ) {
215 // SetEpReportLevel(3); // report non conservation with model etc.
216 // G4double relativeLevel = 1*perCent;
217 // G4double absoluteLevel = 2*MeV;
218 // SetEnergyMomentumCheckLevels(relativeLevel,absoluteLevel);
219 // }
220
221 //if(eventcounter == 100*(eventcounter/100) )
222 eventcounter++;
223 if(getenv("BCDEBUG") ) G4cerr << " ######### Binary Cascade Reaction number starts ######### "<<eventcounter<<G4endl;
224
225 G4LorentzVector initial4Momentum = aTrack.Get4Momentum();
226 G4ParticleDefinition * definition = const_cast<G4ParticleDefinition *>(aTrack.GetDefinition());
227
228 if(initial4Momentum.e()-initial4Momentum.m()<theBCminP &&
229 ( definition==G4Neutron::NeutronDefinition() || definition==G4Proton::ProtonDefinition() ) )
230 {
231 return theDeExcitation->ApplyYourself(aTrack, aNucleus);
232 }
233
235 // initialize the G4V3DNucleus from G4Nucleus
237
238 // Build a KineticTrackVector with the G4Track
239 G4KineticTrackVector * secondaries;// = new G4KineticTrackVector;
240 G4ThreeVector initialPosition(0., 0., 0.); // will be set later
241
242 if(!getenv("I_Am_G4BinaryCascade_Developer") )
243 {
244 if(definition!=G4Neutron::NeutronDefinition() &&
245 definition!=G4Proton::ProtonDefinition() &&
246 definition!=G4PionPlus::PionPlusDefinition() &&
248 {
249 G4cerr << "You are using G4BinaryCascade for projectiles other than nucleons or pions."<<G4endl;
250 G4cerr << "If you want to continue, please switch on the developer environment: "<<G4endl;
251 G4cerr << "setenv I_Am_G4BinaryCascade_Developer 1 "<<G4endl<<G4endl;
252 throw G4HadronicException(__FILE__, __LINE__, "G4BinaryCascade - used for unvalid particle type - Fatal");
253 }
254 }
255
256 // keep primary
257 thePrimaryType = definition;
258 thePrimaryEscape = false;
259
260 // try until an interaction will happen
261 G4ReactionProductVector * products = 0;
262 G4int interactionCounter = 0;
263 do
264 {
265 // reset status that could be changed in previous loop event
266 theCollisionMgr->ClearAndDestroy();
267
268 if(products != 0)
269 { // free memory from previous loop event
270 ClearAndDestroy(products);
271 delete products;
272 products=0;
273 }
274
275 G4int massNumber=aNucleus.GetA_asInt();
276 the3DNucleus->Init(massNumber, aNucleus.GetZ_asInt());
277 thePropagator->Init(the3DNucleus);
278 // GF Leak on kt??? but where to delete?
279 G4KineticTrack * kt;// = new G4KineticTrack(definition, 0., initialPosition, initial4Momentum);
280 do // sample impact parameter until collisions are found
281 {
282 theCurrentTime=0;
283 G4double radius = the3DNucleus->GetOuterRadius()+3*fermi;
284 initialPosition=GetSpherePoint(1.1*radius, initial4Momentum); // get random position
285 kt = new G4KineticTrack(definition, 0., initialPosition, initial4Momentum);
287 // secondaries has been cleared by Propagate() in the previous loop event
288 secondaries= new G4KineticTrackVector;
289 secondaries->push_back(kt);
290 if(massNumber > 1) // 1H1 is special case
291 {
292 products = Propagate(secondaries, the3DNucleus);
293 } else {
294 products = Propagate1H1(secondaries,the3DNucleus);
295 }
296 } while(! products ); // until we FIND a collision...
297
298 if(++interactionCounter>99) break;
299 } while(products->size() == 0); // ...until we find an ALLOWED collision
300
301 if(products->size()>0)
302 {
303 // G4cout << "BIC Applyyourself: number of products " << products->size() << G4endl;
304
305 // Fill the G4ParticleChange * with products
307 G4ReactionProductVector::iterator iter;
308
309 for(iter = products->begin(); iter != products->end(); ++iter)
310 {
311 G4DynamicParticle * aNew =
312 new G4DynamicParticle((*iter)->GetDefinition(),
313 (*iter)->GetTotalEnergy(),
314 (*iter)->GetMomentum());
316 }
317
318 // DebugEpConservation(track, products);
319
320
321 } else { // no interaction, return primary
322 if(getenv("BCDEBUG") ) G4cerr << " ######### Binary Cascade Reaction number void ######### "<<eventcounter<<G4endl;
326 }
327
328 ClearAndDestroy(products);
329 delete products;
330
331 delete the3DNucleus;
332 the3DNucleus = NULL;
333
334 if(getenv("BCDEBUG") ) G4cerr << " ######### Binary Cascade Reaction number ends ######### "<<eventcounter<<G4endl;
335
336 return &theParticleChange;
337}
@ isAlive
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
Hep3Vector unit() const
Hep3Vector vect() const
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *, G4V3DNucleus *)
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
CascadeState SetState(const CascadeState new_state)
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
virtual G4double GetOuterRadius()=0
virtual void Init(G4int theA, G4int theZ)=0
virtual void Init(G4V3DNucleus *theNucleus)=0
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)=0

◆ ModelDescription()

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

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 172 of file G4BinaryCascade.cc.

173{
174 outFile << "G4BinaryCascade is an intra-nuclear cascade model in which\n"
175 << "an incident hadron collides with a nucleon, forming two\n"
176 << "final-state particles, one or both of which may be resonances.\n"
177 << "The resonances then decay hadronically and the decay products\n"
178 << "are then propagated through the nuclear potential along curved\n"
179 << "trajectories until they re-interact or leave the nucleus.\n"
180 << "This model is valid for incident pions up to 1.5 GeV and\n"
181 << "nucleons up to 10 GeV.\n";
182}

◆ Propagate()

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

Implements G4VIntraNuclearTransportModel.

Definition at line 339 of file G4BinaryCascade.cc.

342{
343 G4ping debug("debug_G4BinaryCascade");
344#ifdef debug_BIC_Propagate
345 G4cout << "G4BinaryCascade Propagate starting -------------------------------------------------------" <<G4endl;
346#endif
347
348 the3DNucleus=aNucleus;
350 theOuterRadius = the3DNucleus->GetOuterRadius();
351 theCurrentTime=0;
352 theProjectile4Momentum=G4LorentzVector(0,0,0,0);
353 // build theSecondaryList, theProjectileList and theCapturedList
354 ClearAndDestroy(&theCapturedList);
355 ClearAndDestroy(&theSecondaryList);
356 theSecondaryList.clear();
357 ClearAndDestroy(&theFinalState);
358 std::vector<G4KineticTrack *>::iterator iter;
359 theCollisionMgr->ClearAndDestroy();
360
361 theCutOnP=90*MeV;
362 if(the3DNucleus->GetMass()>30) theCutOnP = 70*MeV;
363 if(the3DNucleus->GetMass()>60) theCutOnP = 50*MeV;
364 if(the3DNucleus->GetMass()>120) theCutOnP = 45*MeV;
365
366
367 BuildTargetList();
368
369#ifdef debug_BIC_GetExcitationEnergy
370 G4cout << "ExcitationEnergy0 " << GetExcitationEnergy() << G4endl;
371#endif
372
373 thePropagator->Init(the3DNucleus);
374
375 G4bool success = BuildLateParticleCollisions(secondaries);
376 if (! success ) // fails if no excitation energy left....
377 {
378 products=HighEnergyModelFSProducts(products, secondaries);
379 ClearAndDestroy(secondaries);
380 delete secondaries;
381
382#ifdef debug_G4BinaryCascade
383 G4cout << "G4BinaryCascade::Propagate: warning - high energy model failed energy conservation, returning unchanged high energy final state" << G4endl;
384#endif
385
386 return products;
387 }
388 // check baryon and charge ...
389
390 _CheckChargeAndBaryonNumber_("lateparticles");
391
392 // if called stand alone find first collisions
393 FindCollisions(&theSecondaryList);
394
395
396 if(theCollisionMgr->Entries() == 0 ) //late particles ALWAYS create Entries
397 {
398 //G4cout << " no collsions -> return 0" << G4endl;
399 delete products;
400#ifdef debug_BIC_return
401 G4cout << "return @ begin2, no collisions "<< G4endl;
402#endif
403 return 0;
404 }
405
406 // end of initialization: do the job now
407 // loop until there are no more collisions
408
409
410 G4bool haveProducts = false;
411 G4int collisionCount=0;
412 theMomentumTransfer=G4ThreeVector(0,0,0);
413 while(theCollisionMgr->Entries() > 0 && currentZ)
414 {
415 if(Absorb()) { // absorb secondaries, pions only
416 haveProducts = true;
417 }
419 if(Capture()) { // capture secondaries, nucleons only
420 haveProducts = true;
421 }
423
424 // propagate to the next collision if any (collisions could have been deleted
425 // by previous absorption or capture)
426 if(theCollisionMgr->Entries() > 0)
427 {
429 nextCollision = theCollisionMgr->GetNextCollision();
430#ifdef debug_BIC_Propagate_Collisions
431 G4cout << " NextCollision * , Time, curtime = " << nextCollision << " "
432 <<nextCollision->GetCollisionTime()<< " " <<
433 theCurrentTime<< G4endl;
434#endif
435 if (!DoTimeStep(nextCollision->GetCollisionTime()-theCurrentTime) )
436 {
437 // Check if nextCollision is still valid, ie. particle did not leave nucleus
438 if (theCollisionMgr->GetNextCollision() != nextCollision )
439 {
440 nextCollision = 0;
441 }
442 }
443
444 if( nextCollision )
445 {
446 if (ApplyCollision(nextCollision))
447 {
448 //G4cerr << "ApplyCollision success " << G4endl;
449 haveProducts = true;
450 collisionCount++;
451 _CheckChargeAndBaryonNumber_("ApplyCollision");
452
453 } else {
454 //G4cerr << "ApplyCollision failure " << G4endl;
455 theCollisionMgr->RemoveCollision(nextCollision);
456 }
457 }
458 }
459 }
460
461 //--------- end of while on Collisions
462 //G4cout << "currentZ @ end loop " << currentZ << G4endl;
463 if ( ! currentZ ){
464 // nucleus completely destroyed, fill in ReactionProductVector
465 products = FillVoidNucleusProducts(products);
466#ifdef debug_BIC_return
467 G4cout << "return @ Z=0 after collision loop "<< G4endl;
468 PrintKTVector(&theSecondaryList,std::string(" theSecondaryList"));
469 G4cout << "theTargetList size: " << theTargetList.size() << G4endl;
470 PrintKTVector(&theTargetList,std::string(" theTargetList"));
471 PrintKTVector(&theCapturedList,std::string(" theCapturedList"));
472
473 G4cout << " ExcitE be4 Correct : " <<GetExcitationEnergy() << G4endl;
474 G4cout << " Mom Transfered to nucleus : " << theMomentumTransfer << " " << theMomentumTransfer.mag() << G4endl;
475 PrintKTVector(&theFinalState,std::string(" FinalState uncorrected"));
476 G4cout << "returned products: " << products->size() << G4endl;
477#endif
478 // @@GF FixMe cannot return here.
479 return products;
480 }
481
482 // No more collisions: absorb, capture and propagate the secondaries out of the nucleus
483 if(Absorb()) {
484 haveProducts = true;
485 // G4cout << "Absorb sucess " << G4endl;
486 }
487
488 if(Capture()) {
489 haveProducts = true;
490 // G4cout << "Capture sucess " << G4endl;
491 }
492
493 if(!haveProducts) // no collisions happened. Return an empty vector.
494 {
495#ifdef debug_BIC_return
496 G4cout << "return 3, no products "<< G4endl;
497#endif
498 return products;
499 }
500
501
502#ifdef debug_BIC_Propagate
503 G4cout << " Momentum transfer to Nucleus " << theMomentumTransfer << " " << theMomentumTransfer.mag() << G4endl;
504 G4cout << " Stepping particles out...... " << G4endl;
505#endif
506
507 StepParticlesOut();
508
509
510 if ( theSecondaryList.size() > 0 )
511 {
512#ifdef debug_G4BinaryCascade
513 G4cerr << "G4BinaryCascade: Warning, have active particles at end" << G4endl;
514 PrintKTVector(&theSecondaryList, "active particles @ end added to theFinalState");
515#endif
516 // add left secondaries to FinalSate
517 for ( iter =theSecondaryList.begin(); iter != theSecondaryList.end(); ++iter)
518 {
519 theFinalState.push_back(*iter);
520 }
521 theSecondaryList.clear();
522
523 }
524 while ( theCollisionMgr->Entries() > 0 )
525 {
526#ifdef debug_G4BinaryCascade
527 G4cerr << " Warning: remove left over collision(s) " << G4endl;
528#endif
529 theCollisionMgr->RemoveCollision(theCollisionMgr->GetNextCollision());
530 }
531
532#ifdef debug_BIC_Propagate_Excitation
533
534 PrintKTVector(&theSecondaryList,std::string(" theSecondaryList"));
535 G4cout << "theTargetList size: " << theTargetList.size() << G4endl;
536 // PrintKTVector(&theTargetList,std::string(" theTargetList"));
537 PrintKTVector(&theCapturedList,std::string(" theCapturedList"));
538
539 G4cout << " ExcitE be4 Correct : " <<GetExcitationEnergy() << G4endl;
540 G4cout << " Mom Transfered to nucleus : " << theMomentumTransfer << " " << theMomentumTransfer.mag() << G4endl;
541 PrintKTVector(&theFinalState,std::string(" FinalState uncorrected"));
542#endif
543
544 //
545
546
547 G4double ExcitationEnergy=GetExcitationEnergy();
548
549#ifdef debug_BIC_Propagate_finals
550 PrintKTVector(&theFinalState,std::string(" FinalState be4 corr"));
551 G4cout << " Excitation Energy prefinal, #collisions:, out, captured "
552 << ExcitationEnergy << " "
553 << collisionCount << " "
554 << theFinalState.size() << " "
555 << theCapturedList.size()<<G4endl;
556#endif
557
558 if (ExcitationEnergy < 0 )
559 {
560 G4int maxtry=5, ntry=0;
561 do {
562 CorrectFinalPandE();
563 ExcitationEnergy=GetExcitationEnergy();
564 } while ( ++ntry < maxtry && ExcitationEnergy < 0 );
565 }
566
567#ifdef debug_BIC_Propagate_finals
568 PrintKTVector(&theFinalState,std::string(" FinalState corrected"));
569 G4cout << " Excitation Energy final, #collisions:, out, captured "
570 << ExcitationEnergy << " "
571 << collisionCount << " "
572 << theFinalState.size() << " "
573 << theCapturedList.size()<<G4endl;
574#endif
575
576
577 if ( ExcitationEnergy < 0. )
578 {
579 // if ( ExcitationEnergy < 0. )
580 {
581 //#ifdef debug_G4BinaryCascade
582 // G4cerr << "G4BinaryCascade-Warning: negative excitation energy ";
583 // G4cerr <<ExcitationEnergy<<G4endl;
584 // PrintKTVector(&theFinalState,std::string("FinalState"));
585 // PrintKTVector(&theCapturedList,std::string("captured"));
586 // G4cout << "negative ExE:Final 4Momentum .mag: " << GetFinal4Momentum()
587 // << " "<< GetFinal4Momentum().mag()<< G4endl
588 // << "negative ExE:FinalNucleusMom .mag: " << GetFinalNucleusMomentum()
589 // << " "<< GetFinalNucleusMomentum().mag()<< G4endl;
590 //#endif
591 }
592 ClearAndDestroy(products);
593 //G4cout << " negative Excitation E return empty products " << products << G4endl;
594#ifdef debug_BIC_return
595 G4cout << "return 4, excit < 0 "<< G4endl;
596#endif
597 return products; // return empty products- FIXME
598 }
599
600 G4ReactionProductVector * precompoundProducts=DeExcite();
601
602
603 G4DecayKineticTracks decay(&theFinalState);
604
605 products= ProductsAddFinalState(products, theFinalState);
606
607 products= ProductsAddPrecompound(products, precompoundProducts);
608
609// products=ProductsAddFakeGamma(products);
610
611
612 thePrimaryEscape = true;
613
614 #ifdef debug_BIC_return
615 G4cout << "return @end, all ok "<< G4endl;
616 //G4cout << " return products " << products << G4endl;
617 #endif
618
619 return products;
620}
#define _CheckChargeAndBaryonNumber_(val)
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition: G4Types.hh:67
G4DLLIMPORT std::ostream G4cout
double mag() const
void RemoveCollision(G4CollisionInitialState *collision)
G4CollisionInitialState * GetNextCollision()
virtual G4double GetMass()=0
Definition: G4ping.hh:34

Referenced by ApplyYourself().

◆ PropagateModelDescription()

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

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 183 of file G4BinaryCascade.cc.

184{
185 outFile << "G4BinaryCascade propagtes secondaries produced by a high\n"
186 << "energy model through the wounded nucleus.\n"
187 << "Secondaries are followed after the formation time and if\n"
188 << "within the nucleus are propagated through the nuclear\n"
189 << "potential along curved trajectories until they interact\n"
190 << "with a nucleon, decay, or leave the nucleus.\n"
191 << "An interaction of a secondary with a nucleon produces two\n"
192 << "final-state particles, one or both of which may be resonances.\n"
193 << "Resonances decay hadronically and the decay products\n"
194 << "are in turn propagated through the nuclear potential along curved\n"
195 << "trajectories until they re-interact or leave the nucleus.\n"
196 << "This model is valid for pions up to 1.5 GeV and\n"
197 << "nucleons up to about 3.5 GeV.\n";
198}

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