#include <G4INCLNucleus.hh>
Inheritance diagram for G4INCL::Nucleus:
Public Member Functions | |
Nucleus (G4int mass, G4int charge, Config const *const conf, const G4double universeRadius=-1.) | |
virtual | ~Nucleus () |
void | initializeParticles () |
void | insertParticle (Particle *p) |
Insert a new particle (e.g. a projectile) in the nucleus. | |
void | applyFinalState (FinalState *) |
G4int | getInitialA () const |
G4int | getInitialZ () const |
ParticleList const & | getCreatedParticles () const |
ParticleList const & | getUpdatedParticles () const |
Particle * | getBlockedDelta () const |
Get the delta that could not decay. | |
void | propagateParticles (G4double step) |
G4int | getNumberOfEnteringProtons () const |
G4int | getNumberOfEnteringNeutrons () const |
G4double | computeSeparationEnergyBalance () const |
Outgoing - incoming separation energies. | |
G4bool | decayOutgoingDeltas () |
Force the decay of outgoing deltas. | |
G4bool | decayInsideDeltas () |
Force the decay of deltas inside the nucleus. | |
G4bool | decayOutgoingClusters () |
Force the decay of unstable outgoing clusters. | |
G4bool | decayMe () |
Force the phase-space decay of the Nucleus. | |
void | emitInsidePions () |
Force emission of all pions inside the nucleus. | |
void | computeRecoilKinematics () |
Compute the recoil momentum and spin of the nucleus. | |
ThreeVector | computeCenterOfMass () const |
Compute the current center-of-mass position. | |
G4double | computeTotalEnergy () const |
Compute the current total energy. | |
G4double | computeExcitationEnergy () const |
Compute the current excitation energy. | |
void | setIncomingAngularMomentum (const ThreeVector &j) |
Set the incoming angular-momentum vector. | |
const ThreeVector & | getIncomingAngularMomentum () const |
Get the incoming angular-momentum vector. | |
void | setIncomingMomentum (const ThreeVector &p) |
Set the incoming momentum vector. | |
const ThreeVector & | getIncomingMomentum () const |
Get the incoming momentum vector. | |
void | setInitialEnergy (const G4double e) |
Set the initial energy. | |
G4double | getInitialEnergy () const |
Get the initial energy. | |
G4double | getExcitationEnergy () const |
Get the excitation energy of the nucleus. | |
G4bool | containsDeltas () |
Returns true if the nucleus contains any deltas. | |
std::string | print () |
std::string | dump () |
Store * | getStore () const |
void | setStore (Store *s) |
G4double | getInitialInternalEnergy () const |
G4bool | isEventTransparent () const |
Is the event transparent? | |
G4bool | hasRemnant () const |
Does the nucleus give a cascade remnant? | |
void | fillEventInfo (EventInfo *eventInfo) |
G4bool | getTryCompoundNucleus () |
G4int | getProjectileChargeNumber () const |
Return the charge number of the projectile. | |
G4int | getProjectileMassNumber () const |
Return the mass number of the projectile. | |
void | setProjectileChargeNumber (G4int n) |
Set the charge number of the projectile. | |
void | setProjectileMassNumber (G4int n) |
Set the mass number of the projectile. | |
G4bool | isForcedTransparent () |
Returns true if a transparent event should be forced. | |
G4double | getTransmissionBarrier (Particle const *const p) |
Get the transmission barrier. | |
ConservationBalance | getConservationBalance (EventInfo const &theEventInfo, const G4bool afterRecoil) const |
Compute charge, mass, energy and momentum balance. | |
void | useFusionKinematics () |
Adjust the kinematics for complete-fusion events. | |
G4double | getSurfaceRadius (Particle const *const particle) const |
Get the maximum allowed radius for a given particle. | |
G4double | getCoulombRadius (ParticleSpecies const &p) const |
Get the Coulomb radius for a given particle. | |
G4double | getUniverseRadius () const |
Getter for theUniverseRadius. | |
void | setUniverseRadius (const G4double universeRadius) |
Setter for theUniverseRadius. | |
G4bool | isNucleusNucleusCollision () const |
Is it a nucleus-nucleus collision? | |
void | setNucleusNucleusCollision () |
Set a nucleus-nucleus collision. | |
void | setParticleNucleusCollision () |
Set a particle-nucleus collision. | |
void | setProjectileRemnant (ProjectileRemnant *const c) |
Set the projectile remnant. | |
ProjectileRemnant * | getProjectileRemnant () const |
Get the projectile remnant. | |
void | deleteProjectileRemnant () |
Delete the projectile remnant. | |
void | moveProjectileRemnantComponentsToOutgoing () |
Move the components of the projectile remnant to the outgoing list. | |
void | finalizeProjectileRemnant (const G4double emissionTime) |
Finalise the projectile remnant. | |
void | updatePotentialEnergy (Particle *p) |
Update the particle potential energy. | |
NuclearDensity * | getDensity () const |
Getter for theDensity. | |
NuclearPotential::INuclearPotential * | getPotential () const |
Getter for thePotential. | |
Data Structures | |
struct | ConservationBalance |
Struct for conservation laws. More... |
Definition at line 64 of file G4INCLNucleus.hh.
G4INCL::Nucleus::Nucleus | ( | G4int | mass, | |
G4int | charge, | |||
Config const *const | conf, | |||
const G4double | universeRadius = -1. | |||
) |
Definition at line 68 of file G4INCLNucleus.cc.
References G4INCL::ConstantPotential, G4INCL::NuclearDensityFactory::createDensity(), FATAL, G4INCL::Config::getPionPotential(), G4INCL::Config::getPotentialType(), G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy(), G4INCL::IsospinEnergyPotential, G4INCL::IsospinEnergySmoothPotential, G4INCL::IsospinPotential, G4INCL::Neutron, G4INCL::Proton, G4INCL::ParticleSampler::setDensity(), G4INCL::ParticleTable::setNeutronSeparationEnergy(), G4INCL::ParticleSampler::setPotential(), G4INCL::ParticleTable::setProtonSeparationEnergy(), G4INCL::Particle::theA, G4INCL::Cluster::theParticleSampler, and G4INCL::Particle::theZ.
00069 : Cluster(charge,mass), 00070 theInitialZ(charge), theInitialA(mass), 00071 theNpInitial(0), theNnInitial(0), 00072 initialInternalEnergy(0.), 00073 incomingAngularMomentum(0.,0.,0.), incomingMomentum(0.,0.,0.), 00074 initialCenterOfMass(0.,0.,0.), 00075 remnant(true), 00076 blockedDelta(NULL), 00077 initialEnergy(0.), 00078 tryCN(false), 00079 forceTransparent(false), 00080 projectileZ(0), 00081 projectileA(0), 00082 theUniverseRadius(universeRadius), 00083 isNucleusNucleus(false), 00084 theProjectileRemnant(NULL), 00085 theDensity(NULL), 00086 thePotential(NULL) 00087 { 00088 PotentialType potentialType; 00089 G4bool pionPotential; 00090 if(conf) { 00091 potentialType = conf->getPotentialType(); 00092 pionPotential = conf->getPionPotential(); 00093 } else { // By default we don't use energy dependent 00094 // potential. This is convenient for some tests. 00095 potentialType = IsospinPotential; 00096 pionPotential = true; 00097 } 00098 switch(potentialType) { 00099 case IsospinEnergySmoothPotential: 00100 thePotential = new NuclearPotential::NuclearPotentialEnergyIsospinSmooth(theA, theZ, pionPotential); 00101 break; 00102 case IsospinEnergyPotential: 00103 thePotential = new NuclearPotential::NuclearPotentialEnergyIsospin(theA, theZ, pionPotential); 00104 break; 00105 case IsospinPotential: 00106 thePotential = new NuclearPotential::NuclearPotentialIsospin(theA, theZ, pionPotential); 00107 break; 00108 case ConstantPotential: 00109 thePotential = new NuclearPotential::NuclearPotentialConstant(theA, theZ, pionPotential); 00110 break; 00111 default: 00112 FATAL("Unrecognized potential type at Nucleus creation." << std::endl); 00113 break; 00114 } 00115 00116 ParticleTable::setProtonSeparationEnergy(thePotential->getSeparationEnergy(Proton)); 00117 ParticleTable::setNeutronSeparationEnergy(thePotential->getSeparationEnergy(Neutron)); 00118 00119 theDensity = NuclearDensityFactory::createDensity(theA, theZ); 00120 00121 theParticleSampler->setPotential(thePotential); 00122 theParticleSampler->setDensity(theDensity); 00123 00124 if(theUniverseRadius<0) 00125 theUniverseRadius = theDensity->getMaximumRadius(); 00126 theStore = new Store(conf); 00127 toBeUpdated.clear(); 00128 }
G4INCL::Nucleus::~Nucleus | ( | ) | [virtual] |
Definition at line 130 of file G4INCLNucleus.cc.
00130 { 00131 delete theStore; 00132 delete thePotential; 00133 /* We don't delete the density here any more -- the Factory is caching them 00134 delete theDensity;*/ 00135 }
void G4INCL::Nucleus::applyFinalState | ( | FinalState * | ) |
Apply reaction final state information to the nucleus.
Definition at line 166 of file G4INCLNucleus.cc.
References G4INCL::Store::add(), G4INCL::Store::addToOutgoing(), DEBUG, deleted, ERROR, G4INCL::FinalState::getBlockedDelta(), G4INCL::Store::getBook(), G4INCL::FinalState::getCreatedParticles(), G4INCL::Book::getCurrentTime(), G4INCL::FinalState::getDestroyedParticles(), G4INCL::FinalState::getEnteringParticles(), G4INCL::FinalState::getModifiedParticles(), G4INCL::FinalState::getOutgoingParticles(), G4INCL::FinalState::getTotalEnergyBeforeInteraction(), G4INCL::FinalState::getValidity(), insertParticle(), G4INCL::ParticleBelowFermiFS, G4INCL::ParticleBelowZeroFS, G4INCL::Store::particleHasBeenDestroyed(), G4INCL::Store::particleHasBeenEjected(), G4INCL::Store::particleHasBeenUpdated(), G4INCL::PauliBlockedFS, G4INCL::FinalState::print(), G4INCL::Particle::theA, G4INCL::Particle::theZ, and G4INCL::ValidFS.
Referenced by decayInsideDeltas().
00166 { 00167 justCreated.clear(); 00168 toBeUpdated.clear(); // Clear the list of particles to be updated by the propagation model. 00169 blockedDelta = NULL; 00170 G4double totalEnergy = 0.0; 00171 00172 FinalStateValidity const validity = finalstate->getValidity(); 00173 if(validity == ValidFS) { 00174 00175 ParticleList const &created = finalstate->getCreatedParticles(); 00176 for(ParticleIter iter = created.begin(); iter != created.end(); ++iter) { 00177 theStore->add((*iter)); 00178 if(!(*iter)->isOutOfWell()) { 00179 totalEnergy += (*iter)->getEnergy() - (*iter)->getPotentialEnergy(); 00180 justCreated.push_back((*iter)); // New particle, so we must create avatars for it 00181 } 00182 } 00183 00184 ParticleList const &deleted = finalstate->getDestroyedParticles(); 00185 for(ParticleIter iter = deleted.begin(); iter != deleted.end(); ++iter) { 00186 theStore->particleHasBeenDestroyed((*iter)->getID()); 00187 } 00188 00189 ParticleList const &modified = finalstate->getModifiedParticles(); 00190 for(ParticleIter iter = modified.begin(); iter != modified.end(); ++iter) { 00191 theStore->particleHasBeenUpdated((*iter)->getID()); 00192 totalEnergy += (*iter)->getEnergy() - (*iter)->getPotentialEnergy(); 00193 toBeUpdated.push_back((*iter)); // Particle is modified so we have to create new avatars for it. 00194 } 00195 00196 ParticleList const &out = finalstate->getOutgoingParticles(); 00197 for(ParticleIter iter = out.begin(); iter != out.end(); ++iter) { 00198 if((*iter)->isCluster()) { 00199 Cluster *clusterOut = dynamic_cast<Cluster*>((*iter)); 00200 ParticleList const components = clusterOut->getParticles(); 00201 for(ParticleIter in = components.begin(); in != components.end(); ++in) 00202 theStore->particleHasBeenEjected((*in)->getID()); 00203 } else { 00204 theStore->particleHasBeenEjected((*iter)->getID()); 00205 } 00206 totalEnergy += (*iter)->getEnergy(); // No potential here because the particle is gone 00207 theA -= (*iter)->getA(); 00208 theZ -= (*iter)->getZ(); 00209 theStore->addToOutgoing(*iter); 00210 (*iter)->setEmissionTime(theStore->getBook()->getCurrentTime()); 00211 } 00212 00213 ParticleList const &entering = finalstate->getEnteringParticles(); 00214 for(ParticleIter iter = entering.begin(); iter != entering.end(); ++iter) { 00215 insertParticle(*iter); 00216 totalEnergy += (*iter)->getEnergy() - (*iter)->getPotentialEnergy(); 00217 toBeUpdated.push_back((*iter)); // Particle is modified so we have to create new avatars for it. 00218 } 00219 } else if(validity == PauliBlockedFS) { 00220 blockedDelta = finalstate->getBlockedDelta(); 00221 } else if(validity == ParticleBelowFermiFS) { 00222 DEBUG("A Particle is entering below the Fermi sea:" << std::endl << finalstate->print() << std::endl); 00223 tryCN = true; 00224 ParticleList const &entering = finalstate->getEnteringParticles(); 00225 for(ParticleIter iter = entering.begin(); iter != entering.end(); ++iter) { 00226 insertParticle(*iter); 00227 } 00228 } else if(validity == ParticleBelowZeroFS) { 00229 DEBUG("A Particle is entering below zero energy:" << std::endl << finalstate->print() << std::endl); 00230 forceTransparent = true; 00231 ParticleList const &entering = finalstate->getEnteringParticles(); 00232 for(ParticleIter iter = entering.begin(); iter != entering.end(); ++iter) { 00233 insertParticle(*iter); 00234 } 00235 } 00236 00237 if(validity==ValidFS && 00238 std::abs(totalEnergy - finalstate->getTotalEnergyBeforeInteraction()) > 0.1) { 00239 ERROR("Energy nonconservation! Energy at the beginning of the event = " 00240 << finalstate->getTotalEnergyBeforeInteraction() 00241 <<" and after interaction = " 00242 << totalEnergy << std::endl 00243 << finalstate->print()); 00244 } 00245 }
ThreeVector G4INCL::Nucleus::computeCenterOfMass | ( | ) | const |
Compute the current center-of-mass position.
Definition at line 295 of file G4INCLNucleus.cc.
References G4INCL::Store::getParticles().
Referenced by computeRecoilKinematics().
00295 { 00296 ThreeVector cm(0.,0.,0.); 00297 G4double totalMass = 0.0; 00298 ParticleList inside = theStore->getParticles(); 00299 for(ParticleIter p=inside.begin(); p!=inside.end(); ++p) { 00300 const G4double mass = (*p)->getMass(); 00301 cm += (*p)->getPosition() * mass; 00302 totalMass += mass; 00303 } 00304 cm /= totalMass; 00305 return cm; 00306 }
G4double G4INCL::Nucleus::computeExcitationEnergy | ( | ) | const |
Compute the current excitation energy.
Definition at line 308 of file G4INCLNucleus.cc.
References computeSeparationEnergyBalance(), and computeTotalEnergy().
00308 { 00309 const G4double totalEnergy = computeTotalEnergy(); 00310 const G4double separationEnergies = computeSeparationEnergyBalance(); 00311 00312 return totalEnergy - initialInternalEnergy - separationEnergies; 00313 }
void G4INCL::Nucleus::computeRecoilKinematics | ( | ) |
Compute the recoil momentum and spin of the nucleus.
Definition at line 265 of file G4INCLNucleus.cc.
References G4INCL::Particle::adjustEnergyFromMomentum(), computeCenterOfMass(), emitInsidePions(), G4INCL::Store::getOutgoingParticles(), G4INCL::ParticleTable::getTableMass, G4INCL::Particle::setMass(), G4INCL::Particle::theA, G4INCL::Cluster::theExcitationEnergy, G4INCL::Particle::theMomentum, G4INCL::Particle::thePosition, G4INCL::Cluster::theSpin, and G4INCL::Particle::theZ.
00265 { 00266 // If the remnant consists of only one nucleon, we need to apply a special 00267 // procedure to put it on mass shell. 00268 if(theA==1) { 00269 emitInsidePions(); 00270 computeOneNucleonRecoilKinematics(); 00271 remnant=false; 00272 return; 00273 } 00274 00275 // Compute the recoil momentum and angular momentum 00276 theMomentum = incomingMomentum; 00277 theSpin = incomingAngularMomentum; 00278 00279 ParticleList outgoing = theStore->getOutgoingParticles(); 00280 for(ParticleIter p=outgoing.begin(); p!=outgoing.end(); ++p) 00281 { 00282 theMomentum -= (*p)->getMomentum(); 00283 theSpin -= (*p)->getAngularMomentum(); 00284 } 00285 00286 // Subtract orbital angular momentum 00287 thePosition = computeCenterOfMass(); 00288 theSpin -= (thePosition-initialCenterOfMass).vector(theMomentum); 00289 00290 setMass(ParticleTable::getTableMass(theA,theZ) + theExcitationEnergy); 00291 adjustEnergyFromMomentum(); 00292 remnant=true; 00293 }
G4double G4INCL::Nucleus::computeSeparationEnergyBalance | ( | ) | const [inline] |
Outgoing - incoming separation energies.
Used by CDPP.
Definition at line 122 of file G4INCLNucleus.hh.
References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::Store::getOutgoingParticles(), G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy(), G4INCL::Neutron, G4INCL::PiMinus, G4INCL::PiPlus, and G4INCL::Proton.
Referenced by computeExcitationEnergy(), and G4INCL::CDPP::isBlocked().
00122 { 00123 G4double S = 0.0; 00124 ParticleList outgoing = theStore->getOutgoingParticles(); 00125 for(ParticleIter i = outgoing.begin(); i != outgoing.end(); ++i) { 00126 const ParticleType t = (*i)->getType(); 00127 switch(t) { 00128 case Proton: 00129 case Neutron: 00130 case DeltaPlusPlus: 00131 case DeltaPlus: 00132 case DeltaZero: 00133 case DeltaMinus: 00134 S += thePotential->getSeparationEnergy(*i); 00135 break; 00136 case Composite: 00137 S += (*i)->getZ() * thePotential->getSeparationEnergy(Proton) 00138 + ((*i)->getA() - (*i)->getZ()) * thePotential->getSeparationEnergy(Neutron); 00139 break; 00140 case PiPlus: 00141 S += thePotential->getSeparationEnergy(Proton) - thePotential->getSeparationEnergy(Neutron); 00142 break; 00143 case PiMinus: 00144 S += thePotential->getSeparationEnergy(Neutron) - thePotential->getSeparationEnergy(Proton); 00145 break; 00146 default: 00147 break; 00148 } 00149 } 00150 00151 S -= theNpInitial * thePotential->getSeparationEnergy(Proton); 00152 S -= theNnInitial * thePotential->getSeparationEnergy(Neutron); 00153 return S; 00154 }
G4double G4INCL::Nucleus::computeTotalEnergy | ( | ) | const |
Compute the current total energy.
Definition at line 251 of file G4INCLNucleus.cc.
References G4INCL::ParticleTable::effectiveNucleonMass, and G4INCL::Store::getParticles().
Referenced by computeExcitationEnergy(), and initializeParticles().
00251 { 00252 G4double totalEnergy = 0.0; 00253 ParticleList inside = theStore->getParticles(); 00254 for(ParticleIter p=inside.begin(); p!=inside.end(); ++p) { 00255 if((*p)->isNucleon()) // Ugly: we should calculate everything using total energies! 00256 totalEnergy += (*p)->getKineticEnergy() - (*p)->getPotentialEnergy(); 00257 else if((*p)->isResonance()) 00258 totalEnergy += (*p)->getEnergy() - (*p)->getPotentialEnergy() - ParticleTable::effectiveNucleonMass; 00259 else 00260 totalEnergy += (*p)->getEnergy() - (*p)->getPotentialEnergy(); 00261 } 00262 return totalEnergy; 00263 }
G4bool G4INCL::Nucleus::containsDeltas | ( | ) | [inline] |
Returns true if the nucleus contains any deltas.
Definition at line 237 of file G4INCLNucleus.hh.
References G4INCL::Store::getParticles().
00237 { 00238 ParticleList inside = theStore->getParticles(); 00239 for(ParticleIter i=inside.begin(); i!=inside.end(); ++i) 00240 if((*i)->isDelta()) return true; 00241 return false; 00242 }
G4bool G4INCL::Nucleus::decayInsideDeltas | ( | ) |
Force the decay of deltas inside the nucleus.
Definition at line 392 of file G4INCLNucleus.cc.
References applyFinalState(), DEBUG, emitInsidePions(), G4INCL::IAvatar::getFinalState(), G4INCL::Store::getParticles(), G4INCL::FinalState::getValidity(), G4INCL::NuclearPotential::INuclearPotential::hasPionPotential(), G4INCL::Particle::theA, G4INCL::Particle::theZ, and G4INCL::ValidFS.
00392 { 00393 /* If there is a pion potential, do nothing (deltas will be counted as 00394 * excitation energy). 00395 * If, however, the remnant is unphysical (Z<0 or Z>A), force the deltas to 00396 * decay and get rid of all the pions. In case you're wondering, you can 00397 * end up with Z<0 or Z>A if the remnant contains more pi- than protons or 00398 * more pi+ than neutrons, respectively. 00399 */ 00400 const G4bool unphysicalRemnant = (theZ<0 || theZ>theA); 00401 if(thePotential->hasPionPotential() && !unphysicalRemnant) 00402 return false; 00403 00404 // Build a list of deltas (avoid modifying the list you are iterating on). 00405 ParticleList inside = theStore->getParticles(); 00406 ParticleList deltas; 00407 for(ParticleIter i = inside.begin(); i != inside.end(); ++i) 00408 if((*i)->isDelta()) deltas.push_back((*i)); 00409 00410 // Loop over the deltas, make them decay 00411 for(ParticleIter i = deltas.begin(); i != deltas.end(); ++i) { 00412 DEBUG("Decay inside delta particle:" << std::endl 00413 << (*i)->print() << std::endl); 00414 // Create a forced-decay avatar. Note the last boolean parameter. Note 00415 // also that if the remnant is unphysical we more or less explicitly give 00416 // up energy conservation and CDPP by passing a NULL pointer for the 00417 // nucleus. 00418 IAvatar *decay; 00419 if(unphysicalRemnant) 00420 decay = new DecayAvatar((*i), 0.0, NULL, true); 00421 else 00422 decay = new DecayAvatar((*i), 0.0, this, true); 00423 FinalState *fs = decay->getFinalState(); 00424 00425 // The pion can be ejected only if we managed to satisfy energy 00426 // conservation and if pion emission does not lead to negative excitation 00427 // energies. 00428 if(fs->getValidity()==ValidFS) { 00429 // Apply the final state to the nucleus 00430 applyFinalState(fs); 00431 } 00432 delete fs; 00433 delete decay; 00434 } 00435 00436 // If the remnant is unphysical, emit all the pions 00437 if(unphysicalRemnant) { 00438 DEBUG("Remnant is unphysical: Z=" << theZ << ", A=" << theA << std::endl); 00439 emitInsidePions(); 00440 } 00441 00442 return true; 00443 }
G4bool G4INCL::Nucleus::decayMe | ( | ) |
Force the phase-space decay of the Nucleus.
Only applied if Z==0 or Z==A.
Definition at line 463 of file G4INCLNucleus.cc.
References G4INCL::Store::addToOutgoing(), G4INCL::ClusterDecay::decay(), G4INCL::Particle::theA, and G4INCL::Particle::theZ.
00463 { 00464 // Do the phase-space decay only if Z=0 or Z=A 00465 if(theA<=1 || (theZ!=0 && theA!=theZ)) 00466 return false; 00467 00468 ParticleList decayProducts = ClusterDecay::decay(this); 00469 for(ParticleIter j = decayProducts.begin(); j!=decayProducts.end(); ++j) 00470 theStore->addToOutgoing(*j); 00471 00472 return true; 00473 }
G4bool G4INCL::Nucleus::decayOutgoingClusters | ( | ) |
Force the decay of unstable outgoing clusters.
Definition at line 445 of file G4INCLNucleus.cc.
References G4INCL::Store::addToOutgoing(), G4INCL::ClusterDecay::decay(), and G4INCL::Store::getOutgoingParticles().
00445 { 00446 ParticleList out = theStore->getOutgoingParticles(); 00447 ParticleList clusters; 00448 for(ParticleIter i = out.begin(); i != out.end(); ++i) { 00449 if((*i)->isCluster()) clusters.push_back((*i)); 00450 } 00451 if(clusters.empty()) return false; 00452 00453 for(ParticleIter i = clusters.begin(); i != clusters.end(); ++i) { 00454 Cluster *cluster = dynamic_cast<Cluster*>(*i); // Can't avoid using a cast here 00455 cluster->deleteParticles(); // Don't need them 00456 ParticleList decayProducts = ClusterDecay::decay(cluster); 00457 for(ParticleIter j = decayProducts.begin(); j!=decayProducts.end(); ++j) 00458 theStore->addToOutgoing(*j); 00459 } 00460 return true; 00461 }
G4bool G4INCL::Nucleus::decayOutgoingDeltas | ( | ) |
Force the decay of outgoing deltas.
Definition at line 339 of file G4INCLNucleus.cc.
References G4INCL::Store::addToOutgoing(), G4INCL::Particle::adjustEnergyFromMomentum(), G4INCL::Particle::boost(), DEBUG, G4INCL::Store::getBook(), G4INCL::FinalState::getCreatedParticles(), G4INCL::Book::getCurrentTime(), G4INCL::IAvatar::getFinalState(), G4INCL::FinalState::getModifiedParticles(), G4INCL::Particle::getMomentum(), G4INCL::Store::getOutgoingParticles(), G4INCL::Particle::getTableMass(), G4INCL::ThreeVector::mag(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::Particle::setEmissionTime(), G4INCL::Particle::setMomentum(), and G4INCL::Particle::setTableMass().
00339 { 00340 ParticleList out = theStore->getOutgoingParticles(); 00341 ParticleList deltas; 00342 for(ParticleIter i = out.begin(); i != out.end(); ++i) { 00343 if((*i)->isDelta()) deltas.push_back((*i)); 00344 } 00345 if(deltas.empty()) return false; 00346 00347 for(ParticleIter i = deltas.begin(); i != deltas.end(); ++i) { 00348 DEBUG("Decay outgoing delta particle:" << std::endl 00349 << (*i)->print() << std::endl); 00350 const ThreeVector beta = -(*i)->boostVector(); 00351 const G4double deltaMass = (*i)->getMass(); 00352 00353 // Set the delta momentum to zero and sample the decay in the CM frame. 00354 // This makes life simpler if we are using real particle masses. 00355 (*i)->setMomentum(ThreeVector()); 00356 (*i)->setEnergy((*i)->getMass()); 00357 00358 // Use a DecayAvatar 00359 IAvatar *decay = new DecayAvatar((*i), 0.0, NULL); 00360 FinalState *fs = decay->getFinalState(); 00361 Particle * const pion = fs->getCreatedParticles().front(); 00362 Particle * const nucleon = fs->getModifiedParticles().front(); 00363 00364 // Adjust the decay momentum if we are using the real masses 00365 const G4double decayMomentum = KinematicsUtils::momentumInCM(deltaMass, 00366 nucleon->getTableMass(), 00367 pion->getTableMass()); 00368 ThreeVector newMomentum = pion->getMomentum(); 00369 newMomentum *= decayMomentum / newMomentum.mag(); 00370 00371 pion->setTableMass(); 00372 pion->setMomentum(newMomentum); 00373 pion->adjustEnergyFromMomentum(); 00374 pion->setEmissionTime(theStore->getBook()->getCurrentTime()); 00375 pion->boost(beta); 00376 00377 nucleon->setTableMass(); 00378 nucleon->setMomentum(-newMomentum); 00379 nucleon->adjustEnergyFromMomentum(); 00380 nucleon->setEmissionTime(theStore->getBook()->getCurrentTime()); 00381 nucleon->boost(beta); 00382 00383 theStore->addToOutgoing(pion); 00384 00385 delete fs; 00386 delete decay; 00387 } 00388 00389 return true; 00390 }
void G4INCL::Nucleus::deleteProjectileRemnant | ( | ) | [inline] |
Delete the projectile remnant.
Definition at line 401 of file G4INCLNucleus.hh.
Referenced by finalizeProjectileRemnant().
std::string G4INCL::Nucleus::dump | ( | ) |
Definition at line 152 of file G4INCLNucleus.cc.
References G4INCL::Store::getParticipants().
00152 { 00153 std::stringstream ss; 00154 ss <<"(list ;; List of participants " << std::endl; 00155 ParticleList participants = theStore->getParticipants(); 00156 for(ParticleIter i = participants.begin(); i != participants.end(); ++i) { 00157 ss <<"(make-particle-avatar-map " << std::endl 00158 << (*i)->dump() 00159 << "(list ;; List of avatars in this particle" << std::endl 00160 << ")) ;; Close the list of avatars and the particle-avatar-map" << std::endl; 00161 } 00162 ss << ")" << std::endl; 00163 return ss.str(); 00164 }
void G4INCL::Nucleus::emitInsidePions | ( | ) |
Force emission of all pions inside the nucleus.
Definition at line 475 of file G4INCLNucleus.cc.
References G4INCL::Store::addToOutgoing(), G4INCL::Store::getBook(), G4INCL::Book::getCurrentTime(), G4INCL::Store::getParticles(), G4INCL::Store::particleHasBeenEjected(), G4INCL::Particle::theA, G4INCL::Particle::theZ, and WARN.
Referenced by computeRecoilKinematics(), and decayInsideDeltas().
00475 { 00476 /* Forcing emissions of all pions in the nucleus. This probably violates 00477 * energy conservation (although the computation of the recoil kinematics 00478 * might sweep this under the carpet). 00479 */ 00480 WARN("Forcing emissions of all pions in the nucleus." << std::endl); 00481 00482 // Emit the pions with this kinetic energy 00483 const G4double tinyPionEnergy = 0.1; // MeV 00484 00485 // Push out the emitted pions 00486 ParticleList inside = theStore->getParticles(); 00487 for(ParticleIter i = inside.begin(); i != inside.end(); ++i) { 00488 if((*i)->isPion()) { 00489 (*i)->setEmissionTime(theStore->getBook()->getCurrentTime()); 00490 // Correction for real masses 00491 const G4double theQValueCorrection = (*i)->getEmissionQValueCorrection(theA,theZ); 00492 const G4double kineticEnergyOutside = (*i)->getKineticEnergy() - (*i)->getPotentialEnergy() + theQValueCorrection; 00493 (*i)->setTableMass(); 00494 if(kineticEnergyOutside > 0.0) 00495 (*i)->setEnergy((*i)->getMass()+kineticEnergyOutside); 00496 else 00497 (*i)->setEnergy((*i)->getMass()+tinyPionEnergy); 00498 (*i)->adjustMomentumFromEnergy(); 00499 (*i)->setPotentialEnergy(0.); 00500 theZ -= (*i)->getZ(); 00501 theStore->particleHasBeenEjected((*i)->getID()); 00502 theStore->addToOutgoing(*i); 00503 } 00504 } 00505 }
void G4INCL::Nucleus::fillEventInfo | ( | EventInfo * | eventInfo | ) |
Fill the event info which contains INCL output data
Definition at line 695 of file G4INCLNucleus.cc.
References G4INCL::EventInfo::A, G4INCL::EventInfo::ARem, G4INCL::CollisionAvatarType, G4INCL::DecayAvatarType, G4INCL::EventInfo::EKin, G4INCL::EventInfo::EKinRem, G4INCL::EventInfo::emissionTime, G4INCL::EventInfo::EStarRem, G4INCL::EventInfo::firstCollisionTime, G4INCL::EventInfo::firstCollisionXSec, G4INCL::EventInfo::forcedCompoundNucleus, G4INCL::Particle::getA(), G4INCL::Book::getAcceptedCollisions(), G4INCL::Book::getAcceptedDecays(), G4INCL::Book::getAvatars(), G4INCL::Book::getBlockedCollisions(), G4INCL::Book::getBlockedDecays(), G4INCL::Store::getBook(), getExcitationEnergy(), G4INCL::Book::getFirstCollisionTime(), G4INCL::Book::getFirstCollisionXSec(), G4INCL::Particle::getKineticEnergy(), G4INCL::Particle::getMomentum(), G4INCL::Store::getOutgoingParticles(), G4INCL::Cluster::getSpin(), getStore(), G4INCL::ThreeVector::getX(), G4INCL::ThreeVector::getY(), G4INCL::Particle::getZ(), G4INCL::ThreeVector::getZ(), hasRemnant(), G4INCL::PhysicalConstants::hc, G4INCL::EventInfo::history, G4INCL::EventInfo::JRem, G4INCL::EventInfo::jxRem, G4INCL::EventInfo::jyRem, G4INCL::EventInfo::jzRem, G4INCL::ThreeVector::mag(), G4INCL::EventInfo::nBlockedCollisions, G4INCL::EventInfo::nBlockedDecays, G4INCL::EventInfo::nCascadeParticles, G4INCL::EventInfo::nCollisionAvatars, G4INCL::EventInfo::nCollisions, G4INCL::EventInfo::nDecayAvatars, G4INCL::EventInfo::nDecays, G4INCL::Neutron, G4INCL::EventInfo::nParticles, G4INCL::EventInfo::nReflectionAvatars, G4INCL::EventInfo::nRemnants, G4INCL::EventInfo::nucleonAbsorption, G4INCL::EventInfo::origin, G4INCL::EventInfo::phi, G4INCL::ThreeVector::phi(), G4INCL::EventInfo::phiRem, G4INCL::PiMinus, G4INCL::EventInfo::pionAbsorption, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::EventInfo::projectileType, G4INCL::Proton, G4INCL::EventInfo::px, G4INCL::EventInfo::pxRem, G4INCL::EventInfo::py, G4INCL::EventInfo::pyRem, G4INCL::EventInfo::pz, G4INCL::EventInfo::pzRem, G4INCL::SurfaceAvatarType, G4INCL::EventInfo::theta, G4INCL::ThreeVector::theta(), G4INCL::EventInfo::thetaRem, G4INCL::Math::toDegrees(), WARN, G4INCL::EventInfo::Z, and G4INCL::EventInfo::ZRem.
00695 { 00696 eventInfo->nParticles = 0; 00697 G4bool isNucleonAbsorption = false; 00698 00699 G4bool isPionAbsorption = false; 00700 // It is possible to have pion absorption event only if the 00701 // projectile is pion. 00702 if(eventInfo->projectileType == PiPlus || 00703 eventInfo->projectileType == PiMinus || 00704 eventInfo->projectileType == PiZero) { 00705 isPionAbsorption = true; 00706 } 00707 00708 // Forced CN 00709 eventInfo->forcedCompoundNucleus = tryCN; 00710 00711 // Outgoing particles 00712 ParticleList outgoingParticles = getStore()->getOutgoingParticles(); 00713 00714 // Check if we have a nucleon absorption event: nucleon projectile 00715 // and no ejected particles. 00716 if(outgoingParticles.size() == 0 && 00717 (eventInfo->projectileType == Proton || 00718 eventInfo->projectileType == Neutron)) { 00719 isNucleonAbsorption = true; 00720 } 00721 00722 // Reset the remnant counter 00723 eventInfo->nRemnants = 0; 00724 eventInfo->history.clear(); 00725 00726 for( ParticleIter i = outgoingParticles.begin(); i != outgoingParticles.end(); ++i ) { 00727 // If the particle is a cluster and has excitation energy, treat it as a cluster 00728 if((*i)->isCluster()) { 00729 Cluster const * const c = dynamic_cast<Cluster *>(*i); 00730 // assert(c); 00731 #ifdef INCLXX_IN_GEANT4_MODE 00732 if(!c) 00733 continue; 00734 #endif 00735 const G4double eStar = c->getExcitationEnergy(); 00736 if(std::abs(eStar)>1E-10) { 00737 if(eStar<0.) { 00738 WARN("Negative excitation energy in outgoing cluster! EStar = " << eStar << std::endl); 00739 } 00740 eventInfo->ARem[eventInfo->nRemnants] = c->getA(); 00741 eventInfo->ZRem[eventInfo->nRemnants] = c->getZ(); 00742 eventInfo->EStarRem[eventInfo->nRemnants] = eStar; 00743 ThreeVector remnantSpin = c->getSpin(); 00744 Float_t remnantSpinMag; 00745 if(eventInfo->ARem[eventInfo->nRemnants]%2==0) { // even-A nucleus 00746 remnantSpinMag = (G4int) (remnantSpin.mag()/PhysicalConstants::hc + 0.5); 00747 } else { // odd-A nucleus 00748 remnantSpinMag = ((G4int) (remnantSpin.mag()/PhysicalConstants::hc)) + 0.5; 00749 } 00750 remnantSpin *= remnantSpinMag/remnantSpin.mag(); 00751 eventInfo->JRem[eventInfo->nRemnants] = remnantSpinMag; 00752 eventInfo->jxRem[eventInfo->nRemnants] = remnantSpin.getX(); 00753 eventInfo->jyRem[eventInfo->nRemnants] = remnantSpin.getY(); 00754 eventInfo->jzRem[eventInfo->nRemnants] = remnantSpin.getZ(); 00755 eventInfo->EKinRem[eventInfo->nRemnants] = c->getKineticEnergy(); 00756 ThreeVector mom = c->getMomentum(); 00757 eventInfo->pxRem[eventInfo->nRemnants] = mom.getX(); 00758 eventInfo->pyRem[eventInfo->nRemnants] = mom.getY(); 00759 eventInfo->pzRem[eventInfo->nRemnants] = mom.getZ(); 00760 eventInfo->thetaRem[eventInfo->nRemnants] = Math::toDegrees(mom.theta()); 00761 eventInfo->phiRem[eventInfo->nRemnants] = Math::toDegrees(mom.phi()); 00762 eventInfo->nRemnants++; 00763 continue; // don't add it as a particle 00764 } 00765 } 00766 00767 // We have a pion absorption event only if the projectile is 00768 // pion and there are no ejected pions. 00769 if(isPionAbsorption) { 00770 if((*i)->isPion()) { 00771 isPionAbsorption = false; 00772 } 00773 } 00774 eventInfo->A[eventInfo->nParticles] = (*i)->getA(); 00775 eventInfo->Z[eventInfo->nParticles] = (*i)->getZ(); 00776 eventInfo->emissionTime[eventInfo->nParticles] = (*i)->getEmissionTime(); 00777 eventInfo->EKin[eventInfo->nParticles] = (*i)->getKineticEnergy(); 00778 ThreeVector mom = (*i)->getMomentum(); 00779 eventInfo->px[eventInfo->nParticles] = mom.getX(); 00780 eventInfo->py[eventInfo->nParticles] = mom.getY(); 00781 eventInfo->pz[eventInfo->nParticles] = mom.getZ(); 00782 eventInfo->theta[eventInfo->nParticles] = Math::toDegrees(mom.theta()); 00783 eventInfo->phi[eventInfo->nParticles] = Math::toDegrees(mom.phi()); 00784 eventInfo->origin[eventInfo->nParticles] = -1; 00785 eventInfo->history.push_back(""); 00786 eventInfo->nParticles++; 00787 } 00788 eventInfo->nucleonAbsorption = isNucleonAbsorption; 00789 eventInfo->pionAbsorption = isPionAbsorption; 00790 eventInfo->nCascadeParticles = eventInfo->nParticles; 00791 00792 // Remnant characteristics 00793 if(hasRemnant()) { 00794 eventInfo->ARem[eventInfo->nRemnants] = getA(); 00795 eventInfo->ZRem[eventInfo->nRemnants] = getZ(); 00796 eventInfo->EStarRem[eventInfo->nRemnants] = getExcitationEnergy(); 00797 if(eventInfo->EStarRem[eventInfo->nRemnants]<0.) { 00798 WARN("Negative excitation energy! EStarRem = " << eventInfo->EStarRem[eventInfo->nRemnants] << std::endl); 00799 } 00800 if(eventInfo->ARem[eventInfo->nRemnants]%2==0) { // even-A nucleus 00801 eventInfo->JRem[eventInfo->nRemnants] = (G4int) (getSpin().mag()/PhysicalConstants::hc + 0.5); 00802 } else { // odd-A nucleus 00803 eventInfo->JRem[eventInfo->nRemnants] = ((G4int) (getSpin().mag()/PhysicalConstants::hc)) + 0.5; 00804 } 00805 eventInfo->EKinRem[eventInfo->nRemnants] = getKineticEnergy(); 00806 ThreeVector mom = getMomentum(); 00807 eventInfo->pxRem[eventInfo->nRemnants] = mom.getX(); 00808 eventInfo->pyRem[eventInfo->nRemnants] = mom.getY(); 00809 eventInfo->pzRem[eventInfo->nRemnants] = mom.getZ(); 00810 eventInfo->thetaRem[eventInfo->nRemnants] = Math::toDegrees(mom.theta()); 00811 eventInfo->phiRem[eventInfo->nRemnants] = Math::toDegrees(mom.phi()); 00812 eventInfo->nRemnants++; 00813 } 00814 00815 // Global counters, flags, etc. 00816 eventInfo->nCollisions = getStore()->getBook()->getAcceptedCollisions(); 00817 eventInfo->nBlockedCollisions = getStore()->getBook()->getBlockedCollisions(); 00818 eventInfo->nDecays = getStore()->getBook()->getAcceptedDecays(); 00819 eventInfo->nBlockedDecays = getStore()->getBook()->getBlockedDecays(); 00820 eventInfo->firstCollisionTime = getStore()->getBook()->getFirstCollisionTime(); 00821 eventInfo->firstCollisionXSec = getStore()->getBook()->getFirstCollisionXSec(); 00822 eventInfo->nReflectionAvatars = getStore()->getBook()->getAvatars(SurfaceAvatarType); 00823 eventInfo->nCollisionAvatars = getStore()->getBook()->getAvatars(CollisionAvatarType); 00824 eventInfo->nDecayAvatars = getStore()->getBook()->getAvatars(DecayAvatarType); 00825 }
void G4INCL::Nucleus::finalizeProjectileRemnant | ( | const G4double | emissionTime | ) |
Finalise the projectile remnant.
Complete the treatment of the projectile remnant. If it contains nucleons, assign its excitation energy and spin. Move stuff to the outgoing list, if appropriate.
emissionTime | the emission time of the projectile remnant |
Definition at line 869 of file G4INCLNucleus.cc.
References G4INCL::Store::addToOutgoing(), deleteProjectileRemnant(), G4INCL::Particle::getA(), G4INCL::Particle::getInvariantMass(), G4INCL::ProjectileRemnant::getNumberStoredComponents(), G4INCL::Cluster::getParticles(), G4INCL::ParticleTable::getTableMass, G4INCL::Particle::getZ(), G4INCL::Particle::setEmissionTime(), G4INCL::Cluster::setExcitationEnergy(), G4INCL::Particle::setMass(), G4INCL::Cluster::setSpin(), and G4INCL::DeJongSpin::shoot().
00869 { 00870 // Deal with the projectile remnant 00871 if(theProjectileRemnant->getA()>1) { 00872 // Set the mass 00873 const G4double aMass = theProjectileRemnant->getInvariantMass(); 00874 theProjectileRemnant->setMass(aMass); 00875 00876 // Compute the excitation energy from the invariant mass 00877 const G4double anExcitationEnergy = aMass 00878 - ParticleTable::getTableMass(theProjectileRemnant->getA(), theProjectileRemnant->getZ()); 00879 00880 // Set the excitation energy 00881 theProjectileRemnant->setExcitationEnergy(anExcitationEnergy); 00882 00883 // Set the spin 00884 theProjectileRemnant->setSpin(DeJongSpin::shoot(theProjectileRemnant->getNumberStoredComponents(), theProjectileRemnant->getA())); 00885 00886 // Set the emission time 00887 theProjectileRemnant->setEmissionTime(anEmissionTime); 00888 00889 // Put it in the outgoing list 00890 theStore->addToOutgoing(theProjectileRemnant); 00891 00892 // NULL theProjectileRemnant 00893 theProjectileRemnant = NULL; 00894 } else if(theProjectileRemnant->getA()==1) { 00895 // Put the nucleon in the outgoing list 00896 Particle *theNucleon = theProjectileRemnant->getParticles().front(); 00897 theStore->addToOutgoing(theNucleon); 00898 // Delete the remnant 00899 deleteProjectileRemnant(); 00900 } else 00901 deleteProjectileRemnant(); 00902 }
Particle* G4INCL::Nucleus::getBlockedDelta | ( | ) | const [inline] |
Get the delta that could not decay.
Definition at line 106 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::propagate().
Nucleus::ConservationBalance G4INCL::Nucleus::getConservationBalance | ( | EventInfo const & | theEventInfo, | |
const G4bool | afterRecoil | |||
) | const |
Compute charge, mass, energy and momentum balance.
Definition at line 827 of file G4INCLNucleus.cc.
References G4INCL::Nucleus::ConservationBalance::A, G4INCL::EventInfo::Ap, G4INCL::EventInfo::At, G4INCL::Nucleus::ConservationBalance::energy, G4INCL::Particle::getA(), getExcitationEnergy(), getIncomingMomentum(), getInitialEnergy(), G4INCL::Particle::getKineticEnergy(), G4INCL::Particle::getMomentum(), G4INCL::Store::getOutgoingParticles(), G4INCL::ParticleTable::getTableMass, G4INCL::Particle::getZ(), hasRemnant(), G4INCL::Nucleus::ConservationBalance::momentum, G4INCL::Nucleus::ConservationBalance::Z, G4INCL::EventInfo::Zp, and G4INCL::EventInfo::Zt.
00827 { 00828 ConservationBalance theBalance; 00829 // Initialise balance variables with the incoming values 00830 theBalance.Z = theEventInfo.Zp + theEventInfo.Zt; 00831 theBalance.A = theEventInfo.Ap + theEventInfo.At; 00832 00833 theBalance.energy = getInitialEnergy(); 00834 theBalance.momentum = getIncomingMomentum(); 00835 00836 // Process outgoing particles 00837 ParticleList outgoingParticles = theStore->getOutgoingParticles(); 00838 for( ParticleIter i = outgoingParticles.begin(); i != outgoingParticles.end(); ++i ) { 00839 theBalance.Z -= (*i)->getZ(); 00840 theBalance.A -= (*i)->getA(); 00841 // For outgoing clusters, the total energy automatically includes the 00842 // excitation energy 00843 theBalance.energy -= (*i)->getEnergy(); // Note that outgoing particles should have the real mass 00844 theBalance.momentum -= (*i)->getMomentum(); 00845 } 00846 00847 // Remnant contribution, if present 00848 if(hasRemnant()) { 00849 theBalance.Z -= getZ(); 00850 theBalance.A -= getA(); 00851 theBalance.energy -= ParticleTable::getTableMass(getA(),getZ()) + 00852 getExcitationEnergy(); 00853 if(afterRecoil) 00854 theBalance.energy -= getKineticEnergy(); 00855 theBalance.momentum -= getMomentum(); 00856 } 00857 00858 return theBalance; 00859 }
G4double G4INCL::Nucleus::getCoulombRadius | ( | ParticleSpecies const & | p | ) | const [inline] |
Get the Coulomb radius for a given particle.
That's the radius of the sphere that the Coulomb trajectory of the incoming particle should intersect. The intersection point is used to determine the effective impact parameter of the trajectory and the new entrance angle.
If the particle is not a Cluster, the Coulomb radius reduces to the surface radius. We use a parametrisation for d, t, He3 and alphas. For heavier clusters we fall back to the surface radius.
p | the particle species |
Definition at line 350 of file G4INCLNucleus.hh.
References G4InuclParticleNames::ap, G4INCL::Composite, ERROR, G4INCL::ParticleTable::eSquared, getUniverseRadius(), G4INCL::Math::pow23(), G4INCL::Particle::theA, G4INCL::ParticleSpecies::theA, G4INCL::ParticleSpecies::theType, G4INCL::Particle::theZ, and G4INCL::ParticleSpecies::theZ.
00350 { 00351 if(p.theType == Composite) { 00352 const G4int zp = p.theZ; 00353 const G4int ap = p.theA; 00354 G4double barr, radius = 0.; 00355 if(zp==1 && ap==2) { // d 00356 barr = 0.2565*Math::pow23((G4double)theA)-0.78; 00357 radius = ParticleTable::eSquared*zp*theZ/barr - 2.5; 00358 } else if((zp==1 || zp==2) && ap==3) { // t, He3 00359 barr = 0.5*(0.5009*Math::pow23((G4double)theA)-1.16); 00360 radius = ParticleTable::eSquared*theZ/barr - 0.5; 00361 } else if(zp==2 && ap==4) { // alpha 00362 barr = 0.5939*Math::pow23((G4double)theA)-1.64; 00363 radius = ParticleTable::eSquared*zp*theZ/barr - 0.5; 00364 } else if(zp>2) { 00365 radius = getUniverseRadius(); 00366 } 00367 if(radius<=0.) { 00368 ERROR("Negative Coulomb radius! Using the universe radius" << std::endl); 00369 radius = getUniverseRadius(); 00370 } 00371 return radius; 00372 } else 00373 return getUniverseRadius(); 00374 }
ParticleList const& G4INCL::Nucleus::getCreatedParticles | ( | ) | const [inline] |
Get the list of particles that were created by the last applied final state
Definition at line 98 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::propagate().
NuclearDensity* G4INCL::Nucleus::getDensity | ( | ) | const [inline] |
Getter for theDensity.
Definition at line 428 of file G4INCLNucleus.hh.
Referenced by G4INCL::CoulombNonRelativistic::distortOut(), G4INCL::PauliStandard::getBlockingProbability(), and G4INCL::ClusteringModelIntercomparison::getCluster().
G4double G4INCL::Nucleus::getExcitationEnergy | ( | ) | const [inline] |
Get the excitation energy of the nucleus.
Method computeRecoilKinematics() should be called first.
Reimplemented from G4INCL::Cluster.
Definition at line 234 of file G4INCLNucleus.hh.
References G4INCL::Cluster::theExcitationEnergy.
Referenced by fillEventInfo(), and getConservationBalance().
00234 { return theExcitationEnergy; }
const ThreeVector& G4INCL::Nucleus::getIncomingAngularMomentum | ( | ) | const [inline] |
const ThreeVector& G4INCL::Nucleus::getIncomingMomentum | ( | ) | const [inline] |
Get the incoming momentum vector.
Definition at line 220 of file G4INCLNucleus.hh.
Referenced by getConservationBalance().
G4int G4INCL::Nucleus::getInitialA | ( | ) | const [inline] |
G4double G4INCL::Nucleus::getInitialEnergy | ( | ) | const [inline] |
Get the initial energy.
Definition at line 228 of file G4INCLNucleus.hh.
Referenced by getConservationBalance().
G4double G4INCL::Nucleus::getInitialInternalEnergy | ( | ) | const [inline] |
G4int G4INCL::Nucleus::getInitialZ | ( | ) | const [inline] |
G4int G4INCL::Nucleus::getNumberOfEnteringNeutrons | ( | ) | const [inline] |
G4int G4INCL::Nucleus::getNumberOfEnteringProtons | ( | ) | const [inline] |
NuclearPotential::INuclearPotential* G4INCL::Nucleus::getPotential | ( | ) | const [inline] |
Getter for thePotential.
Definition at line 431 of file G4INCLNucleus.hh.
Referenced by G4INCL::PauliStandard::getBlockingProbability(), G4INCL::SurfaceAvatar::getChannel(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ParticleEntryChannel::getFinalState(), and G4INCL::CDPP::isBlocked().
G4int G4INCL::Nucleus::getProjectileChargeNumber | ( | ) | const [inline] |
G4int G4INCL::Nucleus::getProjectileMassNumber | ( | ) | const [inline] |
ProjectileRemnant* G4INCL::Nucleus::getProjectileRemnant | ( | ) | const [inline] |
Get the projectile remnant.
Definition at line 398 of file G4INCLNucleus.hh.
Referenced by G4INCL::ParticleEntryChannel::getFinalState().
Store* G4INCL::Nucleus::getStore | ( | ) | const [inline] |
Definition at line 251 of file G4INCLNucleus.hh.
Referenced by fillEventInfo(), G4INCL::StandardPropagationModel::generateAllAvatars(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::INCL::initializeTarget(), G4INCL::PauliStrictStandard::isBlocked(), G4INCL::CDPP::isBlocked(), G4INCL::SurfaceAvatar::postInteraction(), G4INCL::InteractionAvatar::postInteraction(), G4INCL::DecayAvatar::postInteraction(), G4INCL::BinaryCollisionAvatar::postInteraction(), G4INCL::StandardPropagationModel::propagate(), G4INCL::StandardPropagationModel::registerAvatar(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::InteractionAvatar::shouldUseLocalEnergy(), and G4INCL::StandardPropagationModel::updateAvatars().
Get the maximum allowed radius for a given particle.
Calls the NuclearDensity::getMaxRFromP() method for nucleons and deltas, and the NuclearDensity::getTrasmissionRadius() method for pions.
particle | pointer to a particle |
Definition at line 322 of file G4INCLNucleus.hh.
References G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::NuclearDensity::getMaxRFromP(), G4INCL::Particle::getMomentum(), getUniverseRadius(), G4INCL::Particle::isPion(), and G4INCL::ThreeVector::mag().
Referenced by G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::StandardPropagationModel::getReflectionTime(), and G4INCL::InteractionAvatar::postInteraction().
00322 { 00323 if(particle->isPion()) 00324 // Temporarily set RPION = RMAX 00325 return getUniverseRadius(); 00326 //return 0.5*(theDensity->getTransmissionRadius(particle)+getUniverseRadius()); 00327 else { 00328 const G4double pr = particle->getMomentum().mag()/thePotential->getFermiMomentum(particle); 00329 if(pr>=1.) 00330 return getUniverseRadius(); 00331 else 00332 return theDensity->getMaxRFromP(pr); 00333 } 00334 }
Get the transmission barrier.
Definition at line 295 of file G4INCLNucleus.hh.
References G4INCL::PhysicalConstants::eSquared, G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::Particle::getZ(), and G4INCL::Particle::theZ.
Referenced by G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::InteractionAvatar::postInteraction().
00295 { 00296 const G4double theTransmissionRadius = theDensity->getTransmissionRadius(p); 00297 const G4double theParticleZ = p->getZ(); 00298 return PhysicalConstants::eSquared*(theZ-theParticleZ)*theParticleZ/theTransmissionRadius; 00299 }
G4bool G4INCL::Nucleus::getTryCompoundNucleus | ( | ) | [inline] |
G4double G4INCL::Nucleus::getUniverseRadius | ( | ) | const [inline] |
Getter for theUniverseRadius.
Definition at line 377 of file G4INCLNucleus.hh.
Referenced by G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), getCoulombRadius(), getSurfaceRadius(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
ParticleList const& G4INCL::Nucleus::getUpdatedParticles | ( | ) | const [inline] |
Get the list of particles that were updated by the last applied final state
Definition at line 103 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::generateAllAvatars(), and G4INCL::StandardPropagationModel::propagate().
G4bool G4INCL::Nucleus::hasRemnant | ( | ) | const [inline] |
Does the nucleus give a cascade remnant?
To be called after computeRecoilKinematics().
Definition at line 269 of file G4INCLNucleus.hh.
Referenced by fillEventInfo(), and getConservationBalance().
void G4INCL::Nucleus::initializeParticles | ( | ) | [virtual] |
Call the Cluster method to generate the initial distribution of particles. At the beginning all particles are assigned as spectators.
Reimplemented from G4INCL::Cluster.
Definition at line 137 of file G4INCLNucleus.cc.
References G4INCL::Store::add(), computeTotalEnergy(), G4INCL::Cluster::initializeParticles(), G4INCL::Cluster::particles, G4INCL::Particle::thePosition, and updatePotentialEnergy().
Referenced by G4INCL::INCL::initializeTarget().
00137 { 00138 // Reset the variables connected with the projectile remnant 00139 delete theProjectileRemnant; 00140 theProjectileRemnant = NULL; 00141 00142 Cluster::initializeParticles(); 00143 for(ParticleIter i = particles.begin(); i != particles.end(); ++i) { 00144 updatePotentialEnergy(*i); 00145 theStore->add(*i); 00146 } 00147 particles.clear(); 00148 initialInternalEnergy = computeTotalEnergy(); 00149 initialCenterOfMass = thePosition; 00150 }
void G4INCL::Nucleus::insertParticle | ( | Particle * | p | ) | [inline] |
Insert a new particle (e.g. a projectile) in the nucleus.
Definition at line 76 of file G4INCLNucleus.hh.
References G4INCL::Particle::getA(), G4INCL::Store::getBook(), G4INCL::ParticleTable::getIsospin(), G4INCL::Particle::getType(), G4INCL::Particle::getZ(), G4INCL::Math::heaviside(), G4INCL::Book::incrementCascading(), G4INCL::Particle::isNucleon(), G4INCL::Particle::isTargetSpectator(), G4INCL::Store::particleHasEntered(), G4INCL::Particle::theA, and G4INCL::Particle::theZ.
Referenced by applyFinalState().
00076 { 00077 theZ += p->getZ(); 00078 theA += p->getA(); 00079 theStore->particleHasEntered(p); 00080 if(p->isNucleon()) { 00081 theNpInitial += Math::heaviside(ParticleTable::getIsospin(p->getType())); 00082 theNnInitial += Math::heaviside(-ParticleTable::getIsospin(p->getType())); 00083 } 00084 if(!p->isTargetSpectator()) theStore->getBook()->incrementCascading(); 00085 };
G4bool G4INCL::Nucleus::isEventTransparent | ( | ) | const |
Is the event transparent?
To be called at the end of the cascade.
Definition at line 507 of file G4INCLNucleus.cc.
References G4INCL::Particle::getA(), G4INCL::Store::getOutgoingParticles(), and G4INCL::Particle::getZ().
00507 { 00508 00509 // Forced transparent 00510 if(forceTransparent) 00511 return true; 00512 00513 ParticleList const &pL = theStore->getOutgoingParticles(); 00514 G4int outZ = 0, outA = 0; 00515 00516 // If any of the particles has undergone a collision, the event is not a 00517 // transparent. 00518 for(ParticleIter p = pL.begin(); p != pL.end(); ++p ) { 00519 if( (*p)->getNumberOfCollisions() != 0 ) return false; 00520 if( (*p)->getNumberOfDecays() != 0 ) return false; 00521 outZ += (*p)->getZ(); 00522 outA += (*p)->getA(); 00523 } 00524 00525 // Add the geometrical spectators to the Z and A count 00526 if(theProjectileRemnant) { 00527 outZ += theProjectileRemnant->getZ(); 00528 outA += theProjectileRemnant->getA(); 00529 } 00530 00531 if(outZ!=projectileZ || outA!=projectileA) return false; 00532 00533 return true; 00534 00535 }
G4bool G4INCL::Nucleus::isForcedTransparent | ( | ) | [inline] |
Returns true if a transparent event should be forced.
Definition at line 292 of file G4INCLNucleus.hh.
G4bool G4INCL::Nucleus::isNucleusNucleusCollision | ( | ) | const [inline] |
Is it a nucleus-nucleus collision?
Definition at line 383 of file G4INCLNucleus.hh.
Referenced by G4INCL::SurfaceAvatar::getChannel(), and G4INCL::ParticleEntryChannel::getFinalState().
void G4INCL::Nucleus::moveProjectileRemnantComponentsToOutgoing | ( | ) | [inline] |
Move the components of the projectile remnant to the outgoing list.
Definition at line 407 of file G4INCLNucleus.hh.
References G4INCL::Store::addToOutgoing(), G4INCL::Cluster::clearParticles(), and G4INCL::Cluster::getParticles().
00407 { 00408 theStore->addToOutgoing(theProjectileRemnant->getParticles()); 00409 theProjectileRemnant->clearParticles(); 00410 }
std::string G4INCL::Nucleus::print | ( | ) |
Print the nucleus info
Definition at line 315 of file G4INCLNucleus.cc.
References G4INCL::Store::getOutgoingParticles(), G4INCL::Store::getParticipants(), and G4INCL::Store::getSpectators().
00316 { 00317 std::stringstream ss; 00318 ss << "Particles in the nucleus:" << std::endl 00319 << "Participants:" << std::endl; 00320 G4int counter = 1; 00321 ParticleList participants = theStore->getParticipants(); 00322 for(ParticleIter p = participants.begin(); p != participants.end(); ++p) { 00323 ss << "index = " << counter << std::endl 00324 << (*p)->print(); 00325 counter++; 00326 } 00327 ss <<"Spectators:" << std::endl; 00328 ParticleList spectators = theStore->getSpectators(); 00329 for(ParticleIter p = spectators.begin(); p != spectators.end(); ++p) 00330 ss << (*p)->print(); 00331 ss <<"Outgoing:" << std::endl; 00332 ParticleList outgoing = theStore->getOutgoingParticles(); 00333 for(ParticleIter p = outgoing.begin(); p != outgoing.end(); ++p) 00334 ss << (*p)->print(); 00335 00336 return ss.str(); 00337 }
void G4INCL::Nucleus::propagateParticles | ( | G4double | step | ) |
Propagate the particles one time step.
step | length of the time step |
Definition at line 247 of file G4INCLNucleus.cc.
References WARN.
00247 { 00248 WARN("Useless Nucleus::propagateParticles -method called." << std::endl); 00249 }
void G4INCL::Nucleus::setIncomingAngularMomentum | ( | const ThreeVector & | j | ) | [inline] |
Set the incoming angular-momentum vector.
Definition at line 207 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
void G4INCL::Nucleus::setIncomingMomentum | ( | const ThreeVector & | p | ) | [inline] |
Set the incoming momentum vector.
Definition at line 215 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
void G4INCL::Nucleus::setInitialEnergy | ( | const G4double | e | ) | [inline] |
Set the initial energy.
Definition at line 225 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
void G4INCL::Nucleus::setNucleusNucleusCollision | ( | ) | [inline] |
Set a nucleus-nucleus collision.
Definition at line 386 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite().
void G4INCL::Nucleus::setParticleNucleusCollision | ( | ) | [inline] |
Set a particle-nucleus collision.
Definition at line 389 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootParticle().
void G4INCL::Nucleus::setProjectileChargeNumber | ( | G4int | n | ) | [inline] |
Set the charge number of the projectile.
Definition at line 286 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
00286 { projectileZ=n; }
void G4INCL::Nucleus::setProjectileMassNumber | ( | G4int | n | ) | [inline] |
Set the mass number of the projectile.
Definition at line 289 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
00289 { projectileA=n; }
void G4INCL::Nucleus::setProjectileRemnant | ( | ProjectileRemnant *const | c | ) | [inline] |
Set the projectile remnant.
Definition at line 392 of file G4INCLNucleus.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite().
void G4INCL::Nucleus::setStore | ( | Store * | s | ) | [inline] |
void G4INCL::Nucleus::setUniverseRadius | ( | const G4double | universeRadius | ) | [inline] |
void G4INCL::Nucleus::updatePotentialEnergy | ( | Particle * | p | ) | [inline] |
Update the particle potential energy.
Definition at line 423 of file G4INCLNucleus.hh.
References G4INCL::NuclearPotential::INuclearPotential::computePotentialEnergy(), and G4INCL::Particle::setPotentialEnergy().
Referenced by G4INCL::ReflectionChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), and initializeParticles().
void G4INCL::Nucleus::useFusionKinematics | ( | ) |
Adjust the kinematics for complete-fusion events.
Definition at line 861 of file G4INCLNucleus.cc.
References G4INCL::Cluster::getTableMass(), G4INCL::ThreeVector::mag2(), G4INCL::Particle::setEnergy(), G4INCL::Particle::setMass(), G4INCL::Particle::setMomentum(), G4INCL::Cluster::setSpin(), G4INCL::Particle::theEnergy, G4INCL::Cluster::theExcitationEnergy, and G4INCL::Particle::theMomentum.
00861 { 00862 setEnergy(initialEnergy); 00863 setMomentum(incomingMomentum); 00864 setSpin(incomingAngularMomentum); 00865 theExcitationEnergy = std::sqrt(theEnergy*theEnergy-theMomentum.mag2()) - getTableMass(); 00866 setMass(getTableMass() + theExcitationEnergy); 00867 }