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

#include <G4CascadeInterface.hh>

+ Inheritance diagram for G4CascadeInterface:

Public Member Functions

 G4CascadeInterface (const G4String &name="BertiniCascade")
 
virtual ~G4CascadeInterface ()
 
G4ReactionProductVectorPropagate (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
void SetVerboseLevel (G4int verbose)
 
G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
G4bool IsApplicable (const G4ParticleDefinition *aPD) const
 
void useCascadeDeexcitation ()
 
void usePreCompoundDeexcitation ()
 
void useAblaDeexcitation ()
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void DumpConfiguration (std::ostream &outFile) const
 
- Public Member Functions inherited from G4VIntraNuclearTransportModel
 G4VIntraNuclearTransportModel (const G4String &mName="CascadeModel", G4VPreCompoundModel *ptr=nullptr)
 
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
 
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 G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
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 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
 

Static Public Member Functions

static void Initialize ()
 

Protected Member Functions

void clear ()
 
G4bool createBullet (const G4HadProjectile &aTrack)
 
G4bool createTarget (G4Nucleus &theNucleus)
 
G4bool createTarget (G4V3DNucleus *theNucleus)
 
G4bool createTarget (G4int A, G4int Z)
 
G4bool coulombBarrierViolation () const
 
G4bool retryInelasticProton () const
 
G4bool retryInelasticNucleus () const
 
void copyOutputToHadronicResult ()
 
G4ReactionProductVectorcopyOutputToReactionProducts ()
 
G4HadFinalStateNoInteraction (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
void checkFinalResult ()
 
void throwNonConservationFailure ()
 
G4DynamicParticlemakeDynamicParticle (const G4InuclElementaryParticle &iep) const
 
G4DynamicParticlemakeDynamicParticle (const G4InuclNuclei &inuc) const
 
- 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 ()
 

Additional Inherited Members

- 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 88 of file G4CascadeInterface.hh.

Constructor & Destructor Documentation

◆ G4CascadeInterface()

G4CascadeInterface::G4CascadeInterface ( const G4String & name = "BertiniCascade")

Definition at line 147 of file G4CascadeInterface.cc.

150 maximumTries(20), numberOfTries(0),
151 collider(new G4InuclCollider), balance(new G4CascadeCheckBalance(name)),
152 ltcollider(new G4LightTargetCollider),
153 bullet(0), target(0), output(new G4CollisionOutput), secID(-1)
154{
155 // Set up global objects for master thread or sequential build
157
158 SetEnergyMomentumCheckLevels(5*perCent, 10*MeV);
159 balance->setLimits(5*perCent, 10*MeV/GeV); // Bertini internal units
161
164 } else if ( G4CascadeParameters::useAbla() ) {
166 } else {
168 }
169
170 secID = G4PhysicsModelCatalog::GetModelID( "model_BertiniCascade" );
171}
void setLimits(G4double relative, G4double absolute)
void SetVerboseLevel(G4int verbose)
static G4bool usePreCompound()
static const G4String & randomFile()
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
static G4int GetModelID(const G4int modelIndex)
G4VIntraNuclearTransportModel(const G4String &mName="CascadeModel", G4VPreCompoundModel *ptr=nullptr)
const char * name(G4int ptype)
G4bool IsMasterThread()

◆ ~G4CascadeInterface()

G4CascadeInterface::~G4CascadeInterface ( )
virtual

Definition at line 173 of file G4CascadeInterface.cc.

173 {
174 clear();
175 delete collider; collider=0;
176 delete ltcollider; ltcollider = 0;
177 delete balance; balance=0;
178 delete output; output=0;
179}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 261 of file G4CascadeInterface.cc.

262 {
263 if (verboseLevel)
264 G4cout << " >>> G4CascadeInterface::ApplyYourself" << G4endl;
265
266 if (aTrack.GetKineticEnergy() < 0.) {
267 G4cerr << " >>> G4CascadeInterface got negative-energy track: "
268 << aTrack.GetDefinition()->GetParticleName() << " Ekin = "
269 << aTrack.GetKineticEnergy() << G4endl;
270 }
271
272#ifdef G4CASCADE_DEBUG_INTERFACE
273 static G4int counter(0);
274 counter++;
275 G4cerr << "Reaction number "<< counter << " "
276 << aTrack.GetDefinition()->GetParticleName() << " "
277 << aTrack.GetKineticEnergy() << " MeV" << G4endl;
278#endif
279
280 if (!randomFile.empty()) { // User requested random-seed capture
281 if (verboseLevel>1)
282 G4cout << " Saving random engine state to " << randomFile << G4endl;
284 }
285
287 clear();
288
289 // Abort processing if no interaction is possible
290 if (!IsApplicable(aTrack, theNucleus)) {
291 if (verboseLevel) G4cerr << " No interaction possible " << G4endl;
292 return NoInteraction(aTrack, theNucleus);
293 }
294
295 // If target A < 3 skip all cascade machinery and do scattering on
296 // nucleons
297
298 if (aTrack.GetDefinition() == G4Gamma::Gamma() &&
299 theNucleus.GetA_asInt() < 3) {
300 output->reset();
301 createBullet(aTrack);
302 createTarget(theNucleus);
303 // Due to binning, gamma-p cross sections between 130 MeV and the inelastic threshold
304 // (144 for pi0 p, 152 for pi+ n) are non-zero, causing energy non-conservation
305 // So, if Egamma is between 144 and 152, only pi0 p is allowed.
306
307 // Also, inelastic gamma-p cross section from G4PhotoNuclearCrossSection seems to be
308 // non-zero below between pi0 mass (135 MeV) and threshold (144 MeV)
309 ltcollider->collide(bullet, target, *output);
310
311 } else {
312
313 // Make conversion between native Geant4 and Bertini cascade classes.
314 if (!createBullet(aTrack)) {
315 if (verboseLevel) G4cerr << " Unable to create usable bullet" << G4endl;
316 return NoInteraction(aTrack, theNucleus);
317 }
318
319 if (!createTarget(theNucleus)) {
320 if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
321 return NoInteraction(aTrack, theNucleus);
322 }
323
324 // Different retry conditions for proton target vs. nucleus
325 const G4bool isHydrogen = (theNucleus.GetA_asInt() == 1);
326
327 numberOfTries = 0;
328 do { // we try to create inelastic interaction
329 if (verboseLevel > 1)
330 G4cout << " Generating cascade attempt " << numberOfTries << G4endl;
331
332 output->reset();
333 collider->collide(bullet, target, *output);
334 balance->collide(bullet, target, *output);
335
336 numberOfTries++;
337 /* Loop checking 08.06.2015 MHK */
338 } while ( isHydrogen ? retryInelasticProton() : retryInelasticNucleus() );
339
340 // Null event if unsuccessful
341 if (numberOfTries >= maximumTries) {
342 if (verboseLevel)
343 G4cout << " Cascade aborted after trials " << numberOfTries << G4endl;
344 return NoInteraction(aTrack, theNucleus);
345 }
346
347 // Abort job if energy or momentum are not conserved
348 if (!balance->okay()) {
350 return NoInteraction(aTrack, theNucleus);
351 }
352
353 // Successful cascade -- clean up and return
354 if (verboseLevel) {
355 G4cout << " Cascade output after trials " << numberOfTries << G4endl;
356 if (verboseLevel > 1) output->printCollisionOutput();
357 }
358
359 } // end cascade-style collisions
360
362
363 // Report violations of conservation laws in original frame
365
366 // Clean up and return final result;
367 clear();
368/*
369 G4int nSec = theParticleChange.GetNumberOfSecondaries();
370 for (G4int i = 0; i < nSec; i++) {
371 G4HadSecondary* sec = theParticleChange.GetSecondary(i);
372 G4DynamicParticle* dp = sec->GetParticle();
373 if (dp->GetDefinition()->GetParticleName() == "neutron")
374 G4cout << dp->GetDefinition()->GetParticleName() << " has "
375 << dp->GetKineticEnergy()/MeV << " MeV " << G4endl;
376 }
377*/
378 return &theParticleChange;
379}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static void saveEngineStatus(const char filename[]="Config.conf")
Definition Random.cc:278
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4bool retryInelasticNucleus() const
G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4HadFinalState * NoInteraction(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4bool createBullet(const G4HadProjectile &aTrack)
G4bool createTarget(G4Nucleus &theNucleus)
G4bool retryInelasticProton() const
void printCollisionOutput(std::ostream &os=G4cout) const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
const G4String & GetParticleName() const

◆ checkFinalResult()

void G4CascadeInterface::checkFinalResult ( )
protected

Definition at line 660 of file G4CascadeInterface.cc.

660 {
661 balance->collide(bullet, target, *output);
662
663 if (verboseLevel > 2) {
664 if (!balance->baryonOkay()) {
665 G4cerr << "ERROR: no baryon number conservation, sum of baryons = "
666 << balance->deltaB() << G4endl;
667 }
668
669 if (!balance->chargeOkay()) {
670 G4cerr << "ERROR: no charge conservation, sum of charges = "
671 << balance->deltaQ() << G4endl;
672 }
673
674 if (std::abs(balance->deltaKE()) > 0.01 ) { // GeV
675 G4cerr << "Kinetic energy conservation violated by "
676 << balance->deltaKE() << " GeV" << G4endl;
677 }
678
679 G4double eInit = bullet->getEnergy() + target->getEnergy();
680 G4double eFinal = eInit + balance->deltaE();
681
682 G4cout << "Initial energy " << eInit << " final energy " << eFinal
683 << "\nTotal energy conservation at level "
684 << balance->deltaE() * GeV << " MeV" << G4endl;
685
686 if (balance->deltaKE() > 5.0e-5 ) { // 0.05 MeV
687 G4cerr << "FATAL ERROR: kinetic energy created "
688 << balance->deltaKE() * GeV << " MeV" << G4endl;
689 }
690 }
691}
double G4double
Definition G4Types.hh:83
G4double getEnergy() const

Referenced by ApplyYourself().

◆ clear()

void G4CascadeInterface::clear ( )
protected

Definition at line 199 of file G4CascadeInterface.cc.

199 {
200 bullet=0;
201 target=0;
202}

Referenced by ApplyYourself(), Propagate(), and ~G4CascadeInterface().

◆ copyOutputToHadronicResult()

void G4CascadeInterface::copyOutputToHadronicResult ( )
protected

Definition at line 589 of file G4CascadeInterface.cc.

589 {
590 if (verboseLevel > 1)
591 G4cout << " >>> G4CascadeInterface::copyOutputToHadronicResult" << G4endl;
592
593 const std::vector<G4InuclNuclei>& outgoingNuclei = output->getOutgoingNuclei();
594 const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
595
598
599 // Get outcoming particles
600 if (!particles.empty()) {
601 particleIterator ipart = particles.begin();
602 for (; ipart != particles.end(); ipart++) {
604 }
605 }
606
607 // get nuclei fragments
608 if (!outgoingNuclei.empty()) {
609 nucleiIterator ifrag = outgoingNuclei.begin();
610 for (; ifrag != outgoingNuclei.end(); ifrag++) {
612 }
613 }
614}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
@ stopAndKill
G4DynamicParticle * makeDynamicParticle(const G4InuclElementaryParticle &iep) const
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)

Referenced by ApplyYourself().

◆ copyOutputToReactionProducts()

G4ReactionProductVector * G4CascadeInterface::copyOutputToReactionProducts ( )
protected

Definition at line 616 of file G4CascadeInterface.cc.

616 {
617 if (verboseLevel > 1)
618 G4cout << " >>> G4CascadeInterface::copyOutputToReactionProducts" << G4endl;
619
620 const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
621 const std::vector<G4InuclNuclei>& fragments = output->getOutgoingNuclei();
622
624
625 G4ReactionProduct* rp = 0; // Buffers to create outgoing tracks
626 G4DynamicParticle* dp = 0;
627
628 // Get outcoming particles
629 if (!particles.empty()) {
630 particleIterator ipart = particles.begin();
631 for (; ipart != particles.end(); ipart++) {
632 rp = new G4ReactionProduct;
633 dp = makeDynamicParticle(*ipart);
634 (*rp) = (*dp); // This does all the necessary copying
635 rp->SetCreatorModelID(secID);
636 propResult->push_back(rp);
637 delete dp;
638 }
639 }
640
641 // get nuclei fragments
642 if (!fragments.empty()) {
643 nucleiIterator ifrag = fragments.begin();
644 for (; ifrag != fragments.end(); ifrag++) {
645 rp = new G4ReactionProduct;
646 dp = makeDynamicParticle(*ifrag);
647 (*rp) = (*dp); // This does all the necessary copying
648 rp->SetCreatorModelID(secID);
649 propResult->push_back(rp);
650 delete dp;
651 }
652 }
653
654 return propResult;
655}
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetCreatorModelID(const G4int mod)

Referenced by Propagate().

◆ coulombBarrierViolation()

G4bool G4CascadeInterface::coulombBarrierViolation ( ) const
protected

Definition at line 696 of file G4CascadeInterface.cc.

696 {
697 G4bool violated = false; // by default coulomb analysis is OK
698
699 const G4double coulumbBarrier = 8.7 * MeV/GeV; // Bertini uses GeV
700
701 const std::vector<G4InuclElementaryParticle>& p =
702 output->getOutgoingParticles();
703
704 for (particleIterator ipart=p.begin(); ipart != p.end(); ipart++) {
705 if (ipart->type() == proton) {
706 violated |= (ipart->getKineticEnergy() < coulumbBarrier);
707 }
708 }
709
710 return violated;
711}

Referenced by retryInelasticNucleus().

◆ createBullet()

G4bool G4CascadeInterface::createBullet ( const G4HadProjectile & aTrack)
protected

Definition at line 474 of file G4CascadeInterface.cc.

474 {
475 const G4ParticleDefinition* trkDef = aTrack.GetDefinition();
476 G4int bulletType = 0; // For elementary particles
477 G4int bulletA = 0, bulletZ = 0; // For nucleus projectile
478
479 if (trkDef->GetAtomicMass() <= 1) {
480 bulletType = G4InuclElementaryParticle::type(trkDef);
481 } else {
482 bulletA = trkDef->GetAtomicMass();
483 bulletZ = trkDef->GetAtomicNumber();
484 }
485
486 if (0 == bulletType && 0 == bulletA*bulletZ) {
487 if (verboseLevel) {
488 G4cerr << " G4CascadeInterface: " << trkDef->GetParticleName()
489 << " not usable as bullet." << G4endl;
490 }
491 bullet = 0;
492 return false;
493 }
494
495 // Code momentum and energy -- Bertini wants z-axis and GeV units
496 G4LorentzVector projectileMomentum = aTrack.Get4Momentum()/GeV;
497
498 // Rotation/boost to get from z-axis back to original frame
499 // According to bug report #1990 this rotation is unnecessary and causes
500 // irreproducibility. Verifed and fixed by DHW 27 Nov 2017
501 // bulletInLabFrame = G4LorentzRotation::IDENTITY; // Initialize
502 // bulletInLabFrame.rotateZ(-projectileMomentum.phi());
503 // bulletInLabFrame.rotateY(-projectileMomentum.theta());
504 // bulletInLabFrame.invert();
505
506 G4LorentzVector momentumBullet(0., 0., projectileMomentum.rho(),
507 projectileMomentum.e());
508
509 if (G4InuclElementaryParticle::valid(bulletType)) {
510 hadronBullet.fill(momentumBullet, bulletType);
511 bullet = &hadronBullet;
512 } else {
513 nucleusBullet.fill(momentumBullet, bulletA, bulletZ);
514 bullet = &nucleusBullet;
515 }
516
517 if (verboseLevel > 2) G4cout << "Bullet: \n" << *bullet << G4endl;
518
519 return true;
520}
const G4LorentzVector & Get4Momentum() const
void fill(G4int ityp, Model model=DefaultModel)
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
G4int GetAtomicNumber() const
G4int GetAtomicMass() const

Referenced by ApplyYourself(), and Propagate().

◆ createTarget() [1/3]

G4bool G4CascadeInterface::createTarget ( G4int A,
G4int Z )
protected

Definition at line 533 of file G4CascadeInterface.cc.

533 {
534 if (A > 1) {
535 nucleusTarget.fill(A, Z);
536 target = &nucleusTarget;
537 } else {
538 hadronTarget.fill(0., (Z==1?proton:neutron));
539 target = &hadronTarget;
540 }
541
542 if (verboseLevel > 2) G4cout << "Target: \n" << *target << G4endl;
543
544 return true; // Right now, target never fails
545}
const G4double A[17]

◆ createTarget() [2/3]

G4bool G4CascadeInterface::createTarget ( G4Nucleus & theNucleus)
protected

Definition at line 525 of file G4CascadeInterface.cc.

525 {
526 return createTarget(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt());
527}
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105

Referenced by ApplyYourself(), createTarget(), createTarget(), and Propagate().

◆ createTarget() [3/3]

G4bool G4CascadeInterface::createTarget ( G4V3DNucleus * theNucleus)
protected

Definition at line 529 of file G4CascadeInterface.cc.

529 {
530 return createTarget(theNucleus->GetMassNumber(), theNucleus->GetCharge());
531}
virtual G4int GetCharge()=0
virtual G4int GetMassNumber()=0

◆ DumpConfiguration()

void G4CascadeInterface::DumpConfiguration ( std::ostream & outFile) const
virtual

Definition at line 195 of file G4CascadeInterface.cc.

195 {
197}
static void DumpConfiguration(std::ostream &os)

◆ Initialize()

void G4CascadeInterface::Initialize ( )
static

Definition at line 207 of file G4CascadeInterface.cc.

207 {
212 if (!ch || !pn || !nn || !pp) return; // Avoid "unused variables"
213}
static const G4CascadeChannel * GetTable(G4int initialState)
static G4Dineutron * Definition()
static G4Diproton * Definition()
Definition G4Diproton.cc:67
static G4UnboundPN * Definition()

Referenced by G4CascadeInterface().

◆ IsApplicable() [1/2]

G4bool G4CascadeInterface::IsApplicable ( const G4HadProjectile & aTrack,
G4Nucleus & theNucleus )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 244 of file G4CascadeInterface.cc.

245 {
246 return IsApplicable(aTrack.GetDefinition());
247}

Referenced by ApplyYourself(), IsApplicable(), and G4HadronicAbsorptionBertini::IsApplicable().

◆ IsApplicable() [2/2]

G4bool G4CascadeInterface::IsApplicable ( const G4ParticleDefinition * aPD) const

Definition at line 249 of file G4CascadeInterface.cc.

249 {
250 if (aPD->GetAtomicMass() > 1) return true; // Nuclei are okay
251
252 // Valid particle and have interactions available
255}

◆ makeDynamicParticle() [1/2]

G4DynamicParticle * G4CascadeInterface::makeDynamicParticle ( const G4InuclElementaryParticle & iep) const
protected

Definition at line 550 of file G4CascadeInterface.cc.

551 {
552 G4int outgoingType = iep.type();
553
554 if (iep.quasi_deutron()) {
555 G4cerr << " ERROR: G4CascadeInterface incompatible particle type "
556 << outgoingType << G4endl;
557 return 0;
558 }
559
560 // Copy local G4DynPart to public output (handle kaon mixing specially)
561 if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
562 G4ThreeVector momDir = iep.getMomentum().vect().unit();
563 G4double ekin = iep.getKineticEnergy()*GeV; // Bertini -> G4 units
564
566 if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
567
568 return new G4DynamicParticle(pd, momDir, ekin);
569 } else {
570 return new G4DynamicParticle(iep.getDynamicParticle());
571 }
572
573 return 0; // Should never get here!
574}
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
G4double getKineticEnergy() const
G4LorentzVector getMomentum() const
const G4DynamicParticle & getDynamicParticle() const
static G4KaonZeroLong * Definition()
static G4KaonZeroShort * Definition()

Referenced by copyOutputToHadronicResult(), and copyOutputToReactionProducts().

◆ makeDynamicParticle() [2/2]

G4DynamicParticle * G4CascadeInterface::makeDynamicParticle ( const G4InuclNuclei & inuc) const
protected

Definition at line 577 of file G4CascadeInterface.cc.

577 {
578 if (verboseLevel > 2) {
579 G4cout << " Nuclei fragment: \n" << inuc << G4endl;
580 }
581
582 // Copy local G4DynPart to public output
583 return new G4DynamicParticle(inuc.getDynamicParticle());
584}

◆ ModelDescription()

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

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 181 of file G4CascadeInterface.cc.

182{
183 outFile << "The Bertini-style cascade implements the inelastic scattering\n"
184 << "of hadrons by nuclei. Nucleons, pions, kaons and hyperons\n"
185 << "from 0 to 15 GeV may be used as projectiles in this model.\n"
186 << "Final state hadrons are produced by a classical cascade\n"
187 << "consisting of individual hadron-nucleon scatterings which use\n"
188 << "free-space partial cross sections, corrected for various\n"
189 << "nuclear medium effects. The target nucleus is modeled as a\n"
190 << "set of 1, 3 or 6 spherical shells, in which scattered hadrons\n"
191 << "travel in straight lines until they are reflected from or\n"
192 << "transmitted through shell boundaries.\n";
193}

◆ NoInteraction()

G4HadFinalState * G4CascadeInterface::NoInteraction ( const G4HadProjectile & aTrack,
G4Nucleus & theNucleus )
protected

Definition at line 457 of file G4CascadeInterface.cc.

458 {
459 if (verboseLevel)
460 G4cout << " >>> G4CascadeInterface::NoInteraction" << G4endl;
461
464
465 G4double ekin = aTrack.GetKineticEnergy()>0. ? aTrack.GetKineticEnergy() : 0.;
466 theParticleChange.SetEnergyChange(ekin); // Protect against rounding
467
468 return &theParticleChange;
469}
@ isAlive

Referenced by ApplyYourself().

◆ Propagate()

G4ReactionProductVector * G4CascadeInterface::Propagate ( G4KineticTrackVector * theSecondaries,
G4V3DNucleus * theNucleus )
virtual

Implements G4VIntraNuclearTransportModel.

Definition at line 382 of file G4CascadeInterface.cc.

383 {
384 if (verboseLevel) G4cout << " >>> G4CascadeInterface::Propagate" << G4endl;
385
386#ifdef G4CASCADE_DEBUG_INTERFACE
387 if (verboseLevel>1) {
388 G4cout << " G4V3DNucleus A " << theNucleus->GetMassNumber()
389 << " Z " << theNucleus->GetCharge()
390 << "\n " << theSecondaries->size() << " secondaries:" << G4endl;
391 for (size_t i=0; i<theSecondaries->size(); i++) {
392 G4KineticTrack* kt = (*theSecondaries)[i];
393 G4cout << " " << i << ": " << kt->GetDefinition()->GetParticleName()
394 << " p " << kt->Get4Momentum() << " @ " << kt->GetPosition()
395 << " t " << kt->GetFormationTime() << G4endl;
396 }
397 }
398#endif
399
400 if (!randomFile.empty()) { // User requested random-seed capture
401 if (verboseLevel>1)
402 G4cout << " Saving random engine state to " << randomFile << G4endl;
404 }
405
407 clear();
408
409 // Process input secondaries list to eliminate resonances
410 G4DecayKineticTracks decay(theSecondaries);
411
412 // NOTE: Requires 9.4-ref-03 mods to base class and G4TheoFSGenerator
413 const G4HadProjectile* projectile = GetPrimaryProjectile();
414 if (projectile) createBullet(*projectile);
415
416 if (!createTarget(theNucleus)) {
417 if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
418 return 0; // FIXME: This will cause a segfault later
419 }
420
421 numberOfTries = 0;
422 do {
423 if (verboseLevel > 1)
424 G4cout << " Generating rescatter attempt " << numberOfTries << G4endl;
425
426 output->reset();
427 collider->rescatter(bullet, theSecondaries, theNucleus, *output);
428 balance->collide(bullet, target, *output);
429
430 numberOfTries++;
431 // FIXME: retry checks will SEGFAULT until we can define the bullet!
432 } while (retryInelasticNucleus()); /* Loop checking 08.06.2015 MHK */
433
434 // Check whether repeated attempts have all failed; report and exit
435 if (numberOfTries >= maximumTries && !balance->okay()) {
436 throwNonConservationFailure(); // This terminates the job
437 }
438
439 // Successful cascade -- clean up and return
440 if (verboseLevel) {
441 G4cout << " Cascade rescatter after trials " << numberOfTries << G4endl;
442 if (verboseLevel > 1) output->printCollisionOutput();
443 }
444
445 // Does calling code take ownership? I hope so!
447
448 // Clean up and and return final result
449 clear();
450 return propResult;
451}
G4ReactionProductVector * copyOutputToReactionProducts()
void rescatter(G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
G4double GetFormationTime() const
const G4ThreeVector & GetPosition() const
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
const G4HadProjectile * GetPrimaryProjectile() const
ParticleList decay(Cluster *const c)
Carries out a cluster decay.

◆ retryInelasticNucleus()

G4bool G4CascadeInterface::retryInelasticNucleus ( ) const
protected

Definition at line 747 of file G4CascadeInterface.cc.

747 {
748 // Quantities necessary for conditional block below
749 G4int npart = output->numberOfOutgoingParticles();
750 G4int nfrag = output->numberOfOutgoingNuclei();
751
752 const G4ParticleDefinition* firstOut = (npart == 0) ? 0 :
753 output->getOutgoingParticles().begin()->getDefinition();
754
755#ifdef G4CASCADE_DEBUG_INTERFACE
756 // Report on all retry conditions, in order of return logic
757 G4cout << " retryInelasticNucleus: numberOfTries "
758 << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
759 << "\n retryInelasticNucleus: AND outputParticles "
760 << ((npart != 0) ? "NON-ZERO (t)" : "EMPTY (f)")
761#ifdef G4CASCADE_COULOMB_DEV
762 << "\n retryInelasticNucleus: AND coulombBarrier (COULOMB_DEV) "
763 << (coulombBarrierViolation() ? "VIOLATED (t)" : "PASSED (f)")
764 << "\n retryInelasticNucleus: AND collision type (COULOMB_DEV) "
765 << ((npart+nfrag > 2) ? "INELASTIC (t)" : "ELASTIC (f)")
766#else
767 << "\n retryInelasticNucleus: AND collision type "
768 << ((npart+nfrag < 3) ? "ELASTIC (t)" : "INELASTIC (f)")
769 << "\n retryInelasticNucleus: AND Leading particle bullet "
770 << ((firstOut == bullet->getDefinition()) ? "YES (t)" : "NO (f)")
771#endif
772 << "\n retryInelasticNucleus: OR conservation "
773 << (!balance->okay() ? "FAILED (t)" : "PASSED (f)")
774 << G4endl;
775#endif
776
777 return ( (numberOfTries < maximumTries) &&
778 ( ((npart != 0) &&
779#ifdef G4CASCADE_COULOMB_DEV
780 (coulombBarrierViolation() && npart+nfrag > 2)
781#else
782 (npart+nfrag < 3 && firstOut == bullet->getDefinition())
783#endif
784 )
785#ifndef G4CASCADE_SKIP_ECONS
786 || (!balance->okay())
787#endif
788 )
789 );
790}
G4bool coulombBarrierViolation() const
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
const G4ParticleDefinition * getDefinition() const

Referenced by ApplyYourself(), and Propagate().

◆ retryInelasticProton()

G4bool G4CascadeInterface::retryInelasticProton ( ) const
protected

Definition at line 715 of file G4CascadeInterface.cc.

715 {
716 const std::vector<G4InuclElementaryParticle>& out =
717 output->getOutgoingParticles();
718
719#ifdef G4CASCADE_DEBUG_INTERFACE
720 // Report on all retry conditions, in order of return logic
721 G4cout << " retryInelasticProton: number of Tries "
722 << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
723 << "\n retryInelasticProton: AND collision type ";
724 if (out.empty()) G4cout << "FAILED" << G4endl;
725 else {
726 G4cout << (out.size() == 2 ? "ELASTIC (t)" : "INELASTIC (f)")
727 << "\n retryInelasticProton: AND Leading particles bullet "
728 << (out.size() >= 2 &&
729 (out[0].getDefinition() == bullet->getDefinition() ||
730 out[1].getDefinition() == bullet->getDefinition())
731 ? "YES (t)" : "NO (f)")
732 << G4endl;
733 }
734#endif
735
736 return ( (numberOfTries < maximumTries) &&
737 (out.empty() ||
738 (out.size() == 2 &&
739 (out[0].getDefinition() == bullet->getDefinition() ||
740 out[1].getDefinition() == bullet->getDefinition())))
741 );
742}

Referenced by ApplyYourself().

◆ SetVerboseLevel()

void G4CascadeInterface::SetVerboseLevel ( G4int verbose)

Definition at line 234 of file G4CascadeInterface.cc.

234 {
236 collider->setVerboseLevel(verbose);
237 balance->setVerboseLevel(verbose);
238 output->setVerboseLevel(verbose);
239}
void setVerboseLevel(G4int verbose)
void SetVerboseLevel(G4int value)
void setVerboseLevel(G4int verbose=0)
virtual void setVerboseLevel(G4int verbose=0)

Referenced by G4CascadeInterface().

◆ throwNonConservationFailure()

void G4CascadeInterface::throwNonConservationFailure ( )
protected

Definition at line 796 of file G4CascadeInterface.cc.

796 {
797 // NOTE: Once G4HadronicException is changed, use the following line!
798 // G4ExceptionDescription errInfo;
799 std::ostream& errInfo = G4cerr;
800
801 errInfo << " >>> G4CascadeInterface has non-conserving"
802 << " cascade after " << numberOfTries << " attempts." << G4endl;
803
804 G4String throwMsg = "G4CascadeInterface - ";
805 if (!balance->energyOkay()) {
806 throwMsg += "Energy";
807 errInfo << " Energy conservation violated by " << balance->deltaE()
808 << " GeV (" << balance->relativeE() << ")" << G4endl;
809 }
810
811 if (!balance->momentumOkay()) {
812 throwMsg += "Momentum";
813 errInfo << " Momentum conservation violated by " << balance->deltaP()
814 << " GeV/c (" << balance->relativeP() << ")" << G4endl;
815 }
816
817 if (!balance->baryonOkay()) {
818 throwMsg += "Baryon number";
819 errInfo << " Baryon number violated by " << balance->deltaB() << G4endl;
820 }
821
822 if (!balance->chargeOkay()) {
823 throwMsg += "Charge";
824 errInfo << " Charge conservation violated by " << balance->deltaQ()
825 << G4endl;
826 }
827
828 errInfo << " Final event output, for debugging:\n"
829 << " Bullet: \n" << *bullet << G4endl
830 << " Target: \n" << *target << G4endl;
831 output->printCollisionOutput(errInfo);
832
833 throwMsg += " non-conservation. More info in output.";
834 throw G4HadronicException(__FILE__, __LINE__, throwMsg); // Job ends here!
835}

Referenced by ApplyYourself(), and Propagate().

◆ useAblaDeexcitation()

void G4CascadeInterface::useAblaDeexcitation ( )

Definition at line 227 of file G4CascadeInterface.cc.

227 {
228 collider->useAblaDeexcitation();
229}

Referenced by G4CascadeInterface().

◆ useCascadeDeexcitation()

void G4CascadeInterface::useCascadeDeexcitation ( )

Definition at line 219 of file G4CascadeInterface.cc.

219 {
220 collider->useCascadeDeexcitation();
221}

Referenced by G4CascadeInterface().

◆ usePreCompoundDeexcitation()

void G4CascadeInterface::usePreCompoundDeexcitation ( )

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