Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IntraNucleiCascader Class Reference

#include <G4IntraNucleiCascader.hh>

+ Inheritance diagram for G4IntraNucleiCascader:

Public Member Functions

 G4IntraNucleiCascader ()
 
virtual ~G4IntraNucleiCascader ()
 
void collide (G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
 
void rescatter (G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
 
void setVerboseLevel (G4int verbose=0)
 
- Public Member Functions inherited from G4CascadeColliderBase
 G4CascadeColliderBase (const char *name, G4int verbose=0)
 
virtual ~G4CascadeColliderBase ()
 
virtual void rescatter (G4InuclParticle *, G4KineticTrackVector *, G4V3DNucleus *, G4CollisionOutput &)
 
virtual void setVerboseLevel (G4int verbose=0)
 
virtual void setConservationChecks (G4bool doBalance=true)
 
- Public Member Functions inherited from G4VCascadeCollider
 G4VCascadeCollider (const char *name, G4int verbose=0)
 
virtual ~G4VCascadeCollider ()
 
virtual void collide (G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)=0
 
virtual void setVerboseLevel (G4int verbose=0)
 

Protected Member Functions

G4bool initialize (G4InuclParticle *bullet, G4InuclParticle *target)
 
void newCascade (G4int itry)
 
void setupCascade ()
 
void generateCascade ()
 
G4bool finishCascade ()
 
void finalize (G4int itry, G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
 
G4InuclParticlecreateTarget (G4V3DNucleus *theNucleus)
 
void preloadCascade (G4V3DNucleus *theNucleus, G4KineticTrackVector *theSecondaries)
 
void copyWoundedNucleus (G4V3DNucleus *theNucleus)
 
void copySecondaries (G4KineticTrackVector *theSecondaries)
 
void processSecondary (const G4KineticTrack *aSecondary)
 
void releaseSecondary (const G4KineticTrack *aSecondary)
 
void processTrappedParticle (const G4CascadParticle &trapped)
 
void decayTrappedParticle (const G4CascadParticle &trapped)
 
- Protected Member Functions inherited from G4CascadeColliderBase
virtual G4bool useEPCollider (G4InuclParticle *bullet, G4InuclParticle *target) const
 
virtual G4bool explosion (G4InuclNuclei *target) const
 
virtual G4bool explosion (G4Fragment *target) const
 
virtual G4bool explosion (G4int A, G4int Z, G4double excitation) const
 
virtual G4bool inelasticInteractionPossible (G4InuclParticle *bullet, G4InuclParticle *target, G4double ekin) const
 
virtual G4bool validateOutput (G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
 
virtual G4bool validateOutput (G4InuclParticle *bullet, G4InuclParticle *target, const std::vector< G4InuclElementaryParticle > &particles)
 
virtual G4bool validateOutput (G4InuclParticle *bullet, G4InuclParticle *target, const std::vector< G4InuclNuclei > &fragments)
 
- Protected Member Functions inherited from G4VCascadeCollider
virtual void setName (const char *name)
 

Additional Inherited Members

- Protected Attributes inherited from G4CascadeColliderBase
G4InteractionCase interCase
 
G4bool doConservationChecks
 
G4CascadeCheckBalancebalance
 
- Protected Attributes inherited from G4VCascadeCollider
const char * theName
 
G4int verboseLevel
 

Detailed Description

Definition at line 85 of file G4IntraNucleiCascader.hh.

Constructor & Destructor Documentation

◆ G4IntraNucleiCascader()

G4IntraNucleiCascader::G4IntraNucleiCascader ( )

Definition at line 156 of file G4IntraNucleiCascader.cc.

157 : G4CascadeColliderBase("G4IntraNucleiCascader"), model(new G4NucleiModel),
158 theElementaryParticleCollider(new G4ElementaryParticleCollider),
159 theRecoilMaker(new G4CascadeRecoilMaker), theClusterMaker(0),
160 tnuclei(0), bnuclei(0), bparticle(0),
161 minimum_recoil_A(0.), coulombBarrier(0.),
162 nucleusTarget(new G4InuclNuclei),
163 protonTarget(new G4InuclElementaryParticle) {
165 theClusterMaker = new G4CascadeCoalescence;
166}
static G4bool doCoalescence()

◆ ~G4IntraNucleiCascader()

G4IntraNucleiCascader::~G4IntraNucleiCascader ( )
virtual

Definition at line 168 of file G4IntraNucleiCascader.cc.

168 {
169 delete model;
170 delete theElementaryParticleCollider;
171 delete theRecoilMaker;
172 delete theClusterMaker;
173 delete nucleusTarget;
174 delete protonTarget;
175}

Member Function Documentation

◆ collide()

void G4IntraNucleiCascader::collide ( G4InuclParticle bullet,
G4InuclParticle target,
G4CollisionOutput globalOutput 
)
virtual

Implements G4VCascadeCollider.

Definition at line 186 of file G4IntraNucleiCascader.cc.

188 {
189 if (verboseLevel) G4cout << " >>> G4IntraNucleiCascader::collide " << G4endl;
190
191 if (!initialize(bullet, target)) return; // Load buffers and drivers
192
193 G4int itry = 0;
194 do {
195 newCascade(++itry);
196 setupCascade();
198 } while (!finishCascade() && itry<itry_max);
199
200 finalize(itry, bullet, target, globalOutput);
201}
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4bool initialize(G4InuclParticle *bullet, G4InuclParticle *target)
void finalize(G4int itry, G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)

Referenced by G4InuclCollider::collide().

◆ copySecondaries()

void G4IntraNucleiCascader::copySecondaries ( G4KineticTrackVector theSecondaries)
protected

Definition at line 687 of file G4IntraNucleiCascader.cc.

687 {
688 if (verboseLevel > 1)
689 G4cout << " >>> G4IntraNucleiCascader::copySecondaries" << G4endl;
690
691 for (size_t i=0; i<secondaries->size(); i++) {
692 if (verboseLevel > 3) G4cout << " processing secondary " << i << G4endl;
693
694 processSecondary((*secondaries)[i]); // Copy to cascade or to output
695 }
696
697 // Sort list of secondaries to put leading particle first
698 std::sort(cascad_particles.begin(), cascad_particles.end(),
700
701 if (verboseLevel > 2) {
702 G4cout << " Original list of " << secondaries->size() << " secondaries"
703 << " produced " << cascad_particles.size() << " cascade, "
704 << output.numberOfOutgoingParticles() << " released particles, "
705 << output.numberOfOutgoingNuclei() << " fragments" << G4endl;
706 }
707}
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
void processSecondary(const G4KineticTrack *aSecondary)

Referenced by preloadCascade().

◆ copyWoundedNucleus()

void G4IntraNucleiCascader::copyWoundedNucleus ( G4V3DNucleus theNucleus)
protected

Definition at line 657 of file G4IntraNucleiCascader.cc.

657 {
658 if (verboseLevel > 1)
659 G4cout << " >>> G4IntraNucleiCascader::copyWoundedNucleus" << G4endl;
660
661 // Loop over nucleons and count hits as exciton holes
662 theExitonConfiguration.clear();
663 hitNucleons.clear();
664 if (theNucleus->StartLoop()) {
665 G4Nucleon* nucl = 0;
666 G4int nuclType = 0;
667 while ((nucl = theNucleus->GetNextNucleon())) {
668 if (nucl->AreYouHit()) { // Found previously interacted nucleon
670 theExitonConfiguration.incrementHoles(nuclType);
671 hitNucleons.push_back(nucl->GetPosition());
672 }
673 }
674 }
675
676 if (verboseLevel > 3)
677 G4cout << " nucleus has " << theExitonConfiguration.neutronHoles
678 << " neutrons hit, " << theExitonConfiguration.protonHoles
679 << " protons hit" << G4endl;
680
681 // Preload nuclear model with confirmed hits, including locations
682 model->reset(theExitonConfiguration.neutronHoles,
683 theExitonConfiguration.protonHoles, &hitNucleons);
684}
void reset(G4int nHitNeutrons=0, G4int nHitProtons=0, const std::vector< G4ThreeVector > *hitPoints=0)
G4bool AreYouHit() const
Definition: G4Nucleon.hh:97
G4ParticleDefinition * GetParticleType() const
Definition: G4Nucleon.hh:84
virtual const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:68
virtual G4Nucleon * GetNextNucleon()=0
virtual G4bool StartLoop()=0

Referenced by preloadCascade().

◆ createTarget()

G4InuclParticle * G4IntraNucleiCascader::createTarget ( G4V3DNucleus theNucleus)
protected

Definition at line 628 of file G4IntraNucleiCascader.cc.

628 {
629 G4int theNucleusA = theNucleus->GetMassNumber();
630 G4int theNucleusZ = theNucleus->GetCharge();
631
632 if (theNucleusA > 1) {
633 if (!nucleusTarget) nucleusTarget = new G4InuclNuclei; // Just in case
634 nucleusTarget->fill(0., theNucleusA, theNucleusZ, 0.);
635 return nucleusTarget;
636 } else {
637 if (!protonTarget) protonTarget = new G4InuclElementaryParticle;
638 protonTarget->fill(0., (theNucleusZ==1)?proton:neutron);
639 return protonTarget;
640 }
641
642 return 0; // Can never actually get here
643}
@ neutron
void fill(G4int ityp, Model model=DefaultModel)
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
virtual G4int GetCharge()=0
virtual G4int GetMassNumber()=0

Referenced by rescatter().

◆ decayTrappedParticle()

void G4IntraNucleiCascader::decayTrappedParticle ( const G4CascadParticle trapped)
protected

Definition at line 816 of file G4IntraNucleiCascader.cc.

817 {
818 if (verboseLevel > 3)
819 G4cout << " unstable must be decayed in flight" << G4endl;
820
821 const G4InuclElementaryParticle& trappedP = trapped.getParticle();
822
823 G4DecayTable* unstable = trappedP.getDefinition()->GetDecayTable();
824 if (!unstable) { // No decay table; cannot decay!
825 if (verboseLevel > 3)
826 G4cerr << " no decay table! Releasing trapped particle" << G4endl;
827
828 output.addOutgoingParticle(trappedP);
829 return;
830 }
831
832 // Get secondaries from decay in particle's rest frame
833 G4DecayProducts* daughters = unstable->SelectADecayChannel()->DecayIt( trappedP.getDefinition()->GetPDGMass() );
834 if (!daughters) { // No final state; cannot decay!
835 if (verboseLevel > 3)
836 G4cerr << " no daughters! Releasing trapped particle" << G4endl;
837
838 output.addOutgoingParticle(trappedP);
839 return;
840 }
841
842 if (verboseLevel > 3)
843 G4cout << " " << daughters->entries() << " decay daughters" << G4endl;
844
845 // Convert secondaries to lab frame
846 G4double decayEnergy = trappedP.getEnergy();
847 G4ThreeVector decayDir = trappedP.getMomentum().vect().unit();
848 daughters->Boost(decayEnergy, decayDir);
849
850 // Put all the secondaries onto the list for propagation
851 const G4ThreeVector& decayPos = trapped.getPosition();
852 G4int zone = trapped.getCurrentZone();
853 G4int gen = trapped.getGeneration()+1;
854
855 for (G4int i=0; i<daughters->entries(); i++) {
856 G4DynamicParticle* idaug = (*daughters)[i];
857
859
860 // Propagate hadronic secondaries with known interactions (tables)
861 if (G4CascadeChannelTables::GetTable(idaugEP.type())) {
862 if (verboseLevel > 3) G4cout << " propagating " << idaugEP << G4endl;
863 cascad_particles.push_back(G4CascadParticle(idaugEP,decayPos,zone,0.,gen));
864 } else {
865 if (verboseLevel > 3) G4cout << " releasing " << idaugEP << G4endl;
866 output.addOutgoingParticle(idaugEP);
867 }
868 }
869}
double G4double
Definition: G4Types.hh:64
G4DLLIMPORT std::ostream G4cerr
Hep3Vector unit() const
Hep3Vector vect() const
G4int getGeneration() const
const G4InuclElementaryParticle & getParticle() const
G4int getCurrentZone() const
const G4ThreeVector & getPosition() const
static const G4CascadeChannel * GetTable(G4int initialState)
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
G4int entries() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * SelectADecayChannel()
Definition: G4DecayTable.cc:81
G4ParticleDefinition * getDefinition() const
G4LorentzVector getMomentum() const
G4double getEnergy() const
G4DecayTable * GetDecayTable() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0

Referenced by processTrappedParticle().

◆ finalize()

void G4IntraNucleiCascader::finalize ( G4int  itry,
G4InuclParticle bullet,
G4InuclParticle target,
G4CollisionOutput globalOutput 
)
protected

Definition at line 606 of file G4IntraNucleiCascader.cc.

608 {
609 if (itry >= itry_max) {
610 if (verboseLevel) {
611 G4cout << " IntraNucleiCascader-> no inelastic interaction after "
612 << itry << " attempts " << G4endl;
613 }
614
615 output.trivialise(bullet, target);
616 } else if (verboseLevel) {
617 G4cout << " IntraNucleiCascader output after trials " << itry << G4endl;
618 }
619
620 // Copy final generated cascade to output buffer for return
621 globalOutput.add(output);
622}
void add(const G4CollisionOutput &right)
void trivialise(G4InuclParticle *bullet, G4InuclParticle *target)

Referenced by collide(), and rescatter().

◆ finishCascade()

G4bool G4IntraNucleiCascader::finishCascade ( )
protected

Definition at line 457 of file G4IntraNucleiCascader.cc.

457 {
458 if (verboseLevel > 1)
459 G4cout << " >>> G4IntraNucleiCascader::finishCascade ?" << G4endl;
460
461 // Add left-over cascade particles to output
462 output.addOutgoingParticles(cascad_particles);
463 cascad_particles.clear();
464
465 // Cascade is finished. Check if it's OK.
466 if (verboseLevel>3) {
467 G4cout << " G4IntraNucleiCascader finished" << G4endl;
468 output.printCollisionOutput();
469 }
470
471 // Apply cluster coalesence model to produce light ions
472 if (theClusterMaker) {
473 theClusterMaker->setVerboseLevel(verboseLevel);
474 theClusterMaker->FindClusters(output);
475
476 // Update recoil fragment after generating light ions
477 if (verboseLevel>3) G4cout << " Recomputing recoil fragment" << G4endl;
478 theRecoilMaker->collide(interCase.getBullet(), interCase.getTarget(),
479 output);
480 if (verboseLevel>3) {
481 G4cout << " After cluster coalescence" << G4endl;
482 output.printCollisionOutput();
483 }
484 }
485
486 // Use last created recoil fragment instead of re-constructing
487 G4int afin = theRecoilMaker->getRecoilA();
488 G4int zfin = theRecoilMaker->getRecoilZ();
489
490 // FIXME: Should we deal with unbalanced (0,0) case before rejecting?
491
492 // Sanity check before proceeding
493 if (!theRecoilMaker->goodFragment() && !theRecoilMaker->wholeEvent()) {
494 if (verboseLevel > 1)
495 G4cerr << " Recoil nucleus is not physical: A=" << afin << " Z="
496 << zfin << G4endl;
497 return false; // Discard event and try again
498 }
499
500 const G4LorentzVector& presid = theRecoilMaker->getRecoilMomentum();
501
502 if (verboseLevel > 1) {
503 G4cout << " afin " << afin << " zfin " << zfin << G4endl;
504 }
505
506 if (afin == 0) return true; // Whole event fragmented, exit
507
508 if (afin == 1) { // Add bare nucleon to particle list
509 G4int last_type = (zfin==1) ? 1 : 2; // proton=1, neutron=2
510
512 G4double mres = presid.m();
513
514 // Check for sensible kinematics
515 if (mres-mass < -small_ekin) { // Insufficient recoil energy
516 if (verboseLevel > 2) G4cerr << " unphysical recoil nucleon" << G4endl;
517 return false;
518 }
519
520 if (mres-mass > small_ekin) { // Too much extra energy
521 if (verboseLevel > 2)
522 G4cerr << " extra energy with recoil nucleon" << G4endl;
523
524 // FIXME: For now, we add the nucleon as unbalanced, and let
525 // "SetOnShell" fudge things. This should be abandoned.
526 }
527
528 G4InuclElementaryParticle last_particle(presid, last_type,
530
531 if (verboseLevel > 3) {
532 G4cout << " adding recoiling nucleon to output list\n"
533 << last_particle << G4endl;
534 }
535
536 output.addOutgoingParticle(last_particle);
537
538 // Update recoil to include residual nucleon
539 theRecoilMaker->collide(interCase.getBullet(), interCase.getTarget(),
540 output);
541 }
542
543 // Process recoil fragment for consistency, exit or reject
544 if (output.numberOfOutgoingParticles() == 1) {
545 G4double Eex = theRecoilMaker->getRecoilExcitation();
546 if (std::abs(Eex) < quasielast_cut) {
547 if (verboseLevel > 3) {
548 G4cout << " quasi-elastic scatter with " << Eex << " MeV recoil"
549 << G4endl;
550 }
551
552 theRecoilMaker->setRecoilExcitation(Eex=0.);
553 if (verboseLevel > 3) {
554 G4cout << " Eex reset to " << theRecoilMaker->getRecoilExcitation()
555 << G4endl;
556 }
557 }
558 }
559
560 if (theRecoilMaker->goodNucleus()) {
561 theRecoilMaker->addExcitonConfiguration(theExitonConfiguration);
562
563 G4Fragment* recoilFrag = theRecoilMaker->makeRecoilFragment();
564 if (!recoilFrag) {
565 G4cerr << "Got null pointer for recoil fragment!" << G4endl;
566 return false;
567 }
568
569 if (verboseLevel > 2)
570 G4cout << " adding recoil fragment to output list" << G4endl;
571
572 output.addRecoilFragment(*recoilFrag);
573 }
574
575 // Put final-state particles in "leading order" for return
576 std::vector<G4InuclElementaryParticle>& opart = output.getOutgoingParticles();
577 std::sort(opart.begin(), opart.end(), G4ParticleLargerEkin());
578
579 // Adjust final state to balance momentum and energy if necessary
580 if (theRecoilMaker->wholeEvent() || theRecoilMaker->goodNucleus()) {
583 output.setVerboseLevel(0);
584
585 if (output.acceptable()) return true;
586 else if (verboseLevel>2) G4cerr << " Cascade setOnShell failed." << G4endl;
587 }
588
589 // Cascade not physically reasonable
590 if (afin <= minimum_recoil_A && minimum_recoil_A < tnuclei->getA()) {
591 ++minimum_recoil_A;
592 if (verboseLevel > 3) {
593 G4cout << " minimum recoil fragment increased to A " << minimum_recoil_A
594 << G4endl;
595 }
596 }
597
598 if (verboseLevel>2) G4cerr << " Cascade failed. Retrying..." << G4endl;
599 return false;
600}
void FindClusters(G4CollisionOutput &finalState)
void setVerboseLevel(G4int verbose)
void addExcitonConfiguration(const G4ExitonConfiguration exciton)
const G4LorentzVector & getRecoilMomentum() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4double getRecoilExcitation() const
void setRecoilExcitation(G4double Eexc)
G4Fragment * makeRecoilFragment()
void addRecoilFragment(const G4Fragment *aFragment)
void printCollisionOutput(std::ostream &os=G4cout) const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void setOnShell(G4InuclParticle *bullet, G4InuclParticle *target)
void setVerboseLevel(G4int verbose)
G4bool acceptable() const
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
G4InuclParticle * getBullet() const
G4InuclParticle * getTarget() const
static G4double getParticleMass(G4int type)

Referenced by collide(), and rescatter().

◆ generateCascade()

void G4IntraNucleiCascader::generateCascade ( )
protected

Definition at line 336 of file G4IntraNucleiCascader.cc.

336 {
337 if (verboseLevel>1) G4cout << " generateCascade " << G4endl;
338
339 G4int iloop = 0;
340 while (!cascad_particles.empty() && !model->empty()) {
341 iloop++;
342
343 if (verboseLevel > 2) {
344 G4cout << " Iteration " << iloop << ": Number of cparticles "
345 << cascad_particles.size() << " last one: \n"
346 << cascad_particles.back() << G4endl;
347 }
348
349 model->generateParticleFate(cascad_particles.back(),
350 theElementaryParticleCollider,
351 new_cascad_particles);
352
353 if (verboseLevel > 2) {
354 G4cout << " After generate fate: New particles "
355 << new_cascad_particles.size() << G4endl
356 << " Discarding last cparticle from list " << G4endl;
357 }
358
359 cascad_particles.pop_back();
360
361 // handle the result of a new step
362
363 if (new_cascad_particles.size() == 1) { // last particle goes without interaction
364 const G4CascadParticle& currentCParticle = new_cascad_particles[0];
365
366 if (model->stillInside(currentCParticle)) {
367 if (verboseLevel > 3)
368 G4cout << " particle still inside nucleus " << G4endl;
369
370 if (currentCParticle.getNumberOfReflections() < reflection_cut &&
371 model->worthToPropagate(currentCParticle)) {
372 if (verboseLevel > 3) G4cout << " continue reflections " << G4endl;
373 cascad_particles.push_back(currentCParticle);
374 } else {
375 processTrappedParticle(currentCParticle);
376 } // reflection or exciton
377
378 } else { // particle about to leave nucleus - check for Coulomb barrier
379 if (verboseLevel > 3) G4cout << " possible escape " << G4endl;
380
381 const G4InuclElementaryParticle& currentParticle =
382 currentCParticle.getParticle();
383
384 G4double KE = currentParticle.getKineticEnergy();
385 G4double mass = currentParticle.getMass();
386 G4double Q = currentParticle.getCharge();
387
388 if (verboseLevel > 3)
389 G4cout << " KE " << KE << " barrier " << Q*coulombBarrier << G4endl;
390
391 if (KE < Q*coulombBarrier) {
392 // Calculate barrier penetration
393 G4double CBP = 0.0;
394
395 // if (KE > 0.0001) CBP = std::exp(-0.00126*tnuclei->getZ()*0.25*
396 // (1./KE - 1./coulombBarrier));
397 if (KE > 0.0001) CBP = std::exp(-0.0181*0.5*tnuclei->getZ()*
398 (1./KE - 1./coulombBarrier)*
399 std::sqrt(mass*(coulombBarrier-KE)) );
400
401 if (G4UniformRand() < CBP) {
402 if (verboseLevel > 3)
403 G4cout << " tunneled\n" << currentParticle << G4endl;
404
405 // Tunnelling through barrier leaves KE unchanged
406 output.addOutgoingParticle(currentParticle);
407 } else {
408 processTrappedParticle(currentCParticle);
409 }
410 } else {
411 output.addOutgoingParticle(currentParticle);
412
413 if (verboseLevel > 3)
414 G4cout << " Goes out\n" << output.getOutgoingParticles().back()
415 << G4endl;
416 }
417 }
418 } else { // interaction
419 if (verboseLevel > 3)
420 G4cout << " interacted, adding new to list " << G4endl;
421
422 cascad_particles.insert(cascad_particles.end(),
423 new_cascad_particles.begin(),
424 new_cascad_particles.end());
425
426 std::pair<G4int, G4int> holes = model->getTypesOfNucleonsInvolved();
427 if (verboseLevel > 3)
428 G4cout << " adding new exciton holes " << holes.first << ","
429 << holes.second << G4endl;
430
431 theExitonConfiguration.incrementHoles(holes.first);
432
433 if (holes.second > 0)
434 theExitonConfiguration.incrementHoles(holes.second);
435 } // if (new_cascad_particles ...
436
437 // Evaluate nuclear residue
438 theRecoilMaker->collide(interCase.getBullet(), interCase.getTarget(),
439 output, cascad_particles);
440
441 G4double aresid = theRecoilMaker->getRecoilA();
442 if (verboseLevel > 2) {
443 G4cout << " cparticles remaining " << cascad_particles.size()
444 << " nucleus (model) has "
445 << model->getNumberOfNeutrons() << " n, "
446 << model->getNumberOfProtons() << " p "
447 << " residual fragment A " << aresid << G4endl;
448 }
449
450 if (aresid <= minimum_recoil_A) return; // Must have minimum fragment
451 } // while cascade-list and model
452}
#define G4UniformRand()
Definition: Randomize.hh:53
G4int getNumberOfReflections() const
void processTrappedParticle(const G4CascadParticle &trapped)
G4int getZ() const
G4double getKineticEnergy() const
G4double getMass() const
G4double getCharge() const
G4int getNumberOfNeutrons() const
G4bool worthToPropagate(const G4CascadParticle &cparticle) const
G4int getNumberOfProtons() const
void generateParticleFate(G4CascadParticle &cparticle, G4ElementaryParticleCollider *theEPCollider, std::vector< G4CascadParticle > &cascade)
G4bool empty() const
G4bool stillInside(const G4CascadParticle &cparticle)
std::pair< G4int, G4int > getTypesOfNucleonsInvolved() const

Referenced by collide(), and rescatter().

◆ initialize()

G4bool G4IntraNucleiCascader::initialize ( G4InuclParticle bullet,
G4InuclParticle target 
)
protected

Definition at line 227 of file G4IntraNucleiCascader.cc.

228 {
229 if (verboseLevel>1)
230 G4cout << " >>> G4IntraNucleiCascader::initialize " << G4endl;
231
232 // Configure processing modules
233 theRecoilMaker->setTolerance(small_ekin);
234
235 interCase.set(bullet,target); // Classify collision type
236
237 if (verboseLevel > 3) {
239 << *interCase.getTarget() << G4endl;
240 }
241
242 // Bullet may be nucleus or simple particle
243 bnuclei = dynamic_cast<G4InuclNuclei*>(interCase.getBullet());
244 bparticle = dynamic_cast<G4InuclElementaryParticle*>(interCase.getBullet());
245
246 if (!bnuclei && !bparticle) {
247 G4cerr << " G4IntraNucleiCascader: projectile is not a valid particle."
248 << G4endl;
249 return false;
250 }
251
252 // Target _must_ be nucleus
253 tnuclei = dynamic_cast<G4InuclNuclei*>(interCase.getTarget());
254 if (!tnuclei) {
255 if (verboseLevel)
256 G4cerr << " Target is not a nucleus. Abandoning." << G4endl;
257 return false;
258 }
259
260 model->generateModel(tnuclei);
261 coulombBarrier = 0.00126*tnuclei->getZ() / (1.+G4cbrt(tnuclei->getA()));
262
263 // Energy/momentum conservation usually requires a recoiling nuclear fragment
264 // This cut will be increased on each "itry" if momentum could not balance.
265 minimum_recoil_A = 0.;
266
267 if (verboseLevel > 3) {
268 G4LorentzVector momentum_in = bullet->getMomentum() + target->getMomentum();
269 G4cout << " intitial momentum E " << momentum_in.e() << " Px "
270 << momentum_in.x() << " Py " << momentum_in.y() << " Pz "
271 << momentum_in.z() << G4endl;
272 }
273
274 return true;
275}
void setTolerance(G4double tolerance)
void set(G4InuclParticle *part1, G4InuclParticle *part2)
G4int getA() const
void generateModel(G4InuclNuclei *nuclei)

Referenced by collide(), and rescatter().

◆ newCascade()

void G4IntraNucleiCascader::newCascade ( G4int  itry)
protected

Definition at line 279 of file G4IntraNucleiCascader.cc.

279 {
280 if (verboseLevel > 1) {
281 G4cout << " IntraNucleiCascader itry " << itry << " inter_case "
282 << interCase.code() << G4endl;
283 }
284
285 model->reset(); // Start new cascade process
286 output.reset();
287 new_cascad_particles.clear();
288 theExitonConfiguration.clear();
289
290 cascad_particles.clear(); // List of initial secondaries
291}

Referenced by collide(), and rescatter().

◆ preloadCascade()

void G4IntraNucleiCascader::preloadCascade ( G4V3DNucleus theNucleus,
G4KineticTrackVector theSecondaries 
)
protected

Definition at line 648 of file G4IntraNucleiCascader.cc.

649 {
650 if (verboseLevel > 1)
651 G4cout << " >>> G4IntraNucleiCascader::preloadCascade" << G4endl;
652
653 copyWoundedNucleus(theNucleus); // Update interacted nucleon counts
654 copySecondaries(theSecondaries); // Copy original to internal list
655}
void copyWoundedNucleus(G4V3DNucleus *theNucleus)
void copySecondaries(G4KineticTrackVector *theSecondaries)

Referenced by rescatter().

◆ processSecondary()

void G4IntraNucleiCascader::processSecondary ( const G4KineticTrack aSecondary)
protected

Definition at line 712 of file G4IntraNucleiCascader.cc.

712 {
713 if (!ktrack) return; // Sanity check
714
715 // Get particle type to determine whether to keep or release
716 G4ParticleDefinition* kpd = ktrack->GetDefinition();
717 if (!kpd) return;
718
720 if (!ktype) {
721 releaseSecondary(ktrack);
722 return;
723 }
724
725 if (verboseLevel > 1) {
726 G4cout << " >>> G4IntraNucleiCascader::processSecondary "
727 << kpd->GetParticleName() << G4endl;
728 }
729
730 // Allocate next local particle in buffer and fill
731 cascad_particles.resize(cascad_particles.size()+1); // Like push_back();
732 G4CascadParticle& cpart = cascad_particles.back();
733
734 // Convert momentum to Bertini internal units
735 cpart.getParticle().fill(ktrack->Get4Momentum()/GeV, ktype);
736 cpart.setGeneration(0);
737 cpart.setMovingInsideNuclei();
738 cpart.initializePath(0);
739
740 // Convert position units to Bertini's internal scale
741 G4ThreeVector cpos = ktrack->GetPosition()/model->getRadiusUnits();
742
743 cpart.updatePosition(cpos);
744 cpart.updateZone(model->getZone(cpos.mag()));
745
746 if (verboseLevel > 2)
747 G4cout << " Created cascade particle \n" << cpart << G4endl;
748}
double mag() const
void updateZone(G4int izone)
void setGeneration(G4int gen)
void updatePosition(const G4ThreeVector &pos)
void setMovingInsideNuclei(G4bool isMovingIn=true)
void initializePath(G4double npath)
void releaseSecondary(const G4KineticTrack *aSecondary)
G4int getZone(G4double r) const
G4double getRadiusUnits() const
const G4String & GetParticleName() const

Referenced by copySecondaries().

◆ processTrappedParticle()

void G4IntraNucleiCascader::processTrappedParticle ( const G4CascadParticle trapped)
protected

Definition at line 786 of file G4IntraNucleiCascader.cc.

787 {
788 const G4InuclElementaryParticle& trappedP = trapped.getParticle();
789
790 G4int xtype = trappedP.type();
791 if (verboseLevel > 3) G4cout << " exciton of type " << xtype << G4endl;
792
793 if (trappedP.nucleon()) { // normal exciton (proton or neutron)
794 theExitonConfiguration.incrementQP(xtype);
795 return;
796 }
797
798 if (trappedP.hyperon()) { // Not nucleon, so must be hyperon
799 decayTrappedParticle(trapped);
800 return;
801 }
802
803 // non-standard exciton; release it
804 // FIXME: this is a meson, so need to absorb it
805 if (verboseLevel > 3) {
806 G4cout << " non-standard should be absorbed, now released\n"
807 << trapped << G4endl;
808 }
809
810 output.addOutgoingParticle(trappedP);
811}
void decayTrappedParticle(const G4CascadParticle &trapped)

Referenced by generateCascade().

◆ releaseSecondary()

void G4IntraNucleiCascader::releaseSecondary ( const G4KineticTrack aSecondary)
protected

Definition at line 753 of file G4IntraNucleiCascader.cc.

753 {
754 G4ParticleDefinition* kpd = ktrack->GetDefinition();
755
756 if (verboseLevel > 1) {
757 G4cout << " >>> G4IntraNucleiCascader::releaseSecondary "
758 << kpd->GetParticleName() << G4endl;
759 }
760
761 // Convert light ion into nucleus on fragment list
762 if (dynamic_cast<G4Ions*>(kpd)) {
763 // Use resize() and fill() to avoid memory churn
764 output.getOutgoingNuclei().resize(output.numberOfOutgoingNuclei()+1);
765 G4InuclNuclei& inucl = output.getOutgoingNuclei().back();
766
767 inucl.fill(ktrack->Get4Momentum()/GeV,
768 kpd->GetAtomicMass(), kpd->GetAtomicNumber());
769 if (verboseLevel > 2)
770 G4cout << " Created pre-cascade fragment\n" << inucl << G4endl;
771 } else {
772 // Use resize() and fill() to avoid memory churn
773 output.getOutgoingParticles().resize(output.numberOfOutgoingParticles()+1);
774 G4InuclElementaryParticle& ipart = output.getOutgoingParticles().back();
775
776 // SPECIAL: Use G4PartDef directly, allowing unknown type code
777 ipart.fill(ktrack->Get4Momentum()/GeV, ktrack->GetDefinition());
778 if (verboseLevel > 2)
779 G4cout << " Created invalid pre-cascade particle\n" << ipart << G4endl;
780 }
781}
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
Definition: G4Ions.hh:52
G4int GetAtomicNumber() const
G4int GetAtomicMass() const

Referenced by processSecondary().

◆ rescatter()

void G4IntraNucleiCascader::rescatter ( G4InuclParticle bullet,
G4KineticTrackVector theSecondaries,
G4V3DNucleus theNucleus,
G4CollisionOutput globalOutput 
)
virtual

Reimplemented from G4CascadeColliderBase.

Definition at line 206 of file G4IntraNucleiCascader.cc.

209 {
210 if (verboseLevel)
211 G4cout << " >>> G4IntraNucleiCascader::rescatter " << G4endl;
212
213 G4InuclParticle* target = createTarget(theNucleus);
214 if (!initialize(bullet, target)) return; // Load buffers and drivers
215
216 G4int itry = 0;
217 do {
218 newCascade(++itry);
219 preloadCascade(theNucleus, theSecondaries);
221 } while (!finishCascade() && itry<itry_max);
222
223 finalize(itry, bullet, target, globalOutput);
224}
G4InuclParticle * createTarget(G4V3DNucleus *theNucleus)
void preloadCascade(G4V3DNucleus *theNucleus, G4KineticTrackVector *theSecondaries)

Referenced by G4InuclCollider::rescatter().

◆ setupCascade()

void G4IntraNucleiCascader::setupCascade ( )
protected

Definition at line 296 of file G4IntraNucleiCascader.cc.

296 {
297 if (verboseLevel > 1)
298 G4cout << " >>> G4IntraNucleiCascader::setupCascade" << G4endl;
299
300 if (interCase.hadNucleus()) { // particle with nuclei
301 if (verboseLevel > 3)
302 G4cout << " bparticle charge " << bparticle->getCharge()
303 << " baryon number " << bparticle->baryon() << G4endl;
304
305 cascad_particles.push_back(model->initializeCascad(bparticle));
306 } else { // nuclei with nuclei
307 G4int ab = bnuclei->getA();
308 G4int zb = bnuclei->getZ();
309
310 G4NucleiModel::modelLists all_particles; // Buffer to receive lists
311 model->initializeCascad(bnuclei, tnuclei, all_particles);
312
313 cascad_particles = all_particles.first;
314 output.addOutgoingParticles(all_particles.second);
315
316 if (cascad_particles.size() == 0) { // compound nuclei
317 G4int i;
318
319 for (i = 0; i < ab; i++) {
320 G4int knd = i < zb ? 1 : 2;
321 theExitonConfiguration.incrementQP(knd);
322 };
323
324 G4int ihn = G4int(2 * (ab-zb) * inuclRndm() + 0.5);
325 G4int ihz = G4int(2 * zb * inuclRndm() + 0.5);
326
327 for (i = 0; i < ihn; i++) theExitonConfiguration.incrementHoles(2);
328 for (i = 0; i < ihz; i++) theExitonConfiguration.incrementHoles(1);
329 }
330 } // if (interCase ...
331}
G4bool hadNucleus() const
G4CascadParticle initializeCascad(G4InuclElementaryParticle *particle)
std::pair< std::vector< G4CascadParticle >, std::vector< G4InuclElementaryParticle > > modelLists

Referenced by collide().

◆ setVerboseLevel()

void G4IntraNucleiCascader::setVerboseLevel ( G4int  verbose = 0)
virtual

Reimplemented from G4CascadeColliderBase.

Definition at line 177 of file G4IntraNucleiCascader.cc.

177 {
179 model->setVerboseLevel(verbose);
180 theElementaryParticleCollider->setVerboseLevel(verbose);
181 theRecoilMaker->setVerboseLevel(verbose);
182}
virtual void setVerboseLevel(G4int verbose=0)
void setVerboseLevel(G4int verbose)
virtual void setVerboseLevel(G4int verbose=0)

Referenced by G4InuclCollider::setVerboseLevel().


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