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

#include <G4CascadeInterface.hh>

Inheritance diagram for G4CascadeInterface:
G4VIntraNuclearTransportModel G4HadronicInteraction

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 G4ReactionProductVectorPropagateNuclNucl (G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
 
void SetDeExcitation (G4VPreCompoundModel *ptr)
 
void Set3DNucleus (G4V3DNucleus *const value)
 
void SetPrimaryProjectile (const G4HadProjectile &aPrimary)
 
const G4StringGetModelName () const
 
virtual void PropagateModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 

Static Public Member Functions

static void Initialize ()
 

Protected Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from G4VIntraNuclearTransportModel
G4String theTransportModelName
 
G4V3DNucleusthe3DNucleus
 
G4VPreCompoundModeltheDeExcitation
 
const G4HadProjectilethePrimaryProjectile
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 86 of file G4CascadeInterface.hh.

Constructor & Destructor Documentation

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

Definition at line 149 of file G4CascadeInterface.cc.

References python.hepunit::GeV, python.hepunit::MeV, python.hepunit::perCent, G4HadronicInteraction::SetEnergyMomentumCheckLevels(), G4CascadeCheckBalance::setLimits(), SetVerboseLevel(), G4CascadeParameters::usePreCompound(), usePreCompoundDeexcitation(), and G4CascadeParameters::verbose().

150  : G4VIntraNuclearTransportModel(name), numberOfTries(0),
151  collider(new G4InuclCollider), balance(new G4CascadeCheckBalance(name)),
152  bullet(0), target(0), output(new G4CollisionOutput) {
154  balance->setLimits(5*perCent, 10*MeV/GeV); // Bertini internal units
155  // Description(); // Model description
156 
159 }
G4VIntraNuclearTransportModel(const G4String &modelName="CascadeModel", G4VPreCompoundModel *ptr=0)
const XML_Char * target
void setLimits(G4double relative, G4double absolute)
float perCent
Definition: hepunit.py:239
void SetVerboseLevel(G4int verbose)
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
static G4bool usePreCompound()
G4CascadeInterface::~G4CascadeInterface ( )
virtual

Definition at line 161 of file G4CascadeInterface.cc.

References clear().

161  {
162  clear();
163  delete collider; collider=0;
164  delete balance; balance=0;
165  delete output; output=0;
166 }

Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 244 of file G4CascadeInterface.cc.

References checkFinalResult(), G4HadFinalState::Clear(), clear(), G4InuclCollider::collide(), G4CascadeCheckBalance::collide(), copyOutputToHadronicResult(), createBullet(), createTarget(), G4cerr, G4cout, G4endl, G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4HadProjectile::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), IsApplicable(), NoInteraction(), G4CascadeCheckBalance::okay(), G4CollisionOutput::printCollisionOutput(), G4CollisionOutput::reset(), retryInelasticNucleus(), retryInelasticProton(), G4CollisionOutput::rotateEvent(), CLHEP::HepRandom::saveEngineStatus(), G4HadronicInteraction::theParticleChange, throwNonConservationFailure(), and G4HadronicInteraction::verboseLevel.

245  {
246  if (verboseLevel)
247  G4cout << " >>> G4CascadeInterface::ApplyYourself" << G4endl;
248 
249  if (aTrack.GetKineticEnergy() < 0.) {
250  G4cerr << " >>> G4CascadeInterface got negative-energy track: "
251  << aTrack.GetDefinition()->GetParticleName() << " Ekin = "
252  << aTrack.GetKineticEnergy() << G4endl;
253  }
254 
255 #ifdef G4CASCADE_DEBUG_INTERFACE
256  static G4int counter(0);
257  counter++;
258  G4cerr << "Reaction number "<< counter << " "
259  << aTrack.GetDefinition()->GetParticleName() << " "
260  << aTrack.GetKineticEnergy() << " MeV" << G4endl;
261 #endif
262 
263  if (!randomFile.empty()) { // User requested random-seed capture
264  if (verboseLevel>1)
265  G4cout << " Saving random engine state to " << randomFile << G4endl;
267  }
268 
270  clear();
271 
272  // Abort processing if no interaction is possible
273  if (!IsApplicable(aTrack, theNucleus)) {
274  if (verboseLevel) G4cerr << " No interaction possible " << G4endl;
275  return NoInteraction(aTrack, theNucleus);
276  }
277 
278  // Make conversion between native Geant4 and Bertini cascade classes.
279  if (!createBullet(aTrack)) {
280  if (verboseLevel) G4cerr << " Unable to create usable bullet" << G4endl;
281  return NoInteraction(aTrack, theNucleus);
282  }
283 
284  if (!createTarget(theNucleus)) {
285  if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
286  return NoInteraction(aTrack, theNucleus);
287  }
288 
289  // Different retry conditions for proton target vs. nucleus
290  const G4bool isHydrogen = (theNucleus.GetA_asInt() == 1);
291 
292  numberOfTries = 0;
293  do { // we try to create inelastic interaction
294  if (verboseLevel > 1)
295  G4cout << " Generating cascade attempt " << numberOfTries << G4endl;
296 
297  output->reset();
298  collider->collide(bullet, target, *output);
299  balance->collide(bullet, target, *output);
300 
301  numberOfTries++;
302  } while ( isHydrogen ? retryInelasticProton() : retryInelasticNucleus() );
303 
304  // Null event if unsuccessful
305  if (numberOfTries >= maximumTries) {
306  if (verboseLevel)
307  G4cout << " Cascade aborted after trials " << numberOfTries << G4endl;
308  return NoInteraction(aTrack, theNucleus);
309  }
310 
311  // Abort job if energy or momentum are not conserved
312  if (!balance->okay()) {
314  return NoInteraction(aTrack, theNucleus);
315  }
316 
317  // Successful cascade -- clean up and return
318  if (verboseLevel) {
319  G4cout << " Cascade output after trials " << numberOfTries << G4endl;
320  if (verboseLevel > 1) output->printCollisionOutput();
321  }
322 
323  // Rotate event to put Z axis along original projectile direction
324  output->rotateEvent(bulletInLabFrame);
325 
327 
328  // Report violations of conservation laws in original frame
330 
331  // Clean up and return final result;
332  clear();
333 /*
334  G4int nSec = theParticleChange.GetNumberOfSecondaries();
335  for (G4int i = 0; i < nSec; i++) {
336  G4HadSecondary* sec = theParticleChange.GetSecondary(i);
337  G4DynamicParticle* dp = sec->GetParticle();
338  if (dp->GetDefinition()->GetParticleName() == "neutron")
339  G4cout << dp->GetDefinition()->GetParticleName() << " has "
340  << dp->GetKineticEnergy()/MeV << " MeV " << G4endl;
341  }
342 */
343  return &theParticleChange;
344 }
G4bool retryInelasticNucleus() const
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
void printCollisionOutput(std::ostream &os=G4cout) const
const XML_Char * target
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
G4bool createTarget(G4Nucleus &theNucleus)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
void rotateEvent(const G4LorentzRotation &rotate)
static void saveEngineStatus(const char filename[]="Config.conf")
Definition: Random.cc:175
#define G4endl
Definition: G4ios.hh:61
G4bool createBullet(const G4HadProjectile &aTrack)
G4bool retryInelasticProton() const
G4GLOB_DLL std::ostream G4cerr
G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4HadFinalState * NoInteraction(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
void G4CascadeInterface::checkFinalResult ( )
protected

Definition at line 621 of file G4CascadeInterface.cc.

References G4CascadeCheckBalance::baryonOkay(), G4CascadeCheckBalance::chargeOkay(), G4CascadeCheckBalance::collide(), G4CascadeCheckBalance::deltaB(), G4CascadeCheckBalance::deltaE(), G4CascadeCheckBalance::deltaKE(), G4CascadeCheckBalance::deltaQ(), G4cerr, G4cout, G4endl, G4InuclParticle::getEnergy(), python.hepunit::GeV, and G4HadronicInteraction::verboseLevel.

Referenced by ApplyYourself().

621  {
622  balance->collide(bullet, target, *output);
623 
624  if (verboseLevel > 2) {
625  if (!balance->baryonOkay()) {
626  G4cerr << "ERROR: no baryon number conservation, sum of baryons = "
627  << balance->deltaB() << G4endl;
628  }
629 
630  if (!balance->chargeOkay()) {
631  G4cerr << "ERROR: no charge conservation, sum of charges = "
632  << balance->deltaQ() << G4endl;
633  }
634 
635  if (std::abs(balance->deltaKE()) > 0.01 ) { // GeV
636  G4cerr << "Kinetic energy conservation violated by "
637  << balance->deltaKE() << " GeV" << G4endl;
638  }
639 
640  G4double eInit = bullet->getEnergy() + target->getEnergy();
641  G4double eFinal = eInit + balance->deltaE();
642 
643  G4cout << "Initial energy " << eInit << " final energy " << eFinal
644  << "\nTotal energy conservation at level "
645  << balance->deltaE() * GeV << " MeV" << G4endl;
646 
647  if (balance->deltaKE() > 5.0e-5 ) { // 0.05 MeV
648  G4cerr << "FATAL ERROR: kinetic energy created "
649  << balance->deltaKE() * GeV << " MeV" << G4endl;
650  }
651  }
652 }
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4double getEnergy() const
const XML_Char * target
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr
void G4CascadeInterface::clear ( void  )
protected

Definition at line 186 of file G4CascadeInterface.cc.

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

186  {
187  bullet=0;
188  target=0;
189 }
const XML_Char * target
void G4CascadeInterface::copyOutputToHadronicResult ( )
protected

Definition at line 552 of file G4CascadeInterface.cc.

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

Referenced by ApplyYourself().

552  {
553  if (verboseLevel > 1)
554  G4cout << " >>> G4CascadeInterface::copyOutputToHadronicResult" << G4endl;
555 
556  const std::vector<G4InuclNuclei>& outgoingNuclei = output->getOutgoingNuclei();
557  const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
558 
561 
562  // Get outcoming particles
563  if (!particles.empty()) {
564  particleIterator ipart = particles.begin();
565  for (; ipart != particles.end(); ipart++) {
567  }
568  }
569 
570  // get nuclei fragments
571  if (!outgoingNuclei.empty()) {
572  nucleiIterator ifrag = outgoingNuclei.begin();
573  for (; ifrag != outgoingNuclei.end(); ifrag++) {
575  }
576  }
577 }
void SetStatusChange(G4HadFinalStateStatus aS)
G4GLOB_DLL std::ostream G4cout
G4DynamicParticle * makeDynamicParticle(const G4InuclElementaryParticle &iep) const
void SetEnergyChange(G4double anEnergy)
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
#define G4endl
Definition: G4ios.hh:61
void AddSecondary(G4DynamicParticle *aP)
G4ReactionProductVector * G4CascadeInterface::copyOutputToReactionProducts ( )
protected

Definition at line 579 of file G4CascadeInterface.cc.

References G4cout, G4endl, G4CollisionOutput::getOutgoingNuclei(), G4CollisionOutput::getOutgoingParticles(), makeDynamicParticle(), and G4HadronicInteraction::verboseLevel.

Referenced by Propagate().

579  {
580  if (verboseLevel > 1)
581  G4cout << " >>> G4CascadeInterface::copyOutputToReactionProducts" << G4endl;
582 
583  const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
584  const std::vector<G4InuclNuclei>& fragments = output->getOutgoingNuclei();
585 
587 
588  G4ReactionProduct* rp = 0; // Buffers to create outgoing tracks
589  G4DynamicParticle* dp = 0;
590 
591  // Get outcoming particles
592  if (!particles.empty()) {
593  particleIterator ipart = particles.begin();
594  for (; ipart != particles.end(); ipart++) {
595  rp = new G4ReactionProduct;
596  dp = makeDynamicParticle(*ipart);
597  (*rp) = (*dp); // This does all the necessary copying
598  propResult->push_back(rp);
599  delete dp;
600  }
601  }
602 
603  // get nuclei fragments
604  if (!fragments.empty()) {
605  nucleiIterator ifrag = fragments.begin();
606  for (; ifrag != fragments.end(); ifrag++) {
607  rp = new G4ReactionProduct;
608  dp = makeDynamicParticle(*ifrag);
609  (*rp) = (*dp); // This does all the necessary copying
610  propResult->push_back(rp);
611  delete dp;
612  }
613  }
614 
615  return propResult;
616 }
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4GLOB_DLL std::ostream G4cout
G4DynamicParticle * makeDynamicParticle(const G4InuclElementaryParticle &iep) const
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
#define G4endl
Definition: G4ios.hh:61
G4bool G4CascadeInterface::coulombBarrierViolation ( ) const
protected

Definition at line 657 of file G4CascadeInterface.cc.

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

Referenced by retryInelasticNucleus().

657  {
658  G4bool violated = false; // by default coulomb analysis is OK
659 
660  const G4double coulumbBarrier = 8.7 * MeV/GeV; // Bertini uses GeV
661 
662  const std::vector<G4InuclElementaryParticle>& p =
663  output->getOutgoingParticles();
664 
665  for (particleIterator ipart=p.begin(); ipart != p.end(); ipart++) {
666  if (ipart->type() == proton) {
667  violated |= (ipart->getKineticEnergy() < coulumbBarrier);
668  }
669  }
670 
671  return violated;
672 }
const char * p
Definition: xmltok.h:285
bool G4bool
Definition: G4Types.hh:79
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
double G4double
Definition: G4Types.hh:76
G4bool G4CascadeInterface::createBullet ( const G4HadProjectile aTrack)
protected

Definition at line 439 of file G4CascadeInterface.cc.

References CLHEP::HepLorentzVector::e(), G4InuclElementaryParticle::fill(), G4InuclNuclei::fill(), G4cerr, G4cout, G4endl, G4HadProjectile::Get4Momentum(), G4ParticleDefinition::GetAtomicMass(), G4ParticleDefinition::GetAtomicNumber(), G4HadProjectile::GetDefinition(), G4ParticleDefinition::GetParticleName(), python.hepunit::GeV, CLHEP::HepLorentzRotation::IDENTITY, CLHEP::HepLorentzRotation::invert(), CLHEP::HepLorentzVector::phi(), CLHEP::HepLorentzVector::rho(), CLHEP::HepLorentzRotation::rotateY(), CLHEP::HepLorentzRotation::rotateZ(), CLHEP::HepLorentzVector::theta(), G4InuclElementaryParticle::type(), G4InuclElementaryParticle::valid(), and G4HadronicInteraction::verboseLevel.

Referenced by ApplyYourself(), and Propagate().

439  {
440  const G4ParticleDefinition* trkDef = aTrack.GetDefinition();
441  G4int bulletType = 0; // For elementary particles
442  G4int bulletA = 0, bulletZ = 0; // For nucleus projectile
443 
444  if (trkDef->GetAtomicMass() <= 1) {
445  bulletType = G4InuclElementaryParticle::type(trkDef);
446  } else {
447  bulletA = trkDef->GetAtomicMass();
448  bulletZ = trkDef->GetAtomicNumber();
449  }
450 
451  if (0 == bulletType && 0 == bulletA*bulletZ) {
452  if (verboseLevel) {
453  G4cerr << " G4CascadeInterface: " << trkDef->GetParticleName()
454  << " not usable as bullet." << G4endl;
455  }
456  bullet = 0;
457  return false;
458  }
459 
460  // Code momentum and energy -- Bertini wants z-axis and GeV units
461  G4LorentzVector projectileMomentum = aTrack.Get4Momentum()/GeV;
462 
463  // Rotation/boost to get from z-axis back to original frame
464  bulletInLabFrame = G4LorentzRotation::IDENTITY; // Initialize
465  bulletInLabFrame.rotateZ(-projectileMomentum.phi());
466  bulletInLabFrame.rotateY(-projectileMomentum.theta());
467  bulletInLabFrame.invert();
468 
469  G4LorentzVector momentumBullet(0., 0., projectileMomentum.rho(),
470  projectileMomentum.e());
471 
472  if (G4InuclElementaryParticle::valid(bulletType)) {
473  hadronBullet.fill(momentumBullet, bulletType);
474  bullet = &hadronBullet;
475  } else {
476  nucleusBullet.fill(momentumBullet, bulletA, bulletZ);
477  bullet = &nucleusBullet;
478  }
479 
480  if (verboseLevel > 2) G4cout << "Bullet: \n" << *bullet << G4endl;
481 
482  return true;
483 }
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
HepLorentzRotation & rotateY(double delta)
G4int GetAtomicNumber() const
double phi() const
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
double theta() const
double rho() const
const G4LorentzVector & Get4Momentum() const
G4int GetAtomicMass() const
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation & invert()
void fill(G4int ityp, Model model=DefaultModel)
#define G4endl
Definition: G4ios.hh:61
static DLL_API const HepLorentzRotation IDENTITY
G4GLOB_DLL std::ostream G4cerr
G4bool G4CascadeInterface::createTarget ( G4Nucleus theNucleus)
protected

Definition at line 488 of file G4CascadeInterface.cc.

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

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

488  {
489  return createTarget(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt());
490 }
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4bool createTarget(G4Nucleus &theNucleus)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4bool G4CascadeInterface::createTarget ( G4V3DNucleus theNucleus)
protected

Definition at line 492 of file G4CascadeInterface.cc.

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

492  {
493  return createTarget(theNucleus->GetMassNumber(), theNucleus->GetCharge());
494 }
virtual G4int GetCharge()=0
virtual G4int GetMassNumber()=0
G4bool createTarget(G4Nucleus &theNucleus)
G4bool G4CascadeInterface::createTarget ( G4int  A,
G4int  Z 
)
protected

Definition at line 496 of file G4CascadeInterface.cc.

References G4InuclElementaryParticle::fill(), G4InuclNuclei::fill(), G4cout, G4endl, neutron, G4InuclParticleNames::proton, and G4HadronicInteraction::verboseLevel.

496  {
497  if (A > 1) {
498  nucleusTarget.fill(A, Z);
499  target = &nucleusTarget;
500  } else {
501  hadronTarget.fill(0., (Z==1?proton:neutron));
502  target = &hadronTarget;
503  }
504 
505  if (verboseLevel > 2) G4cout << "Target: \n" << *target << G4endl;
506 
507  return true; // Right now, target never fails
508 }
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
const XML_Char * target
G4GLOB_DLL std::ostream G4cout
void fill(G4int ityp, Model model=DefaultModel)
#define G4endl
Definition: G4ios.hh:61
void G4CascadeInterface::DumpConfiguration ( std::ostream &  outFile) const
virtual

Definition at line 182 of file G4CascadeInterface.cc.

References G4CascadeParameters::DumpConfiguration().

182  {
184 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
static void DumpConfiguration(std::ostream &os)
void G4CascadeInterface::Initialize ( )
static
G4bool G4CascadeInterface::IsApplicable ( const G4HadProjectile aTrack,
G4Nucleus theNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 227 of file G4CascadeInterface.cc.

References G4HadProjectile::GetDefinition().

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

228  {
229  return IsApplicable(aTrack.GetDefinition());
230 }
const G4ParticleDefinition * GetDefinition() const
G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4bool G4CascadeInterface::IsApplicable ( const G4ParticleDefinition aPD) const

Definition at line 232 of file G4CascadeInterface.cc.

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

232  {
233  if (aPD->GetAtomicMass() > 1) return true; // Nuclei are okay
234 
235  // Valid particle and have interactions available
237  return (G4CascadeChannelTables::GetTable(type));
238 }
static const G4CascadeChannel * GetTable(G4int initialState)
int G4int
Definition: G4Types.hh:78
G4int GetAtomicMass() const
G4DynamicParticle * G4CascadeInterface::makeDynamicParticle ( const G4InuclElementaryParticle iep) const
protected

Definition at line 514 of file G4CascadeInterface.cc.

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

Referenced by copyOutputToHadronicResult(), and copyOutputToReactionProducts().

514  {
515  G4int outgoingType = iep.type();
516 
517  if (iep.quasi_deutron()) {
518  G4cerr << " ERROR: G4CascadeInterface incompatible particle type "
519  << outgoingType << G4endl;
520  return 0;
521  }
522 
523  // Copy local G4DynPart to public output (handle kaon mixing specially)
524  if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
525  G4ThreeVector momDir = iep.getMomentum().vect().unit();
526  G4double ekin = iep.getKineticEnergy()*GeV; // Bertini -> G4 units
527 
529  if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
530 
531  return new G4DynamicParticle(pd, momDir, ekin);
532  } else {
533  return new G4DynamicParticle(iep.getDynamicParticle());
534  }
535 
536  return 0; // Should never get here!
537 }
const G4DynamicParticle & getDynamicParticle() const
G4LorentzVector getMomentum() const
int G4int
Definition: G4Types.hh:78
G4double getKineticEnergy() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
static G4KaonZeroLong * Definition()
Hep3Vector unit() const
static G4KaonZeroShort * Definition()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr
G4DynamicParticle * G4CascadeInterface::makeDynamicParticle ( const G4InuclNuclei inuc) const
protected

Definition at line 540 of file G4CascadeInterface.cc.

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

540  {
541  if (verboseLevel > 2) {
542  G4cout << " Nuclei fragment: \n" << inuc << G4endl;
543  }
544 
545  // Copy local G4DynPart to public output
546  return new G4DynamicParticle(inuc.getDynamicParticle());
547 }
const G4DynamicParticle & getDynamicParticle() const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4CascadeInterface::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VIntraNuclearTransportModel.

Definition at line 168 of file G4CascadeInterface.cc.

169 {
170  outFile << "The Bertini-style cascade implements the inelastic scattering\n"
171  << "of hadrons by nuclei. Nucleons, pions, kaons and hyperons\n"
172  << "from 0 to 15 GeV may be used as projectiles in this model.\n"
173  << "Final state hadrons are produced by a classical cascade\n"
174  << "consisting of individual hadron-nucleon scatterings which use\n"
175  << "free-space partial cross sections, corrected for various\n"
176  << "nuclear medium effects. The target nucleus is modeled as a\n"
177  << "set of 1, 3 or 6 spherical shells, in which scattered hadrons\n"
178  << "travel in straight lines until they are reflected from or\n"
179  << "transmitted through shell boundaries.\n";
180 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4HadFinalState * G4CascadeInterface::NoInteraction ( const G4HadProjectile aTrack,
G4Nucleus theNucleus 
)
protected

Definition at line 422 of file G4CascadeInterface.cc.

References G4HadFinalState::Clear(), G4cout, G4endl, G4HadProjectile::GetKineticEnergy(), isAlive, G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetStatusChange(), G4HadronicInteraction::theParticleChange, and G4HadronicInteraction::verboseLevel.

Referenced by ApplyYourself().

423  {
424  if (verboseLevel)
425  G4cout << " >>> G4CascadeInterface::NoInteraction" << G4endl;
426 
429 
430  G4double ekin = aTrack.GetKineticEnergy()>0. ? aTrack.GetKineticEnergy() : 0.;
431  theParticleChange.SetEnergyChange(ekin); // Protect against rounding
432 
433  return &theParticleChange;
434 }
void SetStatusChange(G4HadFinalStateStatus aS)
G4GLOB_DLL std::ostream G4cout
G4double GetKineticEnergy() const
void SetEnergyChange(G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4ReactionProductVector * G4CascadeInterface::Propagate ( G4KineticTrackVector theSecondaries,
G4V3DNucleus theNucleus 
)
virtual

Implements G4VIntraNuclearTransportModel.

Definition at line 347 of file G4CascadeInterface.cc.

References G4HadFinalState::Clear(), clear(), G4CascadeCheckBalance::collide(), 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(), G4CascadeCheckBalance::okay(), G4CollisionOutput::printCollisionOutput(), G4InuclCollider::rescatter(), G4CollisionOutput::reset(), retryInelasticNucleus(), CLHEP::HepRandom::saveEngineStatus(), G4HadronicInteraction::theParticleChange, throwNonConservationFailure(), and G4HadronicInteraction::verboseLevel.

348  {
349  if (verboseLevel) G4cout << " >>> G4CascadeInterface::Propagate" << G4endl;
350 
351 #ifdef G4CASCADE_DEBUG_INTERFACE
352  if (verboseLevel>1) {
353  G4cout << " G4V3DNucleus A " << theNucleus->GetMassNumber()
354  << " Z " << theNucleus->GetCharge()
355  << "\n " << theSecondaries->size() << " secondaries:" << G4endl;
356  for (size_t i=0; i<theSecondaries->size(); i++) {
357  G4KineticTrack* kt = (*theSecondaries)[i];
358  G4cout << " " << i << ": " << kt->GetDefinition()->GetParticleName()
359  << " p " << kt->Get4Momentum() << " @ " << kt->GetPosition()
360  << " t " << kt->GetFormationTime() << G4endl;
361  }
362  }
363 #endif
364 
365  if (!randomFile.empty()) { // User requested random-seed capture
366  if (verboseLevel>1)
367  G4cout << " Saving random engine state to " << randomFile << G4endl;
369  }
370 
372  clear();
373 
374  // Process input secondaries list to eliminate resonances
375  G4DecayKineticTracks decay(theSecondaries);
376 
377  // NOTE: Requires 9.4-ref-03 mods to base class and G4TheoFSGenerator
378  const G4HadProjectile* projectile = GetPrimaryProjectile();
379  if (projectile) createBullet(*projectile);
380 
381  if (!createTarget(theNucleus)) {
382  if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
383  return 0; // FIXME: This will cause a segfault later
384  }
385 
386  numberOfTries = 0;
387  do {
388  if (verboseLevel > 1)
389  G4cout << " Generating rescatter attempt " << numberOfTries << G4endl;
390 
391  output->reset();
392  collider->rescatter(bullet, theSecondaries, theNucleus, *output);
393  balance->collide(bullet, target, *output);
394 
395  numberOfTries++;
396  // FIXME: retry checks will SEGFAULT until we can define the bullet!
397  } while (retryInelasticNucleus());
398 
399  // Check whether repeated attempts have all failed; report and exit
400  if (numberOfTries >= maximumTries && !balance->okay()) {
401  throwNonConservationFailure(); // This terminates the job
402  }
403 
404  // Successful cascade -- clean up and return
405  if (verboseLevel) {
406  G4cout << " Cascade rescatter after trials " << numberOfTries << G4endl;
407  if (verboseLevel > 1) output->printCollisionOutput();
408  }
409 
410  // Does calling code take ownership? I hope so!
412 
413  // Clean up and and return final result
414  clear();
415  return propResult;
416 }
G4bool retryInelasticNucleus() const
void rescatter(G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
virtual G4int GetCharge()=0
const G4HadProjectile * GetPrimaryProjectile() const
const G4ThreeVector & GetPosition() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
void printCollisionOutput(std::ostream &os=G4cout) const
virtual G4int GetMassNumber()=0
const XML_Char * target
const G4String & GetParticleName() const
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4ParticleDefinition * GetDefinition() const
G4GLOB_DLL std::ostream G4cout
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
G4double GetFormationTime() const
G4bool createTarget(G4Nucleus &theNucleus)
G4ReactionProductVector * copyOutputToReactionProducts()
static void saveEngineStatus(const char filename[]="Config.conf")
Definition: Random.cc:175
#define G4endl
Definition: G4ios.hh:61
G4bool createBullet(const G4HadProjectile &aTrack)
const G4LorentzVector & Get4Momentum() const
G4GLOB_DLL std::ostream G4cerr
G4bool G4CascadeInterface::retryInelasticNucleus ( ) const
protected

Definition at line 708 of file G4CascadeInterface.cc.

References coulombBarrierViolation(), G4cout, G4endl, G4InuclParticle::getDefinition(), G4CollisionOutput::getOutgoingParticles(), G4CollisionOutput::numberOfOutgoingNuclei(), G4CollisionOutput::numberOfOutgoingParticles(), and G4CascadeCheckBalance::okay().

Referenced by ApplyYourself(), and Propagate().

708  {
709  // Quantities necessary for conditional block below
710  G4int npart = output->numberOfOutgoingParticles();
711  G4int nfrag = output->numberOfOutgoingNuclei();
712 
713  const G4ParticleDefinition* firstOut = (npart == 0) ? 0 :
714  output->getOutgoingParticles().begin()->getDefinition();
715 
716 #ifdef G4CASCADE_DEBUG_INTERFACE
717  // Report on all retry conditions, in order of return logic
718  G4cout << " retryInelasticNucleus: numberOfTries "
719  << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
720  << "\n retryInelasticNucleus: AND outputParticles "
721  << ((npart != 0) ? "NON-ZERO (t)" : "EMPTY (f)")
722 #ifdef G4CASCADE_COULOMB_DEV
723  << "\n retryInelasticNucleus: AND coulombBarrier (COULOMB_DEV) "
724  << (coulombBarrierViolation() ? "VIOLATED (t)" : "PASSED (f)")
725  << "\n retryInelasticNucleus: AND collision type (COULOMB_DEV) "
726  << ((npart+nfrag > 2) ? "INELASTIC (t)" : "ELASTIC (f)")
727 #else
728  << "\n retryInelasticNucleus: AND collsion type "
729  << ((npart+nfrag < 3) ? "ELASTIC (t)" : "INELASTIC (f)")
730  << "\n retryInelasticNucleus: AND Leading particle bullet "
731  << ((firstOut == bullet->getDefinition()) ? "YES (t)" : "NO (f)")
732 #endif
733  << "\n retryInelasticNucleus: OR conservation "
734  << (!balance->okay() ? "FAILED (t)" : "PASSED (f)")
735  << G4endl;
736 #endif
737 
738  return ( ((numberOfTries < maximumTries) &&
739  (npart != 0) &&
740 #ifdef G4CASCADE_COULOMB_DEV
741  (coulombBarrierViolation() && npart+nfrag > 2)
742 #else
743  (npart+nfrag < 3 && firstOut == bullet->getDefinition())
744 #endif
745  )
746 #ifndef G4CASCADE_SKIP_ECONS
747  || (!balance->okay())
748 #endif
749  );
750 }
G4bool coulombBarrierViolation() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
#define G4endl
Definition: G4ios.hh:61
G4ParticleDefinition * getDefinition() const
G4bool G4CascadeInterface::retryInelasticProton ( ) const
protected

Definition at line 676 of file G4CascadeInterface.cc.

References G4cout, G4endl, G4InuclParticle::getDefinition(), and G4CollisionOutput::getOutgoingParticles().

Referenced by ApplyYourself().

676  {
677  const std::vector<G4InuclElementaryParticle>& out =
678  output->getOutgoingParticles();
679 
680 #ifdef G4CASCADE_DEBUG_INTERFACE
681  // Report on all retry conditions, in order of return logic
682  G4cout << " retryInelasticProton: number of Tries "
683  << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
684  << "\n retryInelasticProton: AND collision type ";
685  if (out.empty()) G4cout << "FAILED" << G4endl;
686  else {
687  G4cout << (out.size() == 2 ? "ELASTIC (t)" : "INELASTIC (f)")
688  << "\n retryInelasticProton: AND Leading particles bullet "
689  << (out.size() >= 2 &&
690  (out[0].getDefinition() == bullet->getDefinition() ||
691  out[1].getDefinition() == bullet->getDefinition())
692  ? "YES (t)" : "NO (f)")
693  << G4endl;
694  }
695 #endif
696 
697  return ( (numberOfTries < maximumTries) &&
698  (out.empty() ||
699  (out.size() == 2 &&
700  (out[0].getDefinition() == bullet->getDefinition() ||
701  out[1].getDefinition() == bullet->getDefinition())))
702  );
703 }
G4GLOB_DLL std::ostream G4cout
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
#define G4endl
Definition: G4ios.hh:61
G4ParticleDefinition * getDefinition() const
void G4CascadeInterface::SetVerboseLevel ( G4int  verbose)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 217 of file G4CascadeInterface.cc.

References G4VCascadeCollider::setVerboseLevel(), G4CollisionOutput::setVerboseLevel(), G4InuclCollider::setVerboseLevel(), and G4HadronicInteraction::SetVerboseLevel().

Referenced by G4CascadeInterface().

217  {
219  collider->setVerboseLevel(verbose);
220  balance->setVerboseLevel(verbose);
221  output->setVerboseLevel(verbose);
222 }
void setVerboseLevel(G4int verbose)
virtual void setVerboseLevel(G4int verbose=0)
void setVerboseLevel(G4int verbose=0)
virtual void SetVerboseLevel(G4int value)
void G4CascadeInterface::throwNonConservationFailure ( )
protected

Definition at line 756 of file G4CascadeInterface.cc.

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

Referenced by ApplyYourself(), and Propagate().

756  {
757  // NOTE: Once G4HadronicException is changed, use the following line!
758  // G4ExceptionDescription errInfo;
759  std::ostream& errInfo = G4cerr;
760 
761  errInfo << " >>> G4CascadeInterface has non-conserving"
762  << " cascade after " << numberOfTries << " attempts." << G4endl;
763 
764  G4String throwMsg = "G4CascadeInterface - ";
765  if (!balance->energyOkay()) {
766  throwMsg += "Energy";
767  errInfo << " Energy conservation violated by " << balance->deltaE()
768  << " GeV (" << balance->relativeE() << ")" << G4endl;
769  }
770 
771  if (!balance->momentumOkay()) {
772  throwMsg += "Momentum";
773  errInfo << " Momentum conservation violated by " << balance->deltaP()
774  << " GeV/c (" << balance->relativeP() << ")" << G4endl;
775  }
776 
777  if (!balance->baryonOkay()) {
778  throwMsg += "Baryon number";
779  errInfo << " Baryon number violated by " << balance->deltaB() << G4endl;
780  }
781 
782  if (!balance->chargeOkay()) {
783  throwMsg += "Charge";
784  errInfo << " Charge conservation violated by " << balance->deltaQ()
785  << G4endl;
786  }
787 
788  errInfo << " Final event output, for debugging:\n"
789  << " Bullet: \n" << *bullet << G4endl
790  << " Target: \n" << *target << G4endl;
791  output->printCollisionOutput(errInfo);
792 
793  throwMsg += " non-conservation. More info in output.";
794  throw G4HadronicException(__FILE__, __LINE__, throwMsg); // Job ends here!
795 }
void printCollisionOutput(std::ostream &os=G4cout) const
const XML_Char * target
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
void G4CascadeInterface::useCascadeDeexcitation ( )

Definition at line 206 of file G4CascadeInterface.cc.

References G4InuclCollider::useCascadeDeexcitation().

206  {
207  collider->useCascadeDeexcitation();
208 }
void useCascadeDeexcitation()
void G4CascadeInterface::usePreCompoundDeexcitation ( )

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