Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4CascadeInterface Class Reference

#include <G4CascadeInterface.hh>

Inheritance diagram for G4CascadeInterface:
G4VIntraNuclearTransportModel G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
virtual void DumpConfiguration (std::ostream &outFile) const
 
 G4CascadeInterface (const G4String &name="BertiniCascade")
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
const G4StringGetModelName () const
 
G4double GetRecoilEnergyThreshold () const
 
G4int GetVerboseLevel () const
 
virtual void InitialiseModel ()
 
G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
G4bool IsApplicable (const G4ParticleDefinition *aPD) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
virtual void ModelDescription (std::ostream &outFile) const
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4VIntraNuclearTransportModel &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4VIntraNuclearTransportModel &right) const =delete
 
G4ReactionProductVectorPropagate (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
 
virtual void PropagateModelDescription (std::ostream &outFile) const
 
virtual G4ReactionProductVectorPropagateNuclNucl (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void Set3DNucleus (G4V3DNucleus *const value)
 
void SetDeExcitation (G4VPreCompoundModel *ptr)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetPrimaryProjectile (const G4HadProjectile &aPrimary)
 
void SetRecoilEnergyThreshold (G4double val)
 
void SetVerboseLevel (G4int verbose)
 
void useCascadeDeexcitation ()
 
void usePreCompoundDeexcitation ()
 
virtual ~G4CascadeInterface ()
 

Static Public Member Functions

static void Initialize ()
 

Protected Member Functions

void Block ()
 
void checkFinalResult ()
 
void clear ()
 
void copyOutputToHadronicResult ()
 
G4ReactionProductVectorcopyOutputToReactionProducts ()
 
G4bool coulombBarrierViolation () const
 
G4bool createBullet (const G4HadProjectile &aTrack)
 
G4bool createTarget (G4int A, G4int Z)
 
G4bool createTarget (G4Nucleus &theNucleus)
 
G4bool createTarget (G4V3DNucleus *theNucleus)
 
G4V3DNucleusGet3DNucleus () const
 
G4VPreCompoundModelGetDeExcitation () const
 
const G4HadProjectileGetPrimaryProjectile () const
 
G4bool IsBlocked () const
 
G4DynamicParticlemakeDynamicParticle (const G4InuclElementaryParticle &iep) const
 
G4DynamicParticlemakeDynamicParticle (const G4InuclNuclei &inuc) const
 
G4HadFinalStateNoInteraction (const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
 
G4bool retryInelasticNucleus () const
 
G4bool retryInelasticProton () const
 
void SetModelName (const G4String &nam)
 
void throwNonConservationFailure ()
 

Protected Attributes

G4bool isBlocked
 
G4V3DNucleusthe3DNucleus
 
G4VPreCompoundModeltheDeExcitation
 
G4double theMaxEnergy
 
G4double theMinEnergy
 
G4HadFinalState theParticleChange
 
const G4HadProjectilethePrimaryProjectile
 
G4String theTransportModelName
 
G4int verboseLevel
 

Private Member Functions

 G4CascadeInterface (const G4CascadeInterface &)
 
G4bool operator!= (const G4CascadeInterface &right) const
 
G4CascadeInterfaceoperator= (const G4CascadeInterface &)
 
G4bool operator== (const G4CascadeInterface &right) const
 

Private Attributes

G4CascadeCheckBalancebalance
 
G4InuclParticlebullet
 
G4InuclCollidercollider
 
std::pair< G4double, G4doubleepCheckLevels
 
G4InuclElementaryParticle hadronBullet
 
G4InuclElementaryParticle hadronTarget
 
G4LightTargetColliderltcollider
 
const G4int maximumTries
 
G4InuclNuclei nucleusBullet
 
G4InuclNuclei nucleusTarget
 
G4int numberOfTries
 
G4CollisionOutputoutput
 
const G4String randomFile
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
G4int secID
 
G4InuclParticletarget
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements
 
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
 
G4String theModelName
 

Detailed Description

Definition at line 88 of file G4CascadeInterface.hh.

Constructor & Destructor Documentation

◆ G4CascadeInterface() [1/2]

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

Definition at line 147 of file G4CascadeInterface.cc.

153 bullet(0), target(0), output(new G4CollisionOutput), secID(-1)
154{
155 // Set up global objects for master thread or sequential build
157
159 balance->setLimits(5*perCent, 10*MeV/GeV); // Bertini internal units
161
164 else
166
167 secID = G4PhysicsModelCatalog::GetModelID( "model_BertiniCascade" );
168}
static constexpr double perCent
Definition: G4SIunits.hh:325
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
void setLimits(G4double relative, G4double absolute)
G4CollisionOutput * output
G4CascadeCheckBalance * balance
G4InuclParticle * bullet
G4LightTargetCollider * ltcollider
const G4String randomFile
G4InuclParticle * target
void SetVerboseLevel(G4int verbose)
G4InuclCollider * collider
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()
Definition: G4Threading.cc:124

References balance, G4PhysicsModelCatalog::GetModelID(), GeV, Initialize(), G4Threading::IsMasterThread(), MeV, perCent, secID, G4HadronicInteraction::SetEnergyMomentumCheckLevels(), G4CascadeCheckBalance::setLimits(), SetVerboseLevel(), useCascadeDeexcitation(), G4CascadeParameters::usePreCompound(), usePreCompoundDeexcitation(), and G4CascadeParameters::verbose().

◆ ~G4CascadeInterface()

G4CascadeInterface::~G4CascadeInterface ( )
virtual

Definition at line 170 of file G4CascadeInterface.cc.

170 {
171 clear();
172 delete collider; collider=0;
173 delete ltcollider; ltcollider = 0;
174 delete balance; balance=0;
175 delete output; output=0;
176}

References balance, clear(), collider, ltcollider, and output.

◆ G4CascadeInterface() [2/2]

G4CascadeInterface::G4CascadeInterface ( const G4CascadeInterface )
private

Member Function Documentation

◆ ActivateFor() [1/2]

void G4HadronicInteraction::ActivateFor ( const G4Element anElement)
inlineinherited

◆ ActivateFor() [2/2]

void G4HadronicInteraction::ActivateFor ( const G4Material aMaterial)
inlineinherited

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 254 of file G4CascadeInterface.cc.

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

References balance, bullet, checkFinalResult(), clear(), G4HadFinalState::Clear(), G4InuclCollider::collide(), G4LightTargetCollider::collide(), G4CascadeCheckBalance::collide(), collider, copyOutputToHadronicResult(), createBullet(), createTarget(), G4cerr, G4cout, G4endl, G4Gamma::Gamma(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4HadProjectile::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), IsApplicable(), ltcollider, maximumTries, NoInteraction(), numberOfTries, G4CascadeCheckBalance::okay(), output, G4CollisionOutput::printCollisionOutput(), randomFile, G4CollisionOutput::reset(), retryInelasticNucleus(), retryInelasticProton(), CLHEP::HepRandom::saveEngineStatus(), target, G4HadronicInteraction::theParticleChange, throwNonConservationFailure(), and G4HadronicInteraction::verboseLevel.

Referenced by G4ElectroVDNuclearModel::CalculateHadronicVertex(), and G4MuonVDNuclearModel::CalculateHadronicVertex().

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4HadronicInteraction::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ checkFinalResult()

void G4CascadeInterface::checkFinalResult ( )
protected

Definition at line 653 of file G4CascadeInterface.cc.

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

References balance, G4CascadeCheckBalance::baryonOkay(), bullet, G4CascadeCheckBalance::chargeOkay(), G4CascadeCheckBalance::collide(), G4CascadeCheckBalance::deltaB(), G4CascadeCheckBalance::deltaE(), G4CascadeCheckBalance::deltaKE(), G4CascadeCheckBalance::deltaQ(), G4cerr, G4cout, G4endl, G4InuclParticle::getEnergy(), GeV, output, target, and G4HadronicInteraction::verboseLevel.

Referenced by ApplyYourself().

◆ clear()

void G4CascadeInterface::clear ( )
protected

Definition at line 196 of file G4CascadeInterface.cc.

196 {
197 bullet=0;
198 target=0;
199}

References bullet, and target.

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

◆ copyOutputToHadronicResult()

void G4CascadeInterface::copyOutputToHadronicResult ( )
protected

Definition at line 582 of file G4CascadeInterface.cc.

582 {
583 if (verboseLevel > 1)
584 G4cout << " >>> G4CascadeInterface::copyOutputToHadronicResult" << G4endl;
585
586 const std::vector<G4InuclNuclei>& outgoingNuclei = output->getOutgoingNuclei();
587 const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
588
591
592 // Get outcoming particles
593 if (!particles.empty()) {
594 particleIterator ipart = particles.begin();
595 for (; ipart != particles.end(); ipart++) {
597 }
598 }
599
600 // get nuclei fragments
601 if (!outgoingNuclei.empty()) {
602 nucleiIterator ifrag = outgoingNuclei.begin();
603 for (; ifrag != outgoingNuclei.end(); ifrag++) {
605 }
606 }
607}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
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)

References G4HadFinalState::AddSecondary(), G4cout, G4endl, G4CollisionOutput::getOutgoingNuclei(), G4CollisionOutput::getOutgoingParticles(), makeDynamicParticle(), output, secID, G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, and G4HadronicInteraction::verboseLevel.

Referenced by ApplyYourself().

◆ copyOutputToReactionProducts()

G4ReactionProductVector * G4CascadeInterface::copyOutputToReactionProducts ( )
protected

Definition at line 609 of file G4CascadeInterface.cc.

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

References G4cout, G4endl, G4CollisionOutput::getOutgoingNuclei(), G4CollisionOutput::getOutgoingParticles(), makeDynamicParticle(), output, secID, G4ReactionProduct::SetCreatorModelID(), and G4HadronicInteraction::verboseLevel.

Referenced by Propagate().

◆ coulombBarrierViolation()

G4bool G4CascadeInterface::coulombBarrierViolation ( ) const
protected

Definition at line 689 of file G4CascadeInterface.cc.

689 {
690 G4bool violated = false; // by default coulomb analysis is OK
691
692 const G4double coulumbBarrier = 8.7 * MeV/GeV; // Bertini uses GeV
693
694 const std::vector<G4InuclElementaryParticle>& p =
696
697 for (particleIterator ipart=p.begin(); ipart != p.end(); ipart++) {
698 if (ipart->type() == proton) {
699 violated |= (ipart->getKineticEnergy() < coulumbBarrier);
700 }
701 }
702
703 return violated;
704}

References G4CollisionOutput::getOutgoingParticles(), GeV, MeV, output, and G4InuclParticleNames::proton.

Referenced by retryInelasticNucleus().

◆ createBullet()

G4bool G4CascadeInterface::createBullet ( const G4HadProjectile aTrack)
protected

Definition at line 467 of file G4CascadeInterface.cc.

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

References bullet, CLHEP::HepLorentzVector::e(), G4InuclNuclei::fill(), G4InuclElementaryParticle::fill(), G4cerr, G4cout, G4endl, G4HadProjectile::Get4Momentum(), G4ParticleDefinition::GetAtomicMass(), G4ParticleDefinition::GetAtomicNumber(), G4HadProjectile::GetDefinition(), G4ParticleDefinition::GetParticleName(), GeV, hadronBullet, nucleusBullet, CLHEP::HepLorentzVector::rho(), G4InuclElementaryParticle::type(), G4InuclElementaryParticle::valid(), and G4HadronicInteraction::verboseLevel.

Referenced by ApplyYourself(), and Propagate().

◆ createTarget() [1/3]

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

Definition at line 526 of file G4CascadeInterface.cc.

526 {
527 if (A > 1) {
530 } else {
531 hadronTarget.fill(0., (Z==1?proton:neutron));
533 }
534
535 if (verboseLevel > 2) G4cout << "Target: \n" << *target << G4endl;
536
537 return true; // Right now, target never fails
538}
const G4int Z[17]
const G4double A[17]
G4InuclElementaryParticle hadronTarget
G4InuclNuclei nucleusTarget

References A, G4InuclNuclei::fill(), G4InuclElementaryParticle::fill(), G4cout, G4endl, hadronTarget, G4InuclParticleNames::neutron, nucleusTarget, G4InuclParticleNames::proton, target, G4HadronicInteraction::verboseLevel, and Z.

◆ createTarget() [2/3]

G4bool G4CascadeInterface::createTarget ( G4Nucleus theNucleus)
protected

Definition at line 518 of file G4CascadeInterface.cc.

518 {
519 return createTarget(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt());
520}
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105

References createTarget(), G4Nucleus::GetA_asInt(), and G4Nucleus::GetZ_asInt().

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

◆ createTarget() [3/3]

G4bool G4CascadeInterface::createTarget ( G4V3DNucleus theNucleus)
protected

Definition at line 522 of file G4CascadeInterface.cc.

522 {
523 return createTarget(theNucleus->GetMassNumber(), theNucleus->GetCharge());
524}
virtual G4int GetCharge()=0
virtual G4int GetMassNumber()=0

References createTarget(), G4V3DNucleus::GetCharge(), and G4V3DNucleus::GetMassNumber().

◆ DeActivateFor() [1/2]

void G4HadronicInteraction::DeActivateFor ( const G4Element anElement)
inherited

Definition at line 186 of file G4HadronicInteraction.cc.

187{
188 Block();
189 theBlockedListElements.push_back(anElement);
190}
std::vector< const G4Element * > theBlockedListElements

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedListElements.

◆ DeActivateFor() [2/2]

void G4HadronicInteraction::DeActivateFor ( const G4Material aMaterial)
inherited

Definition at line 180 of file G4HadronicInteraction.cc.

181{
182 Block();
183 theBlockedList.push_back(aMaterial);
184}
std::vector< const G4Material * > theBlockedList

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedList.

Referenced by G4HadronHElasticPhysics::ConstructProcess().

◆ DumpConfiguration()

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

Definition at line 192 of file G4CascadeInterface.cc.

192 {
194}
static void DumpConfiguration(std::ostream &os)

References G4CascadeParameters::DumpConfiguration().

◆ Get3DNucleus()

G4V3DNucleus * G4VIntraNuclearTransportModel::Get3DNucleus ( ) const
inlineprotectedinherited

◆ GetDeExcitation()

G4VPreCompoundModel * G4VIntraNuclearTransportModel::GetDeExcitation ( ) const
inlineprotectedinherited

◆ GetEnergyMomentumCheckLevels()

std::pair< G4double, G4double > G4HadronicInteraction::GetEnergyMomentumCheckLevels ( ) const
virtualinherited

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4HadronicInteraction::GetFatalEnergyCheckLevels ( ) const
virtualinherited

Reimplemented in G4FissLib, G4LFission, G4LENDFission, G4ParticleHPCapture, G4ParticleHPElastic, G4ParticleHPFission, G4ParticleHPInelastic, and G4ParticleHPThermalScattering.

Definition at line 210 of file G4HadronicInteraction.cc.

211{
212 // default level of Check
213 return std::pair<G4double, G4double>(2.*perCent, 1. * GeV);
214}

References GeV, and perCent.

Referenced by G4HadronicProcess::CheckResult().

◆ GetMaxEnergy() [1/2]

G4double G4HadronicInteraction::GetMaxEnergy ( ) const
inlineinherited

◆ GetMaxEnergy() [2/2]

G4double G4HadronicInteraction::GetMaxEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 131 of file G4HadronicInteraction.cc.

133{
134 if(!IsBlocked()) { return theMaxEnergy; }
135 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return 0.0; }
136 if(!theMaxEnergyListElements.empty()) {
137 for(auto const& elmlist : theMaxEnergyListElements) {
138 if( anElement == elmlist.second )
139 { return elmlist.first; }
140 }
141 }
142 if(!theMaxEnergyList.empty()) {
143 for(auto const& matlist : theMaxEnergyList) {
144 if( aMaterial == matlist.second )
145 { return matlist.first; }
146 }
147 }
148 return theMaxEnergy;
149}
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements

References G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMaxEnergy, G4HadronicInteraction::theMaxEnergyList, and G4HadronicInteraction::theMaxEnergyListElements.

◆ GetMinEnergy() [1/2]

G4double G4HadronicInteraction::GetMinEnergy ( ) const
inlineinherited

◆ GetMinEnergy() [2/2]

G4double G4HadronicInteraction::GetMinEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 81 of file G4HadronicInteraction.cc.

83{
84 if(!IsBlocked()) { return theMinEnergy; }
85 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return DBL_MAX; }
86 if(!theMinEnergyListElements.empty()) {
87 for(auto const& elmlist : theMinEnergyListElements) {
88 if( anElement == elmlist.second )
89 { return elmlist.first; }
90 }
91 }
92 if(!theMinEnergyList.empty()) {
93 for(auto const & matlist : theMinEnergyList) {
94 if( aMaterial == matlist.second )
95 { return matlist.first; }
96 }
97 }
98 return theMinEnergy;
99}
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMinEnergy, G4HadronicInteraction::theMinEnergyList, and G4HadronicInteraction::theMinEnergyListElements.

◆ GetModelName()

const G4String & G4VIntraNuclearTransportModel::GetModelName ( ) const
inlineinherited

◆ GetPrimaryProjectile()

const G4HadProjectile * G4VIntraNuclearTransportModel::GetPrimaryProjectile ( ) const
inlineprotectedinherited

◆ GetRecoilEnergyThreshold()

G4double G4HadronicInteraction::GetRecoilEnergyThreshold ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4HadronicInteraction::GetVerboseLevel ( ) const
inlineinherited

Definition at line 109 of file G4HadronicInteraction.hh.

110 { return verboseLevel; }

References G4HadronicInteraction::verboseLevel.

◆ InitialiseModel()

void G4HadronicInteraction::InitialiseModel ( )
virtualinherited

◆ Initialize()

void G4CascadeInterface::Initialize ( )
static

◆ IsApplicable() [1/2]

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

Reimplemented from G4HadronicInteraction.

Definition at line 237 of file G4CascadeInterface.cc.

238 {
239 return IsApplicable(aTrack.GetDefinition());
240}

References G4HadProjectile::GetDefinition(), and IsApplicable().

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

◆ IsApplicable() [2/2]

G4bool G4CascadeInterface::IsApplicable ( const G4ParticleDefinition aPD) const

Definition at line 242 of file G4CascadeInterface.cc.

242 {
243 if (aPD->GetAtomicMass() > 1) return true; // Nuclei are okay
244
245 // Valid particle and have interactions available
248}

References G4ParticleDefinition::GetAtomicMass(), G4CascadeChannelTables::GetTable(), and G4InuclElementaryParticle::type().

◆ IsBlocked() [1/3]

G4bool G4HadronicInteraction::IsBlocked ( ) const
inlineprotectedinherited

◆ IsBlocked() [2/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Element anElement) const
inherited

Definition at line 202 of file G4HadronicInteraction.cc.

203{
204 for (auto const& elm : theBlockedListElements) {
205 if (anElement == elm) return true;
206 }
207 return false;
208}

References G4HadronicInteraction::theBlockedListElements.

◆ IsBlocked() [3/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Material aMaterial) const
inherited

Definition at line 193 of file G4HadronicInteraction.cc.

194{
195 for (auto const& mat : theBlockedList) {
196 if (aMaterial == mat) return true;
197 }
198 return false;
199}

References G4HadronicInteraction::theBlockedList.

◆ makeDynamicParticle() [1/2]

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

Definition at line 543 of file G4CascadeInterface.cc.

544 {
545 G4int outgoingType = iep.type();
546
547 if (iep.quasi_deutron()) {
548 G4cerr << " ERROR: G4CascadeInterface incompatible particle type "
549 << outgoingType << G4endl;
550 return 0;
551 }
552
553 // Copy local G4DynPart to public output (handle kaon mixing specially)
554 if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
555 G4ThreeVector momDir = iep.getMomentum().vect().unit();
556 G4double ekin = iep.getKineticEnergy()*GeV; // Bertini -> G4 units
557
559 if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
560
561 return new G4DynamicParticle(pd, momDir, ekin);
562 } else {
563 return new G4DynamicParticle(iep.getDynamicParticle());
564 }
565
566 return 0; // Should never get here!
567}
#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()

References G4KaonZeroLong::Definition(), G4KaonZeroShort::Definition(), G4cerr, G4endl, G4UniformRand, G4InuclParticle::getDynamicParticle(), G4InuclParticle::getKineticEnergy(), G4InuclParticle::getMomentum(), GeV, G4InuclParticleNames::kaonZero, G4InuclParticleNames::kaonZeroBar, G4InuclElementaryParticle::quasi_deutron(), G4InuclElementaryParticle::type(), CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

Referenced by copyOutputToHadronicResult(), and copyOutputToReactionProducts().

◆ makeDynamicParticle() [2/2]

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

Definition at line 570 of file G4CascadeInterface.cc.

570 {
571 if (verboseLevel > 2) {
572 G4cout << " Nuclei fragment: \n" << inuc << G4endl;
573 }
574
575 // Copy local G4DynPart to public output
576 return new G4DynamicParticle(inuc.getDynamicParticle());
577}

References G4cout, G4endl, G4InuclParticle::getDynamicParticle(), and G4HadronicInteraction::verboseLevel.

◆ ModelDescription()

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

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 178 of file G4CascadeInterface.cc.

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

◆ NoInteraction()

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

◆ operator!=() [1/3]

G4bool G4CascadeInterface::operator!= ( const G4CascadeInterface right) const
inlineprivate

Definition at line 158 of file G4CascadeInterface.hh.

158 {
159 return (this != &right);
160 }

◆ operator!=() [2/3]

G4bool G4HadronicInteraction::operator!= ( const G4HadronicInteraction right) const
deleteinherited

◆ operator!=() [3/3]

G4bool G4VIntraNuclearTransportModel::operator!= ( const G4VIntraNuclearTransportModel right) const
deleteinherited

◆ operator=()

G4CascadeInterface & G4CascadeInterface::operator= ( const G4CascadeInterface )
private

◆ operator==() [1/3]

G4bool G4CascadeInterface::operator== ( const G4CascadeInterface right) const
inlineprivate

Definition at line 154 of file G4CascadeInterface.hh.

154 {
155 return (this == &right);
156 }

◆ operator==() [2/3]

G4bool G4HadronicInteraction::operator== ( const G4HadronicInteraction right) const
deleteinherited

◆ operator==() [3/3]

G4bool G4VIntraNuclearTransportModel::operator== ( const G4VIntraNuclearTransportModel right) const
deleteinherited

◆ Propagate()

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

Implements G4VIntraNuclearTransportModel.

Definition at line 375 of file G4CascadeInterface.cc.

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

References balance, bullet, clear(), G4HadFinalState::Clear(), G4CascadeCheckBalance::collide(), collider, copyOutputToReactionProducts(), createBullet(), createTarget(), G4INCL::ClusterDecay::decay(), G4cerr, G4cout, G4endl, G4KineticTrack::Get4Momentum(), G4V3DNucleus::GetCharge(), G4KineticTrack::GetDefinition(), G4KineticTrack::GetFormationTime(), G4V3DNucleus::GetMassNumber(), G4ParticleDefinition::GetParticleName(), G4KineticTrack::GetPosition(), G4VIntraNuclearTransportModel::GetPrimaryProjectile(), maximumTries, numberOfTries, G4CascadeCheckBalance::okay(), output, G4CollisionOutput::printCollisionOutput(), randomFile, G4InuclCollider::rescatter(), G4CollisionOutput::reset(), retryInelasticNucleus(), CLHEP::HepRandom::saveEngineStatus(), target, G4HadronicInteraction::theParticleChange, throwNonConservationFailure(), and G4HadronicInteraction::verboseLevel.

◆ PropagateModelDescription()

void G4VIntraNuclearTransportModel::PropagateModelDescription ( std::ostream &  outFile) const
virtualinherited

Reimplemented in G4BinaryCascade, and G4GeneratorPrecompoundInterface.

Definition at line 58 of file G4VIntraNuclearTransportModel.cc.

59{
60 outFile << "G4VIntraNuclearTransportModel is abstract class, missing description.\n";
61}

Referenced by G4TheoFSGenerator::ModelDescription().

◆ PropagateNuclNucl()

G4ReactionProductVector * G4VIntraNuclearTransportModel::PropagateNuclNucl ( G4KineticTrackVector theSecondaries,
G4V3DNucleus theNucleus,
G4V3DNucleus theProjectileNucleus 
)
virtualinherited

Reimplemented in G4GeneratorPrecompoundInterface.

Definition at line 64 of file G4VIntraNuclearTransportModel.cc.

66{
67 G4Exception("G4VIntraNuclearTransportModel::Propagate()","G4VINT02",
69 "Propagate method for nucleus-nucleus interactions not implemented");
70 return nullptr;
71}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

References FatalException, and G4Exception().

Referenced by G4TheoFSGenerator::ApplyYourself().

◆ retryInelasticNucleus()

G4bool G4CascadeInterface::retryInelasticNucleus ( ) const
protected

Definition at line 740 of file G4CascadeInterface.cc.

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

References balance, bullet, coulombBarrierViolation(), G4cout, G4endl, G4InuclParticle::getDefinition(), G4CollisionOutput::getOutgoingParticles(), maximumTries, G4CollisionOutput::numberOfOutgoingNuclei(), G4CollisionOutput::numberOfOutgoingParticles(), numberOfTries, G4CascadeCheckBalance::okay(), and output.

Referenced by ApplyYourself(), and Propagate().

◆ retryInelasticProton()

G4bool G4CascadeInterface::retryInelasticProton ( ) const
protected

Definition at line 708 of file G4CascadeInterface.cc.

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

References bullet, G4cout, G4endl, G4InuclParticle::getDefinition(), G4CollisionOutput::getOutgoingParticles(), maximumTries, numberOfTries, and output.

Referenced by ApplyYourself().

◆ SampleInvariantT()

G4double G4HadronicInteraction::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
virtualinherited

◆ Set3DNucleus()

void G4VIntraNuclearTransportModel::Set3DNucleus ( G4V3DNucleus *const  value)
inlineinherited

Definition at line 124 of file G4VIntraNuclearTransportModel.hh.

125{
126 delete the3DNucleus; the3DNucleus = value;
127}

References G4VIntraNuclearTransportModel::the3DNucleus.

◆ SetDeExcitation()

void G4VIntraNuclearTransportModel::SetDeExcitation ( G4VPreCompoundModel ptr)
inlineinherited

◆ SetEnergyMomentumCheckLevels()

void G4HadronicInteraction::SetEnergyMomentumCheckLevels ( G4double  relativeLevel,
G4double  absoluteLevel 
)
inlineinherited

Definition at line 149 of file G4HadronicInteraction.hh.

150 { epCheckLevels.first = relativeLevel;
151 epCheckLevels.second = absoluteLevel; }

References G4HadronicInteraction::epCheckLevels.

Referenced by G4BinaryCascade::G4BinaryCascade(), G4CascadeInterface(), and G4FTFModel::G4FTFModel().

◆ SetMaxEnergy() [1/3]

void G4HadronicInteraction::SetMaxEnergy ( const G4double  anEnergy)
inlineinherited

Definition at line 102 of file G4HadronicInteraction.hh.

103 { theMaxEnergy = anEnergy; }

References G4HadronicInteraction::theMaxEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4IonINCLXXPhysics::AddProcess(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4HadronicBuilder::BuildFTFP_BERT(), G4HadronicBuilder::BuildFTFQGSP_BERT(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4HadronicBuilder::BuildQGSP_FTFP_BERT(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic::G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMaxEnergy() [2/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 151 of file G4HadronicInteraction.cc.

153{
154 Block();
155 if(!theMaxEnergyListElements.empty()) {
156 for(auto & elmlist : theMaxEnergyListElements) {
157 if( anElement == elmlist.second ) {
158 elmlist.first = anEnergy;
159 return;
160 }
161 }
162 }
163 theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
164}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyListElements.

◆ SetMaxEnergy() [3/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 166 of file G4HadronicInteraction.cc.

167{
168 Block();
169 if(!theMaxEnergyList.empty()) {
170 for(auto & matlist: theMaxEnergyList) {
171 if( aMaterial == matlist.second ) {
172 matlist.first = anEnergy;
173 return;
174 }
175 }
176 }
177 theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
178}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyList.

◆ SetMinEnergy() [1/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy)
inlineinherited

Definition at line 89 of file G4HadronicInteraction.hh.

90 { theMinEnergy = anEnergy; }

References G4HadronicInteraction::theMinEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPiKBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronElasticPhysicsHP::ConstructProcess(), G4HadronElasticPhysicsLEND::ConstructProcess(), G4HadronElasticPhysicsPHP::ConstructProcess(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionBertini::G4HadronicAbsorptionBertini(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic::G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4NeutrinoElectronCcModel::IsApplicable(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMinEnergy() [2/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 101 of file G4HadronicInteraction.cc.

103{
104 Block();
105 if(!theMinEnergyListElements.empty()) {
106 for(auto & elmlist : theMinEnergyListElements) {
107 if( anElement == elmlist.second ) {
108 elmlist.first = anEnergy;
109 return;
110 }
111 }
112 }
113 theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
114}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyListElements.

◆ SetMinEnergy() [3/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 116 of file G4HadronicInteraction.cc.

118{
119 Block();
120 if(!theMinEnergyList.empty()) {
121 for(auto & matlist : theMinEnergyList) {
122 if( aMaterial == matlist.second ) {
123 matlist.first = anEnergy;
124 return;
125 }
126 }
127 }
128 theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
129}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyList.

◆ SetModelName()

void G4HadronicInteraction::SetModelName ( const G4String nam)
inlineprotectedinherited

◆ SetPrimaryProjectile()

void G4VIntraNuclearTransportModel::SetPrimaryProjectile ( const G4HadProjectile aPrimary)
inlineinherited

Definition at line 148 of file G4VIntraNuclearTransportModel.hh.

149{
150 // NOTE: Previous pointer is NOT deleted: passed by reference, no ownership
151 thePrimaryProjectile = &aPrimary;
152}

References G4VIntraNuclearTransportModel::thePrimaryProjectile.

Referenced by G4TheoFSGenerator::ApplyYourself().

◆ SetRecoilEnergyThreshold()

void G4HadronicInteraction::SetRecoilEnergyThreshold ( G4double  val)
inlineinherited

◆ SetVerboseLevel()

void G4CascadeInterface::SetVerboseLevel ( G4int  verbose)

◆ throwNonConservationFailure()

void G4CascadeInterface::throwNonConservationFailure ( )
protected

Definition at line 789 of file G4CascadeInterface.cc.

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

References balance, G4CascadeCheckBalance::baryonOkay(), bullet, G4CascadeCheckBalance::chargeOkay(), G4CascadeCheckBalance::deltaB(), G4CascadeCheckBalance::deltaE(), G4CascadeCheckBalance::deltaP(), G4CascadeCheckBalance::deltaQ(), G4CascadeCheckBalance::energyOkay(), G4cerr, G4endl, G4CascadeCheckBalance::momentumOkay(), numberOfTries, output, G4CollisionOutput::printCollisionOutput(), G4CascadeCheckBalance::relativeE(), G4CascadeCheckBalance::relativeP(), and target.

Referenced by ApplyYourself(), and Propagate().

◆ useCascadeDeexcitation()

void G4CascadeInterface::useCascadeDeexcitation ( )

Definition at line 216 of file G4CascadeInterface.cc.

216 {
218}
void useCascadeDeexcitation()

References collider, and G4InuclCollider::useCascadeDeexcitation().

Referenced by G4CascadeInterface().

◆ usePreCompoundDeexcitation()

void G4CascadeInterface::usePreCompoundDeexcitation ( )

Field Documentation

◆ balance

G4CascadeCheckBalance* G4CascadeInterface::balance
private

◆ bullet

G4InuclParticle* G4CascadeInterface::bullet
private

◆ collider

G4InuclCollider* G4CascadeInterface::collider
private

◆ epCheckLevels

std::pair<G4double, G4double> G4HadronicInteraction::epCheckLevels
privateinherited

◆ hadronBullet

G4InuclElementaryParticle G4CascadeInterface::hadronBullet
private

Definition at line 177 of file G4CascadeInterface.hh.

Referenced by createBullet().

◆ hadronTarget

G4InuclElementaryParticle G4CascadeInterface::hadronTarget
private

Definition at line 180 of file G4CascadeInterface.hh.

Referenced by createTarget().

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ ltcollider

G4LightTargetCollider* G4CascadeInterface::ltcollider
private

Definition at line 170 of file G4CascadeInterface.hh.

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

◆ maximumTries

const G4int G4CascadeInterface::maximumTries
private

◆ nucleusBullet

G4InuclNuclei G4CascadeInterface::nucleusBullet
private

Definition at line 178 of file G4CascadeInterface.hh.

Referenced by createBullet().

◆ nucleusTarget

G4InuclNuclei G4CascadeInterface::nucleusTarget
private

Definition at line 181 of file G4CascadeInterface.hh.

Referenced by createTarget().

◆ numberOfTries

G4int G4CascadeInterface::numberOfTries
private

◆ output

G4CollisionOutput* G4CascadeInterface::output
private

◆ randomFile

const G4String G4CascadeInterface::randomFile
private

Definition at line 162 of file G4CascadeInterface.hh.

Referenced by ApplyYourself(), and Propagate().

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ secID

G4int G4CascadeInterface::secID
private

◆ target

G4InuclParticle* G4CascadeInterface::target
private

◆ the3DNucleus

G4V3DNucleus* G4VIntraNuclearTransportModel::the3DNucleus
protectedinherited

◆ theBlockedList

std::vector<const G4Material *> G4HadronicInteraction::theBlockedList
privateinherited

◆ theBlockedListElements

std::vector<const G4Element *> G4HadronicInteraction::theBlockedListElements
privateinherited

◆ theDeExcitation

G4VPreCompoundModel* G4VIntraNuclearTransportModel::theDeExcitation
protectedinherited

◆ theMaxEnergy

G4double G4HadronicInteraction::theMaxEnergy
protectedinherited

◆ theMaxEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMaxEnergyList
privateinherited

◆ theMaxEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMaxEnergyListElements
privateinherited

◆ theMinEnergy

G4double G4HadronicInteraction::theMinEnergy
protectedinherited

◆ theMinEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMinEnergyList
privateinherited

◆ theMinEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMinEnergyListElements
privateinherited

◆ theModelName

G4String G4HadronicInteraction::theModelName
privateinherited

◆ theParticleChange

G4HadFinalState G4HadronicInteraction::theParticleChange
protectedinherited

Definition at line 172 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4ElectroVDNuclearModel::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4NeutrinoElectronNcModel::ApplyYourself(), G4NeutronElectronElModel::ApplyYourself(), G4LFission::ApplyYourself(), G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4MuonVDNuclearModel::ApplyYourself(), G4NeutrinoElectronCcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4BinaryCascade::ApplyYourself(), ApplyYourself(), G4LMsdGenerator::ApplyYourself(), G4ElectroVDNuclearModel::CalculateEMVertex(), G4MuonVDNuclearModel::CalculateEMVertex(), G4ElectroVDNuclearModel::CalculateHadronicVertex(), G4MuonVDNuclearModel::CalculateHadronicVertex(), G4NeutrinoNucleusModel::CoherentPion(), copyOutputToHadronicResult(), G4BinaryCascade::DebugEpConservation(), G4BinaryCascade::DebugFinalEpConservation(), G4NeutrinoNucleusModel::FinalBarion(), G4NeutrinoNucleusModel::FinalMeson(), G4WilsonAbrasionModel::GetAbradedNucleons(), NoInteraction(), Propagate(), G4NeutrinoNucleusModel::RecoilDeexcitation(), G4LEHadronProtonElastic::~G4LEHadronProtonElastic(), G4LEnp::~G4LEnp(), and G4LFission::~G4LFission().

◆ thePrimaryProjectile

const G4HadProjectile* G4VIntraNuclearTransportModel::thePrimaryProjectile
protectedinherited

◆ theTransportModelName

G4String G4VIntraNuclearTransportModel::theTransportModelName
protectedinherited

◆ verboseLevel

G4int G4HadronicInteraction::verboseLevel
protectedinherited

Definition at line 177 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LFission::ApplyYourself(), G4MuMinusCapturePrecompound::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), ApplyYourself(), checkFinalResult(), copyOutputToHadronicResult(), copyOutputToReactionProducts(), G4LENDModel::create_used_target_map(), createBullet(), createTarget(), G4ElasticHadrNucleusHE::DefineHadronValues(), G4ElasticHadrNucleusHE::FillData(), G4ElasticHadrNucleusHE::FillFq2(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElasticHadrNucleusHE::G4ElasticHadrNucleusHE(), G4EMDissociation::G4EMDissociation(), G4hhElastic::G4hhElastic(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4ElasticHadrNucleusHE::GetFt(), G4ElasticHadrNucleusHE::GetLightFq2(), G4ElasticHadrNucleusHE::GetQ2_2(), G4HadronicInteraction::GetVerboseLevel(), G4ElasticHadrNucleusHE::HadronNucleusQ2_2(), G4ElasticHadrNucleusHE::HadronProtonQ2(), G4LFission::init(), G4DiffuseElastic::Initialise(), G4DiffuseElasticV2::Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4DiffuseElastic::InitialiseOnFly(), G4DiffuseElasticV2::InitialiseOnFly(), G4NuclNuclDiffuseElastic::InitialiseOnFly(), makeDynamicParticle(), NoInteraction(), Propagate(), G4ElasticHadrNucleusHE::SampleInvariantT(), G4AntiNuclElastic::SampleThetaCMS(), G4DiffuseElastic::SampleThetaLab(), G4NuclNuclDiffuseElastic::SampleThetaLab(), G4AntiNuclElastic::SampleThetaLab(), G4WilsonAbrasionModel::SetUseAblation(), G4HadronicInteraction::SetVerboseLevel(), G4WilsonAbrasionModel::SetVerboseLevel(), G4DiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElasticV2::ThetaCMStoThetaLab(), G4NuclNuclDiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElastic::ThetaLabToThetaCMS(), G4DiffuseElasticV2::ThetaLabToThetaCMS(), and G4NuclNuclDiffuseElastic::ThetaLabToThetaCMS().


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