Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4BinaryCascade Class Reference

#include <G4BinaryCascade.hh>

Inheritance diagram for G4BinaryCascade:
G4VIntraNuclearTransportModel G4HadronicInteraction

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 G4ReactionProductVectorPropagateNuclNucl (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
 
void SetDeExcitation (G4VPreCompoundModel *ptr)
 
void Set3DNucleus (G4V3DNucleus *const value)
 
void SetPrimaryProjectile (const G4HadProjectile &aPrimary)
 
const G4StringGetModelName () const
 
- 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 &, 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
 
virtual G4int GetVerboseLevel () const
 
virtual 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, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 

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 ( G4VPreCompoundModel ptr = 0)

Definition at line 102 of file G4BinaryCascade.cc.

References G4ShortLivedConstructor::ConstructParticle(), G4HadronicInteractionRegistry::FindModel(), G4VIntraNuclearTransportModel::GetDeExcitation(), G4VPreCompoundModel::GetExcitationHandler(), python.hepunit::GeV, G4HadronicInteractionRegistry::Instance(), python.hepunit::MeV, python.hepunit::perCent, G4VIntraNuclearTransportModel::SetDeExcitation(), G4HadronicInteraction::SetEnergyMomentumCheckLevels(), G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

102  :
103 G4VIntraNuclearTransportModel("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 
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 }
G4VIntraNuclearTransportModel(const G4String &modelName="CascadeModel", G4VPreCompoundModel *ptr=0)
G4VPreCompoundModel * GetDeExcitation() const
const char * p
Definition: xmltok.h:285
void SetMinEnergy(G4double anEnergy)
G4ExcitationHandler * GetExcitationHandler() const
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMaxEnergy(const G4double anEnergy)
void SetDeExcitation(G4VPreCompoundModel *ptr)
float perCent
Definition: hepunit.py:239
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
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

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

Implements G4HadronicInteraction.

Definition at line 208 of file G4BinaryCascade.cc.

References G4HadFinalState::AddSecondary(), G4VPreCompoundModel::ApplyYourself(), G4HadFinalState::Clear(), G4CollisionManager::ClearAndDestroy(), CLHEP::HepLorentzVector::e(), python.hepunit::fermi, G4cerr, G4endl, G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4HadProjectile::GetKineticEnergy(), G4V3DNucleus::GetOuterRadius(), G4Nucleus::GetZ_asInt(), G4VFieldPropagation::Init(), G4V3DNucleus::Init(), isAlive, CLHEP::HepLorentzVector::m(), G4Neutron::NeutronDefinition(), G4KineticTrack::outside, G4PionMinus::PionMinusDefinition(), G4PionPlus::PionPlusDefinition(), Propagate(), G4Proton::ProtonDefinition(), G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), G4KineticTrack::SetState(), G4HadFinalState::SetStatusChange(), stopAndKill, G4VIntraNuclearTransportModel::the3DNucleus, G4VIntraNuclearTransportModel::theDeExcitation, G4HadronicInteraction::theParticleChange, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

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

Implements G4VIntraNuclearTransportModel.

Definition at line 329 of file G4BinaryCascade.cc.

References _CheckChargeAndBaryonNumber_, G4CollisionManager::ClearAndDestroy(), debug, G4INCL::ClusterDecay::decay(), G4CollisionManager::Entries(), G4cerr, G4cout, G4endl, G4CollisionInitialState::GetCollisionTime(), G4V3DNucleus::GetMass(), G4CollisionManager::GetNextCollision(), G4V3DNucleus::GetOuterRadius(), G4VFieldPropagation::Init(), CLHEP::Hep3Vector::mag(), python.hepunit::MeV, G4CollisionManager::RemoveCollision(), and G4VIntraNuclearTransportModel::the3DNucleus.

Referenced by ApplyYourself().

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

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