Geant4 9.6.0
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 ()
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void DumpConfiguration (std::ostream &outFile) const
 
- Public Member Functions inherited from G4VIntraNuclearTransportModel
 G4VIntraNuclearTransportModel (const G4String &modelName="CascadeModel", G4VPreCompoundModel *ptr=0)
 
virtual ~G4VIntraNuclearTransportModel ()
 
virtual G4ReactionProductVectorPropagate (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)=0
 
void SetDeExcitation (G4VPreCompoundModel *ptr)
 
void Set3DNucleus (G4V3DNucleus *const value)
 
void SetPrimaryProjectile (const G4HadProjectile &aPrimary)
 
const G4StringGetModelName () const
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void PropagateModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

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

Constructor & Destructor Documentation

◆ G4CascadeInterface()

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

Definition at line 144 of file G4CascadeInterface.cc.

145 : G4VIntraNuclearTransportModel(name), numberOfTries(0),
146 collider(new G4InuclCollider), balance(new G4CascadeCheckBalance(name)),
147 bullet(0), target(0), output(new G4CollisionOutput) {
148 SetEnergyMomentumCheckLevels(5*perCent, 10*MeV);
149 balance->setLimits(5*perCent, 10*MeV/GeV); // Bertini internal units
150 // Description(); // Model description
151
154}
void setLimits(G4double relative, G4double absolute)
void SetVerboseLevel(G4int verbose)
static G4bool usePreCompound()
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)

◆ ~G4CascadeInterface()

G4CascadeInterface::~G4CascadeInterface ( )
virtual

Definition at line 157 of file G4CascadeInterface.cc.

157 {
158 clear();
159 delete collider; collider=0;
160 delete balance; balance=0;
161 delete output; output=0;
162}

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 230 of file G4CascadeInterface.cc.

231 {
232 if (verboseLevel)
233 G4cout << " >>> G4CascadeInterface::ApplyYourself" << G4endl;
234
235 if (aTrack.GetKineticEnergy() < 0.) {
236 G4cerr << " >>> G4CascadeInterface got negative-energy track: "
237 << aTrack.GetDefinition()->GetParticleName() << " Ekin = "
238 << aTrack.GetKineticEnergy() << G4endl;
239 }
240
241#ifdef G4CASCADE_DEBUG_INTERFACE
242 static G4int counter(0);
243 counter++;
244 G4cerr << "Reaction number "<< counter << " "
245 << aTrack.GetDefinition()->GetParticleName() << " "
246 << aTrack.GetKineticEnergy() << " MeV" << G4endl;
247#endif
248
249 if (!randomFile.empty()) { // User requested random-seed capture
250 if (verboseLevel>1)
251 G4cout << " Saving random engine state to " << randomFile << G4endl;
253 }
254
256 clear();
257
258 // Abort processing if no interaction is possible
259 if (!IsApplicable(aTrack, theNucleus)) {
260 if (verboseLevel) G4cerr << " No interaction possible " << G4endl;
261 return NoInteraction(aTrack, theNucleus);
262 }
263
264 // Make conversion between native Geant4 and Bertini cascade classes.
265 if (!createBullet(aTrack)) {
266 if (verboseLevel) G4cerr << " Unable to create usable bullet" << G4endl;
267 return NoInteraction(aTrack, theNucleus);
268 }
269
270 if (!createTarget(theNucleus)) {
271 if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
272 return NoInteraction(aTrack, theNucleus);
273 }
274
275 // Different retry conditions for proton target vs. nucleus
276 const G4bool isHydrogen = (theNucleus.GetA_asInt() == 1);
277
278 numberOfTries = 0;
279 do { // we try to create inelastic interaction
280 if (verboseLevel > 1)
281 G4cout << " Generating cascade attempt " << numberOfTries << G4endl;
282
283 output->reset();
284 collider->collide(bullet, target, *output);
285 balance->collide(bullet, target, *output);
286
287 numberOfTries++;
288 } while ( isHydrogen ? retryInelasticProton() : retryInelasticNucleus() );
289
290 // Null event if unsuccessful
291 if (numberOfTries >= maximumTries) {
292 if (verboseLevel)
293 G4cout << " Cascade aborted after trials " << numberOfTries << G4endl;
294 return NoInteraction(aTrack, theNucleus);
295 }
296
297 // Abort job if energy or momentum are not conserved
298 if (!balance->okay()) {
300 return NoInteraction(aTrack, theNucleus);
301 }
302
303 // Successful cascade -- clean up and return
304 if (verboseLevel) {
305 G4cout << " Cascade output after trials " << numberOfTries << G4endl;
306 if (verboseLevel > 1) output->printCollisionOutput();
307 }
308
309 // Rotate event to put Z axis along original projectile direction
310 output->rotateEvent(bulletInLabFrame);
311
313
314 // Report violations of conservation laws in original frame
316
317 // Clean up and return final result;
318 clear();
319 return &theParticleChange;
320}
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
static void saveEngineStatus(const char filename[]="Config.conf")
Definition: Random.cc:175
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 rotateEvent(const G4LorentzRotation &rotate)
void printCollisionOutput(std::ostream &os=G4cout) const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
const G4String & GetParticleName() const

◆ checkFinalResult()

void G4CascadeInterface::checkFinalResult ( )
protected

Definition at line 598 of file G4CascadeInterface.cc.

598 {
599 balance->collide(bullet, target, *output);
600
601 if (verboseLevel > 2) {
602 if (!balance->baryonOkay()) {
603 G4cerr << "ERROR: no baryon number conservation, sum of baryons = "
604 << balance->deltaB() << G4endl;
605 }
606
607 if (!balance->chargeOkay()) {
608 G4cerr << "ERROR: no charge conservation, sum of charges = "
609 << balance->deltaQ() << G4endl;
610 }
611
612 if (std::abs(balance->deltaKE()) > 0.01 ) { // GeV
613 G4cerr << "Kinetic energy conservation violated by "
614 << balance->deltaKE() << " GeV" << G4endl;
615 }
616
617 G4double eInit = bullet->getEnergy() + target->getEnergy();
618 G4double eFinal = eInit + balance->deltaE();
619
620 G4cout << "Initial energy " << eInit << " final energy " << eFinal
621 << "\nTotal energy conservation at level "
622 << balance->deltaE() * GeV << " MeV" << G4endl;
623
624 if (balance->deltaKE() > 5.0e-5 ) { // 0.05 MeV
625 G4cerr << "FATAL ERROR: kinetic energy created "
626 << balance->deltaKE() * GeV << " MeV" << G4endl;
627 }
628 }
629}
double G4double
Definition: G4Types.hh:64
G4double getEnergy() const

Referenced by ApplyYourself().

◆ clear()

void G4CascadeInterface::clear ( )
protected

Definition at line 183 of file G4CascadeInterface.cc.

183 {
184 bullet=0;
185 target=0;
186}

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

◆ copyOutputToHadronicResult()

void G4CascadeInterface::copyOutputToHadronicResult ( )
protected

Definition at line 529 of file G4CascadeInterface.cc.

529 {
530 if (verboseLevel > 1)
531 G4cout << " >>> G4CascadeInterface::copyOutputToHadronicResult" << G4endl;
532
533 const std::vector<G4InuclNuclei>& outgoingNuclei = output->getOutgoingNuclei();
534 const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
535
538
539 // Get outcoming particles
540 if (!particles.empty()) {
541 particleIterator ipart = particles.begin();
542 for (; ipart != particles.end(); ipart++) {
544 }
545 }
546
547 // get nuclei fragments
548 if (!outgoingNuclei.empty()) {
549 nucleiIterator ifrag = outgoingNuclei.begin();
550 for (; ifrag != outgoingNuclei.end(); ifrag++) {
552 }
553 }
554}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:60
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)
void SetEnergyChange(G4double anEnergy)

Referenced by ApplyYourself().

◆ copyOutputToReactionProducts()

G4ReactionProductVector * G4CascadeInterface::copyOutputToReactionProducts ( )
protected

Definition at line 556 of file G4CascadeInterface.cc.

556 {
557 if (verboseLevel > 1)
558 G4cout << " >>> G4CascadeInterface::copyOutputToReactionProducts" << G4endl;
559
560 const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
561 const std::vector<G4InuclNuclei>& fragments = output->getOutgoingNuclei();
562
564
565 G4ReactionProduct* rp = 0; // Buffers to create outgoing tracks
566 G4DynamicParticle* dp = 0;
567
568 // Get outcoming particles
569 if (!particles.empty()) {
570 particleIterator ipart = particles.begin();
571 for (; ipart != particles.end(); ipart++) {
572 rp = new G4ReactionProduct;
573 dp = makeDynamicParticle(*ipart);
574 (*rp) = (*dp); // This does all the necessary copying
575 propResult->push_back(rp);
576 delete dp;
577 }
578 }
579
580 // get nuclei fragments
581 if (!fragments.empty()) {
582 nucleiIterator ifrag = fragments.begin();
583 for (; ifrag != fragments.end(); ifrag++) {
584 rp = new G4ReactionProduct;
585 dp = makeDynamicParticle(*ifrag);
586 (*rp) = (*dp); // This does all the necessary copying
587 propResult->push_back(rp);
588 delete dp;
589 }
590 }
591
592 return propResult;
593}
std::vector< G4ReactionProduct * > G4ReactionProductVector

Referenced by Propagate().

◆ coulombBarrierViolation()

G4bool G4CascadeInterface::coulombBarrierViolation ( ) const
protected

Definition at line 634 of file G4CascadeInterface.cc.

634 {
635 G4bool violated = false; // by default coulomb analysis is OK
636
637 const G4double coulumbBarrier = 8.7 * MeV/GeV; // Bertini uses GeV
638
639 const std::vector<G4InuclElementaryParticle>& p =
640 output->getOutgoingParticles();
641
642 for (particleIterator ipart=p.begin(); ipart != p.end(); ipart++) {
643 if (ipart->type() == proton) {
644 violated |= (ipart->getKineticEnergy() < coulumbBarrier);
645 }
646 }
647
648 return violated;
649}

Referenced by retryInelasticNucleus().

◆ createBullet()

G4bool G4CascadeInterface::createBullet ( const G4HadProjectile aTrack)
protected

Definition at line 415 of file G4CascadeInterface.cc.

415 {
416 const G4ParticleDefinition* trkDef = aTrack.GetDefinition();
417
418 G4int bulletType = 0; // For elementary particles
419 G4int bulletA = 0, bulletZ = 0; // For nucleus projectile
420
421 if (trkDef->GetAtomicMass() <= 1) {
422 bulletType = G4InuclElementaryParticle::type(trkDef);
423 } else {
424 bulletA = trkDef->GetAtomicMass();
425 bulletZ = trkDef->GetAtomicNumber();
426 }
427
428 if (0 == bulletType && 0 == bulletA*bulletZ) {
429 if (verboseLevel) {
430 G4cerr << " G4CascadeInterface: " << trkDef->GetParticleName()
431 << " not usable as bullet." << G4endl;
432 }
433 bullet = 0;
434 return false;
435 }
436
437 // Code momentum and energy -- Bertini wants z-axis and GeV units
438 G4LorentzVector projectileMomentum = aTrack.Get4Momentum()/GeV;
439
440 // Rotation/boost to get from z-axis back to original frame
441 bulletInLabFrame = G4LorentzRotation::IDENTITY; // Initialize
442 bulletInLabFrame.rotateZ(-projectileMomentum.phi());
443 bulletInLabFrame.rotateY(-projectileMomentum.theta());
444 bulletInLabFrame.invert();
445
446 G4LorentzVector momentumBullet(0., 0., projectileMomentum.rho(),
447 projectileMomentum.e());
448
449 if (bulletType > 0) {
450 hadronBullet.fill(momentumBullet, bulletType);
451 bullet = &hadronBullet;
452 } else {
453 nucleusBullet.fill(momentumBullet, bulletA, bulletZ);
454 bullet = &nucleusBullet;
455 }
456
457 if (verboseLevel > 2) G4cout << "Bullet: \n" << *bullet << G4endl;
458
459 return true;
460}
static DLL_API const HepLorentzRotation IDENTITY
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation & invert()
double theta() const
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 473 of file G4CascadeInterface.cc.

473 {
474 if (A > 1) {
475 nucleusTarget.fill(A, Z);
476 target = &nucleusTarget;
477 } else {
478 hadronTarget.fill(0., (Z=1?proton:neutron));
479 target = &hadronTarget;
480 }
481
482 if (verboseLevel > 2) G4cout << "Target: \n" << *target << G4endl;
483
484 return true; // Right now, target never fails
485}
@ neutron

◆ createTarget() [2/3]

G4bool G4CascadeInterface::createTarget ( G4Nucleus theNucleus)
protected

Definition at line 465 of file G4CascadeInterface.cc.

465 {
466 return createTarget(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt());
467}
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115

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

◆ createTarget() [3/3]

G4bool G4CascadeInterface::createTarget ( G4V3DNucleus theNucleus)
protected

Definition at line 469 of file G4CascadeInterface.cc.

469 {
470 return createTarget(theNucleus->GetMassNumber(), theNucleus->GetCharge());
471}
virtual G4int GetCharge()=0
virtual G4int GetMassNumber()=0

◆ DumpConfiguration()

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

Definition at line 178 of file G4CascadeInterface.cc.

178 {
180}
static void DumpConfiguration(std::ostream &os)

◆ IsApplicable() [1/2]

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

Reimplemented from G4HadronicInteraction.

Definition at line 213 of file G4CascadeInterface.cc.

214 {
215 return IsApplicable(aTrack.GetDefinition());
216}

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

◆ IsApplicable() [2/2]

G4bool G4CascadeInterface::IsApplicable ( const G4ParticleDefinition aPD) const

Definition at line 218 of file G4CascadeInterface.cc.

218 {
219 if (aPD->GetAtomicMass() > 1) return true; // Nuclei are okay
220
221 // Valid particle and have interactions available
223 return (type>0 && G4CascadeChannelTables::GetTable(type));
224}
static const G4CascadeChannel * GetTable(G4int initialState)

◆ makeDynamicParticle() [1/2]

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

Definition at line 490 of file G4CascadeInterface.cc.

491 {
492 G4int outgoingType = iep.type();
493
494 if (iep.quasi_deutron()) {
495 G4cerr << " ERROR: G4CascadeInterface incompatible particle type "
496 << outgoingType << G4endl;
497 return 0;
498 }
499
500 // Copy local G4DynPart to public output (handle kaon mixing specially)
501 if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
502 G4ThreeVector momDir = iep.getMomentum().vect().unit();
503 G4double ekin = iep.getKineticEnergy()*GeV; // Bertini -> G4 units
504
506 if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
507
508 return new G4DynamicParticle(pd, momDir, ekin);
509 } else {
510 return new G4DynamicParticle(iep.getDynamicParticle());
511 }
512
513 return 0; // Should never get here!
514}
#define G4UniformRand()
Definition: Randomize.hh:53
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 517 of file G4CascadeInterface.cc.

517 {
518 if (verboseLevel > 2) {
519 G4cout << " Nuclei fragment: \n" << inuc << G4endl;
520 }
521
522 // Copy local G4DynPart to public output
523 return new G4DynamicParticle(inuc.getDynamicParticle());
524}

◆ ModelDescription()

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

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 164 of file G4CascadeInterface.cc.

165{
166 outFile << "The Bertini-style cascade implements the inelastic scattering\n"
167 << "of hadrons by nuclei. Nucleons, pions, kaons and hyperons\n"
168 << "from 0 to 15 GeV may be used as projectiles in this model.\n"
169 << "Final state hadrons are produced by a classical cascade\n"
170 << "consisting of individual hadron-nucleon scatterings which use\n"
171 << "free-space partial cross sections, corrected for various\n"
172 << "nuclear medium effects. The target nucleus is modeled as a\n"
173 << "set of 1, 3 or 6 spherical shells, in which scattered hadrons\n"
174 << "travel in straight lines until they are reflected from or\n"
175 << "transmitted through shell boundaries.\n";
176}

◆ NoInteraction()

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

Definition at line 398 of file G4CascadeInterface.cc.

399 {
400 if (verboseLevel)
401 G4cout << " >>> G4CascadeInterface::NoInteraction" << G4endl;
402
405
406 G4double ekin = aTrack.GetKineticEnergy()>0. ? aTrack.GetKineticEnergy() : 0.;
407 theParticleChange.SetEnergyChange(ekin); // Protect against rounding
408
409 return &theParticleChange;
410}
@ isAlive

Referenced by ApplyYourself().

◆ Propagate()

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

Implements G4VIntraNuclearTransportModel.

Definition at line 323 of file G4CascadeInterface.cc.

324 {
325 if (verboseLevel) G4cout << " >>> G4CascadeInterface::Propagate" << G4endl;
326
327#ifdef G4CASCADE_DEBUG_INTERFACE
328 if (verboseLevel>1) {
329 G4cout << " G4V3DNucleus A " << theNucleus->GetMassNumber()
330 << " Z " << theNucleus->GetCharge()
331 << "\n " << theSecondaries->size() << " secondaries:" << G4endl;
332 for (size_t i=0; i<theSecondaries->size(); i++) {
333 G4KineticTrack* kt = (*theSecondaries)[i];
334 G4cout << " " << i << ": " << kt->GetDefinition()->GetParticleName()
335 << " p " << kt->Get4Momentum() << " @ " << kt->GetPosition()
336 << " t " << kt->GetFormationTime() << G4endl;
337 }
338 }
339#endif
340
341 if (!randomFile.empty()) { // User requested random-seed capture
342 if (verboseLevel>1)
343 G4cout << " Saving random engine state to " << randomFile << G4endl;
345 }
346
348 clear();
349
350 // Process input secondaries list to eliminate resonances
351 G4DecayKineticTracks decay(theSecondaries);
352
353 // NOTE: Requires 9.4-ref-03 mods to base class and G4TheoFSGenerator
354 const G4HadProjectile* projectile = GetPrimaryProjectile();
355 if (projectile) createBullet(*projectile);
356
357 if (!createTarget(theNucleus)) {
358 if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
359 return 0; // FIXME: This will cause a segfault later
360 }
361
362 numberOfTries = 0;
363 do {
364 if (verboseLevel > 1)
365 G4cout << " Generating rescatter attempt " << numberOfTries << G4endl;
366
367 output->reset();
368 collider->rescatter(bullet, theSecondaries, theNucleus, *output);
369 balance->collide(bullet, target, *output);
370
371 numberOfTries++;
372 // FIXME: retry checks will SEGFAULT until we can define the bullet!
373 } while (retryInelasticNucleus());
374
375 // Check whether repeated attempts have all failed; report and exit
376 if (numberOfTries >= maximumTries && !balance->okay()) {
377 throwNonConservationFailure(); // This terminates the job
378 }
379
380 // Successful cascade -- clean up and return
381 if (verboseLevel) {
382 G4cout << " Cascade rescatter after trials " << numberOfTries << G4endl;
383 if (verboseLevel > 1) output->printCollisionOutput();
384 }
385
386 // Does calling code take ownership? I hope so!
388
389 // Clean up and and return final result
390 clear();
391 return propResult;
392}
G4ReactionProductVector * copyOutputToReactionProducts()
void rescatter(G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
G4double GetFormationTime() const
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
const G4HadProjectile * GetPrimaryProjectile() const

◆ retryInelasticNucleus()

G4bool G4CascadeInterface::retryInelasticNucleus ( ) const
protected

Definition at line 685 of file G4CascadeInterface.cc.

685 {
686 // Quantities necessary for conditional block below
687 G4int npart = output->numberOfOutgoingParticles();
688 G4int nfrag = output->numberOfOutgoingNuclei();
689
690 const G4ParticleDefinition* firstOut = (npart == 0) ? 0 :
691 output->getOutgoingParticles().begin()->getDefinition();
692
693#ifdef G4CASCADE_DEBUG_INTERFACE
694 // Report on all retry conditions, in order of return logic
695 G4cout << " retryInelasticNucleus: numberOfTries "
696 << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
697 << "\n retryInelasticNucleus: AND outputParticles "
698 << ((npart != 0) ? "NON-ZERO (t)" : "EMPTY (f)")
699#ifdef G4CASCADE_COULOMB_DEV
700 << "\n retryInelasticNucleus: AND coulombBarrier (COULOMB_DEV) "
701 << (coulombBarrierViolation() ? "VIOLATED (t)" : "PASSED (f)")
702 << "\n retryInelasticNucleus: AND collision type (COULOMB_DEV) "
703 << ((npart+nfrag > 2) ? "INELASTIC (t)" : "ELASTIC (f)")
704#else
705 << "\n retryInelasticNucleus: AND collsion type "
706 << ((npart+nfrag < 3) ? "ELASTIC (t)" : "INELASTIC (f)")
707 << "\n retryInelasticNucleus: AND Leading particle bullet "
708 << ((firstOut == bullet->getDefinition()) ? "YES (t)" : "NO (f)")
709#endif
710 << "\n retryInelasticNucleus: OR conservation "
711 << (!balance->okay() ? "FAILED (t)" : "PASSED (f)")
712 << G4endl;
713#endif
714
715 return ( ((numberOfTries < maximumTries) &&
716 (npart != 0) &&
717#ifdef G4CASCADE_COULOMB_DEV
718 (coulombBarrierViolation() && npart+nfrag > 2)
719#else
720 (npart+nfrag < 3 && firstOut == bullet->getDefinition())
721#endif
722 )
723#ifndef G4CASCADE_SKIP_ECONS
724 || (!balance->okay())
725#endif
726 );
727}
G4bool coulombBarrierViolation() const
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
G4ParticleDefinition * getDefinition() const

Referenced by ApplyYourself(), and Propagate().

◆ retryInelasticProton()

G4bool G4CascadeInterface::retryInelasticProton ( ) const
protected

Definition at line 653 of file G4CascadeInterface.cc.

653 {
654 const std::vector<G4InuclElementaryParticle>& out =
655 output->getOutgoingParticles();
656
657#ifdef G4CASCADE_DEBUG_INTERFACE
658 // Report on all retry conditions, in order of return logic
659 G4cout << " retryInelasticProton: number of Tries "
660 << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
661 << "\n retryInelasticProton: AND collision type ";
662 if (out.empty()) G4cout << "FAILED" << G4endl;
663 else {
664 G4cout << (out.size() == 2 ? "ELASTIC (t)" : "INELASTIC (f)")
665 << "\n retryInelasticProton: AND Leading particles bullet "
666 << (out.size() >= 2 &&
667 (out[0].getDefinition() == bullet->getDefinition() ||
668 out[1].getDefinition() == bullet->getDefinition())
669 ? "YES (t)" : "NO (f)")
670 << G4endl;
671 }
672#endif
673
674 return ( (numberOfTries < maximumTries) &&
675 (out.empty() ||
676 (out.size() == 2 &&
677 (out[0].getDefinition() == bullet->getDefinition() ||
678 out[1].getDefinition() == bullet->getDefinition())))
679 );
680}

Referenced by ApplyYourself().

◆ SetVerboseLevel()

void G4CascadeInterface::SetVerboseLevel ( G4int  verbose)

Definition at line 203 of file G4CascadeInterface.cc.

203 {
205 collider->setVerboseLevel(verbose);
206 balance->setVerboseLevel(verbose);
207 output->setVerboseLevel(verbose);
208}
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 733 of file G4CascadeInterface.cc.

733 {
734 // NOTE: Once G4HadronicException is changed, use the following line!
735 // G4ExceptionDescription errInfo;
736 std::ostream& errInfo = G4cerr;
737
738 errInfo << " >>> G4CascadeInterface has non-conserving"
739 << " cascade after " << numberOfTries << " attempts." << G4endl;
740
741 G4String throwMsg = "G4CascadeInterface - ";
742 if (!balance->energyOkay()) {
743 throwMsg += "Energy";
744 errInfo << " Energy conservation violated by " << balance->deltaE()
745 << " GeV (" << balance->relativeE() << ")" << G4endl;
746 }
747
748 if (!balance->momentumOkay()) {
749 throwMsg += "Momentum";
750 errInfo << " Momentum conservation violated by " << balance->deltaP()
751 << " GeV/c (" << balance->relativeP() << ")" << G4endl;
752 }
753
754 if (!balance->baryonOkay()) {
755 throwMsg += "Baryon number";
756 errInfo << " Baryon number violated by " << balance->deltaB() << G4endl;
757 }
758
759 if (!balance->chargeOkay()) {
760 throwMsg += "Charge";
761 errInfo << " Charge conservation violated by " << balance->deltaQ()
762 << G4endl;
763 }
764
765 errInfo << " Final event output, for debugging:\n"
766 << " Bullet: \n" << *bullet << G4endl
767 << " Target: \n" << *target << G4endl;
768 output->printCollisionOutput(errInfo);
769
770 throwMsg += " non-conservation. More info in output.";
771 throw G4HadronicException(__FILE__, __LINE__, throwMsg); // Job ends here!
772}

Referenced by ApplyYourself(), and Propagate().

◆ useCascadeDeexcitation()

void G4CascadeInterface::useCascadeDeexcitation ( )

Definition at line 192 of file G4CascadeInterface.cc.

192 {
193 collider->useCascadeDeexcitation();
194}
void useCascadeDeexcitation()

◆ usePreCompoundDeexcitation()

void G4CascadeInterface::usePreCompoundDeexcitation ( )

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