#include <G4WilsonAbrasionModel.hh>
Inheritance diagram for G4WilsonAbrasionModel:
Public Member Functions | |
G4WilsonAbrasionModel (G4bool useAblation1=false) | |
G4WilsonAbrasionModel (G4ExcitationHandler *) | |
~G4WilsonAbrasionModel () | |
G4WilsonAbrasionModel (const G4WilsonAbrasionModel &right) | |
const G4WilsonAbrasionModel & | operator= (G4WilsonAbrasionModel &right) |
virtual G4HadFinalState * | ApplyYourself (const G4HadProjectile &, G4Nucleus &) |
void | SetVerboseLevel (G4int) |
void | SetUseAblation (G4bool) |
G4bool | GetUseAblation () |
void | SetConserveMomentum (G4bool) |
G4bool | GetConserveMomentum () |
void | SetExcitationHandler (G4ExcitationHandler *) |
G4ExcitationHandler * | GetExcitationHandler () |
virtual void | ModelDescription (std::ostream &) const |
Definition at line 77 of file G4WilsonAbrasionModel.hh.
G4WilsonAbrasionModel::G4WilsonAbrasionModel | ( | G4bool | useAblation1 = false |
) |
Definition at line 108 of file G4WilsonAbrasionModel.cc.
References G4HadronicInteraction::isBlocked, G4ExcitationHandler::SetEvaporation(), G4ExcitationHandler::SetFermiModel(), G4ExcitationHandler::SetMaxAandZForFermiBreakUp(), G4HadronicInteraction::SetMaxEnergy(), G4ExcitationHandler::SetMinEForMultiFrag(), G4HadronicInteraction::SetMinEnergy(), G4ExcitationHandler::SetMultiFragmentation(), G4WilsonAblationModel::SetVerboseLevel(), and G4HadronicInteraction::verboseLevel.
00109 :G4HadronicInteraction("G4WilsonAbrasion") 00110 { 00111 // Send message to stdout to advise that the G4Abrasion model is being used. 00112 PrintWelcomeMessage(); 00113 00114 // Set the default verbose level to 0 - no output. 00115 verboseLevel = 0; 00116 useAblation = useAblation1; 00117 00118 // No de-excitation handler has been supplied - define the default handler. 00119 00120 theExcitationHandler = new G4ExcitationHandler; 00121 theExcitationHandlerx = new G4ExcitationHandler; 00122 if (useAblation) 00123 { 00124 theAblation = new G4WilsonAblationModel; 00125 theAblation->SetVerboseLevel(verboseLevel); 00126 theExcitationHandler->SetEvaporation(theAblation); 00127 theExcitationHandlerx->SetEvaporation(theAblation); 00128 } 00129 else 00130 { 00131 theAblation = NULL; 00132 G4Evaporation * theEvaporation = new G4Evaporation; 00133 G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp; 00134 G4StatMF * theMF = new G4StatMF; 00135 theExcitationHandler->SetEvaporation(theEvaporation); 00136 theExcitationHandler->SetFermiModel(theFermiBreakUp); 00137 theExcitationHandler->SetMultiFragmentation(theMF); 00138 theExcitationHandler->SetMaxAandZForFermiBreakUp(12, 6); 00139 theExcitationHandler->SetMinEForMultiFrag(5.0*MeV); 00140 00141 theEvaporation = new G4Evaporation; 00142 theFermiBreakUp = new G4FermiBreakUp; 00143 theExcitationHandlerx->SetEvaporation(theEvaporation); 00144 theExcitationHandlerx->SetFermiModel(theFermiBreakUp); 00145 theExcitationHandlerx->SetMaxAandZForFermiBreakUp(12, 6); 00146 } 00147 00148 // Set the minimum and maximum range for the model (despite nomanclature, 00149 // this is in energy per nucleon number). 00150 00151 SetMinEnergy(70.0*MeV); 00152 SetMaxEnergy(10.1*GeV); 00153 isBlocked = false; 00154 00155 // npK, when mutiplied by the nuclear Fermi momentum, determines the range of 00156 // momentum over which the secondary nucleon momentum is sampled. 00157 00158 r0sq = 0.0; 00159 npK = 5.0; 00160 B = 10.0 * MeV; 00161 third = 1.0 / 3.0; 00162 fradius = 0.99; 00163 conserveEnergy = false; 00164 conserveMomentum = true; 00165 }
G4WilsonAbrasionModel::G4WilsonAbrasionModel | ( | G4ExcitationHandler * | ) |
Definition at line 179 of file G4WilsonAbrasionModel.cc.
References G4HadronicInteraction::isBlocked, G4ExcitationHandler::SetEvaporation(), G4ExcitationHandler::SetFermiModel(), G4ExcitationHandler::SetMaxAandZForFermiBreakUp(), G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), and G4HadronicInteraction::verboseLevel.
00180 { 00181 // Send message to stdout to advise that the G4Abrasion model is being used. 00182 00183 PrintWelcomeMessage(); 00184 00185 // Set the default verbose level to 0 - no output. 00186 00187 verboseLevel = 0; 00188 00189 theAblation = NULL; //A.R. 26-Jul-2012 Coverity fix. 00190 useAblation = false; //A.R. 14-Aug-2012 Coverity fix. 00191 00192 // 00193 // The user is able to provide the excitation handler as well as an argument 00194 // which is provided in this instantiation is used to determine 00195 // whether the spectators of the interaction are free following the abrasion. 00196 // 00197 theExcitationHandler = aExcitationHandler; 00198 theExcitationHandlerx = new G4ExcitationHandler; 00199 G4Evaporation * theEvaporation = new G4Evaporation; 00200 G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp; 00201 theExcitationHandlerx->SetEvaporation(theEvaporation); 00202 theExcitationHandlerx->SetFermiModel(theFermiBreakUp); 00203 theExcitationHandlerx->SetMaxAandZForFermiBreakUp(12, 6); 00204 // 00205 // 00206 // Set the minimum and maximum range for the model (despite nomanclature, this 00207 // is in energy per nucleon number). 00208 // 00209 SetMinEnergy(70.0*MeV); 00210 SetMaxEnergy(10.1*GeV); 00211 isBlocked = false; 00212 // 00213 // 00214 // npK, when mutiplied by the nuclear Fermi momentum, determines the range of 00215 // momentum over which the secondary nucleon momentum is sampled. 00216 // 00217 r0sq = 0.0; //A.R. 14-Aug-2012 Coverity fix. 00218 npK = 5.0; 00219 B = 10.0 * MeV; 00220 third = 1.0 / 3.0; 00221 fradius = 0.99; 00222 conserveEnergy = false; 00223 conserveMomentum = true; 00224 }
G4WilsonAbrasionModel::~G4WilsonAbrasionModel | ( | ) |
Definition at line 227 of file G4WilsonAbrasionModel.cc.
00228 { 00229 // 00230 // 00231 // The destructor doesn't have to do a great deal! 00232 // 00233 if (theExcitationHandler) delete theExcitationHandler; 00234 if (theExcitationHandlerx) delete theExcitationHandlerx; 00235 if (useAblation) delete theAblation; 00236 // delete theExcitationHandler; 00237 // delete theExcitationHandlerx; 00238 }
G4WilsonAbrasionModel::G4WilsonAbrasionModel | ( | const G4WilsonAbrasionModel & | right | ) |
G4HadFinalState * G4WilsonAbrasionModel::ApplyYourself | ( | const G4HadProjectile & | , | |
G4Nucleus & | ||||
) | [virtual] |
Implements G4HadronicInteraction.
Definition at line 241 of file G4WilsonAbrasionModel.cc.
References G4HadFinalState::AddSecondary(), G4Nucleus::AtomicMass(), G4ExcitationHandler::BreakItUp(), G4HadFinalState::Clear(), G4DynamicParticle::DumpInfo(), G4NuclearAbrasionGeometry::F(), G4cout, G4endl, G4Poisson(), G4UniformRand, G4DynamicParticle::Get4Momentum(), G4HadProjectile::Get4Momentum(), G4Fragment::GetA(), G4Nucleus::GetA_asInt(), G4ParticleDefinition::GetBaryonNumber(), G4DynamicParticle::GetDefinition(), G4HadProjectile::GetDefinition(), G4Nucleus::GetEnergyDeposit(), G4NuclearAbrasionGeometry::GetExcitationEnergyOfProjectile(), G4NuclearAbrasionGeometry::GetExcitationEnergyOfTarget(), G4Fragment::GetGroundStateMass(), G4HadProjectile::GetKineticEnergy(), G4Fragment::GetMomentum(), G4HadFinalState::GetNumberOfSecondaries(), G4HadSecondary::GetParticle(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4HadFinalState::GetSecondary(), G4HadProjectile::GetTotalEnergy(), G4WilsonRadius::GetWilsonRadius(), G4Fragment::GetZ(), G4Nucleus::GetZ_asInt(), isAlive, G4InuclParticleNames::lambda, CLHEP::detail::n, G4HadFinalState::SetEnergyChange(), G4Fragment::SetMomentum(), G4HadFinalState::SetMomentumChange(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, and G4HadronicInteraction::verboseLevel.
00243 { 00244 // 00245 // 00246 // The secondaries will be returned in G4HadFinalState &theParticleChange - 00247 // initialise this. The original track will always be discontinued and 00248 // secondaries followed. 00249 // 00250 theParticleChange.Clear(); 00251 theParticleChange.SetStatusChange(stopAndKill); 00252 // 00253 // 00254 // Get relevant information about the projectile and target (A, Z, energy/nuc, 00255 // momentum, etc). 00256 // 00257 const G4ParticleDefinition *definitionP = theTrack.GetDefinition(); 00258 const G4double AP = definitionP->GetBaryonNumber(); 00259 const G4double ZP = definitionP->GetPDGCharge(); 00260 G4LorentzVector pP = theTrack.Get4Momentum(); 00261 G4double E = theTrack.GetKineticEnergy()/AP; 00262 G4double AT = theTarget.GetA_asInt(); 00263 G4double ZT = theTarget.GetZ_asInt(); 00264 G4double TotalEPre = theTrack.GetTotalEnergy() + 00265 theTarget.AtomicMass(AT, ZT) + theTarget.GetEnergyDeposit(); 00266 G4double TotalEPost = 0.0; 00267 // 00268 // 00269 // Determine the radii of the projectile and target nuclei. 00270 // 00271 G4WilsonRadius aR; 00272 G4double rP = aR.GetWilsonRadius(AP); 00273 G4double rT = aR.GetWilsonRadius(AT); 00274 G4double rPsq = rP * rP; 00275 G4double rTsq = rT * rT; 00276 if (verboseLevel >= 2) 00277 { 00278 G4cout <<"########################################" 00279 <<"########################################" 00280 <<G4endl; 00281 G4cout.precision(6); 00282 G4cout <<"IN G4WilsonAbrasionModel" <<G4endl; 00283 G4cout <<"Initial projectile A=" <<AP 00284 <<", Z=" <<ZP 00285 <<", radius = " <<rP/fermi <<" fm" 00286 <<G4endl; 00287 G4cout <<"Initial target A=" <<AT 00288 <<", Z=" <<ZT 00289 <<", radius = " <<rT/fermi <<" fm" 00290 <<G4endl; 00291 G4cout <<"Projectile momentum and Energy/nuc = " <<pP <<" ," <<E <<G4endl; 00292 } 00293 // 00294 // 00295 // The following variables are used to determine the impact parameter in the 00296 // near-field (i.e. taking into consideration the electrostatic repulsion). 00297 // 00298 G4double rm = ZP * ZT * elm_coupling / (E * AP); 00299 G4double r = 0.0; 00300 G4double rsq = 0.0; 00301 // 00302 // 00303 // Initialise some of the variables which wll be used to calculate the chord- 00304 // length for nucleons in the projectile and target, and hence calculate the 00305 // number of abraded nucleons and the excitation energy. 00306 // 00307 G4NuclearAbrasionGeometry *theAbrasionGeometry = NULL; 00308 G4double CT = 0.0; 00309 G4double F = 0.0; 00310 G4int Dabr = 0; 00311 // 00312 // 00313 // The following loop is performed until the number of nucleons which are 00314 // abraded by the process is >1, i.e. an interaction MUST occur. 00315 // 00316 while (Dabr == 0) 00317 { 00318 // Added by MHM 20050119 to fix leaking memory on second pass through this loop 00319 if (theAbrasionGeometry) 00320 { 00321 delete theAbrasionGeometry; 00322 theAbrasionGeometry = NULL; 00323 } 00324 // 00325 // 00326 // Sample the impact parameter. For the moment, this class takes account of 00327 // electrostatic effects on the impact parameter, but (like HZETRN AND NUCFRG2) 00328 // does not make any correction for the effects of nuclear-nuclear repulsion. 00329 // 00330 G4double rPT = rP + rT; 00331 G4double rPTsq = rPT * rPT; 00332 // 00333 // 00334 // This is a "catch" to make sure we don't go into an infinite loop because the 00335 // energy is too low to overcome nuclear repulsion. PRT 20091023. If the 00336 // value of rm < fradius * rPT then we're unlikely to sample a small enough 00337 // impact parameter (energy of incident particle is too low). Return primary 00338 // and don't complete nuclear interaction analysis. 00339 // 00340 if (rm >= fradius * rPT) { 00341 theParticleChange.SetStatusChange(isAlive); 00342 theParticleChange.SetEnergyChange(theTrack.GetKineticEnergy()); 00343 theParticleChange.SetMomentumChange(theTrack.Get4Momentum().vect().unit()); 00344 if (verboseLevel >= 2) { 00345 G4cout <<"Particle energy too low to overcome repulsion." <<G4endl; 00346 G4cout <<"Event rejected and original track maintained" <<G4endl; 00347 G4cout <<"########################################" 00348 <<"########################################" 00349 <<G4endl; 00350 } 00351 return &theParticleChange; 00352 } 00353 // 00354 // 00355 // Now sample impact parameter until the criterion is met that projectile 00356 // and target overlap, but repulsion is taken into consideration. 00357 // 00358 G4int evtcnt = 0; 00359 r = 1.1 * rPT; 00360 while (r > rPT && ++evtcnt < 1000) 00361 { 00362 G4double bsq = rPTsq * G4UniformRand(); 00363 r = (rm + std::sqrt(rm*rm + 4.0*bsq)) / 2.0; 00364 } 00365 // 00366 // 00367 // We've tried to sample this 1000 times, but failed. Assume nuclei do not 00368 // collide. 00369 // 00370 if (evtcnt >= 1000) { 00371 theParticleChange.SetStatusChange(isAlive); 00372 theParticleChange.SetEnergyChange(theTrack.GetKineticEnergy()); 00373 theParticleChange.SetMomentumChange(theTrack.Get4Momentum().vect().unit()); 00374 if (verboseLevel >= 2) { 00375 G4cout <<"Particle energy too low to overcome repulsion." <<G4endl; 00376 G4cout <<"Event rejected and original track maintained" <<G4endl; 00377 G4cout <<"########################################" 00378 <<"########################################" 00379 <<G4endl; 00380 } 00381 return &theParticleChange; 00382 } 00383 00384 00385 rsq = r * r; 00386 // 00387 // 00388 // Now determine the chord-length through the target nucleus. 00389 // 00390 if (rT > rP) 00391 { 00392 G4double x = (rPsq + rsq - rTsq) / 2.0 / r; 00393 if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x); 00394 else CT = 2.0 * std::sqrt(rTsq - rsq); 00395 } 00396 else 00397 { 00398 G4double x = (rTsq + rsq - rPsq) / 2.0 / r; 00399 if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x); 00400 else CT = 2.0 * rT; 00401 } 00402 // 00403 // 00404 // Determine the number of abraded nucleons. Note that the mean number of 00405 // abraded nucleons is used to sample the Poisson distribution. The Poisson 00406 // distribution is sampled only ten times with the current impact parameter, 00407 // and if it fails after this to find a case for which the number of abraded 00408 // nucleons >1, the impact parameter is re-sampled. 00409 // 00410 theAbrasionGeometry = new G4NuclearAbrasionGeometry(AP,AT,r); 00411 F = theAbrasionGeometry->F(); 00412 G4double lambda = 16.6*fermi / std::pow(E/MeV,0.26); 00413 G4double Mabr = F * AP * (1.0 - std::exp(-CT/lambda)); 00414 G4long n = 0; 00415 for (G4int i = 0; i<10; i++) 00416 { 00417 n = G4Poisson(Mabr); 00418 if (n > 0) 00419 { 00420 if (n>AP) Dabr = (G4int) AP; 00421 else Dabr = (G4int) n; 00422 break; 00423 } 00424 } 00425 } 00426 if (verboseLevel >= 2) 00427 { 00428 G4cout <<G4endl; 00429 G4cout <<"Impact parameter = " <<r/fermi <<" fm" <<G4endl; 00430 G4cout <<"# Abraded nucleons = " <<Dabr <<G4endl; 00431 } 00432 // 00433 // 00434 // The number of abraded nucleons must be no greater than the number of 00435 // nucleons in either the projectile or the target. If AP - Dabr < 2 or 00436 // AT - Dabr < 2 then either we have only a nucleon left behind in the 00437 // projectile/target or we've tried to abrade too many nucleons - and Dabr 00438 // should be limited. 00439 // 00440 if (AP - (G4double) Dabr < 2.0) Dabr = (G4int) AP; 00441 if (AT - (G4double) Dabr < 2.0) Dabr = (G4int) AT; 00442 // 00443 // 00444 // Determine the abraded secondary nucleons from the projectile. *fragmentP 00445 // is a pointer to the prefragment from the projectile and nSecP is the number 00446 // of nucleons in theParticleChange which have been abraded. The total energy 00447 // from these is determined. 00448 // 00449 G4ThreeVector boost = pP.findBoostToCM(); 00450 G4Fragment *fragmentP = GetAbradedNucleons (Dabr, AP, ZP, rP); 00451 G4int nSecP = theParticleChange.GetNumberOfSecondaries(); 00452 G4int i = 0; 00453 for (i=0; i<nSecP; i++) 00454 { 00455 TotalEPost += theParticleChange.GetSecondary(i)-> 00456 GetParticle()->GetTotalEnergy(); 00457 } 00458 // 00459 // 00460 // Determine the number of spectators in the interaction region for the 00461 // projectile. 00462 // 00463 G4int DspcP = (G4int) (AP*F) - Dabr; 00464 if (DspcP <= 0) DspcP = 0; 00465 else if (DspcP > AP-Dabr) DspcP = ((G4int) AP) - Dabr; 00466 // 00467 // 00468 // Determine excitation energy associated with excess surface area of the 00469 // projectile (EsP) and the excitation due to scattering of nucleons which are 00470 // retained within the projectile (ExP). Add the total energy from the excited 00471 // nucleus to the total energy of the secondaries. 00472 // 00473 G4bool excitationAbsorbedByProjectile = false; 00474 if (fragmentP != NULL) 00475 { 00476 G4double EsP = theAbrasionGeometry->GetExcitationEnergyOfProjectile(); 00477 G4double ExP = 0.0; 00478 if (Dabr < AT) 00479 excitationAbsorbedByProjectile = G4UniformRand() < 0.5; 00480 if (excitationAbsorbedByProjectile) 00481 ExP = GetNucleonInducedExcitation(rP, rT, r); 00482 G4double xP = EsP + ExP; 00483 if (xP > B*(AP-Dabr)) xP = B*(AP-Dabr); 00484 G4LorentzVector lorentzVector = fragmentP->GetMomentum(); 00485 lorentzVector.setE(lorentzVector.e()+xP); 00486 fragmentP->SetMomentum(lorentzVector); 00487 TotalEPost += lorentzVector.e(); 00488 } 00489 G4double EMassP = TotalEPost; 00490 // 00491 // 00492 // Determine the abraded secondary nucleons from the target. Note that it's 00493 // assumed that the same number of nucleons are abraded from the target as for 00494 // the projectile, and obviously no boost is applied to the products. *fragmentT 00495 // is a pointer to the prefragment from the target and nSec is the total number 00496 // of nucleons in theParticleChange which have been abraded. The total energy 00497 // from these is determined. 00498 // 00499 G4Fragment *fragmentT = GetAbradedNucleons (Dabr, AT, ZT, rT); 00500 G4int nSec = theParticleChange.GetNumberOfSecondaries(); 00501 for (i=nSecP; i<nSec; i++) 00502 { 00503 TotalEPost += theParticleChange.GetSecondary(i)-> 00504 GetParticle()->GetTotalEnergy(); 00505 } 00506 // 00507 // 00508 // Determine the number of spectators in the interaction region for the 00509 // target. 00510 // 00511 G4int DspcT = (G4int) (AT*F) - Dabr; 00512 if (DspcT <= 0) DspcT = 0; 00513 else if (DspcT > AP-Dabr) DspcT = ((G4int) AT) - Dabr; 00514 // 00515 // 00516 // Determine excitation energy associated with excess surface area of the 00517 // target (EsT) and the excitation due to scattering of nucleons which are 00518 // retained within the target (ExT). Add the total energy from the excited 00519 // nucleus to the total energy of the secondaries. 00520 // 00521 if (fragmentT != NULL) 00522 { 00523 G4double EsT = theAbrasionGeometry->GetExcitationEnergyOfTarget(); 00524 G4double ExT = 0.0; 00525 if (!excitationAbsorbedByProjectile) 00526 ExT = GetNucleonInducedExcitation(rT, rP, r); 00527 G4double xT = EsT + ExT; 00528 if (xT > B*(AT-Dabr)) xT = B*(AT-Dabr); 00529 G4LorentzVector lorentzVector = fragmentT->GetMomentum(); 00530 lorentzVector.setE(lorentzVector.e()+xT); 00531 fragmentT->SetMomentum(lorentzVector); 00532 TotalEPost += lorentzVector.e(); 00533 } 00534 // 00535 // 00536 // Now determine the difference between the pre and post interaction 00537 // energy - this will be used to determine the Lorentz boost if conservation 00538 // of energy is to be imposed/attempted. 00539 // 00540 G4double deltaE = TotalEPre - TotalEPost; 00541 if (deltaE > 0.0 && conserveEnergy) 00542 { 00543 G4double beta = std::sqrt(1.0 - EMassP*EMassP/std::pow(deltaE+EMassP,2.0)); 00544 boost = boost / boost.mag() * beta; 00545 } 00546 // 00547 // 00548 // Now boost the secondaries from the projectile. 00549 // 00550 G4ThreeVector pBalance = pP.vect(); 00551 for (i=0; i<nSecP; i++) 00552 { 00553 G4DynamicParticle *dynamicP = theParticleChange.GetSecondary(i)-> 00554 GetParticle(); 00555 G4LorentzVector lorentzVector = dynamicP->Get4Momentum(); 00556 lorentzVector.boost(-boost); 00557 dynamicP->Set4Momentum(lorentzVector); 00558 pBalance -= lorentzVector.vect(); 00559 } 00560 // 00561 // 00562 // Set the boost for the projectile prefragment. This is now based on the 00563 // conservation of momentum. However, if the user selected momentum of the 00564 // prefragment is not to be conserved this simply boosted to the velocity of the 00565 // original projectile times the ratio of the unexcited to the excited mass 00566 // of the prefragment (the excitation increases the effective mass of the 00567 // prefragment, and therefore modifying the boost is an attempt to prevent 00568 // the momentum of the prefragment being excessive). 00569 // 00570 if (fragmentP != NULL) 00571 { 00572 G4LorentzVector lorentzVector = fragmentP->GetMomentum(); 00573 G4double fragmentM = lorentzVector.m(); 00574 if (conserveMomentum) 00575 fragmentP->SetMomentum 00576 (G4LorentzVector(pBalance,std::sqrt(pBalance.mag2()+fragmentM*fragmentM+1.0*eV*eV))); 00577 else 00578 { 00579 G4double fragmentGroundStateM = fragmentP->GetGroundStateMass(); 00580 fragmentP->SetMomentum(lorentzVector.boost(-boost * fragmentGroundStateM/fragmentM)); 00581 } 00582 } 00583 // 00584 // 00585 // Output information to user if verbose information requested. 00586 // 00587 if (verboseLevel >= 2) 00588 { 00589 G4cout <<G4endl; 00590 G4cout <<"-----------------------------------" <<G4endl; 00591 G4cout <<"Secondary nucleons from projectile:" <<G4endl; 00592 G4cout <<"-----------------------------------" <<G4endl; 00593 G4cout.precision(7); 00594 for (i=0; i<nSecP; i++) 00595 { 00596 G4cout <<"Particle # " <<i <<G4endl; 00597 theParticleChange.GetSecondary(i)->GetParticle()->DumpInfo(); 00598 G4DynamicParticle *dyn = theParticleChange.GetSecondary(i)->GetParticle(); 00599 G4cout <<"New nucleon (P) " <<dyn->GetDefinition()->GetParticleName() 00600 <<" : " <<dyn->Get4Momentum() 00601 <<G4endl; 00602 } 00603 G4cout <<"---------------------------" <<G4endl; 00604 G4cout <<"The projectile prefragment:" <<G4endl; 00605 G4cout <<"---------------------------" <<G4endl; 00606 if (fragmentP != NULL) 00607 G4cout <<*fragmentP <<G4endl; 00608 else 00609 G4cout <<"(No residual prefragment)" <<G4endl; 00610 G4cout <<G4endl; 00611 G4cout <<"-------------------------------" <<G4endl; 00612 G4cout <<"Secondary nucleons from target:" <<G4endl; 00613 G4cout <<"-------------------------------" <<G4endl; 00614 G4cout.precision(7); 00615 for (i=nSecP; i<nSec; i++) 00616 { 00617 G4cout <<"Particle # " <<i <<G4endl; 00618 theParticleChange.GetSecondary(i)->GetParticle()->DumpInfo(); 00619 G4DynamicParticle *dyn = theParticleChange.GetSecondary(i)->GetParticle(); 00620 G4cout <<"New nucleon (T) " <<dyn->GetDefinition()->GetParticleName() 00621 <<" : " <<dyn->Get4Momentum() 00622 <<G4endl; 00623 } 00624 G4cout <<"-----------------------" <<G4endl; 00625 G4cout <<"The target prefragment:" <<G4endl; 00626 G4cout <<"-----------------------" <<G4endl; 00627 if (fragmentT != NULL) 00628 G4cout <<*fragmentT <<G4endl; 00629 else 00630 G4cout <<"(No residual prefragment)" <<G4endl; 00631 } 00632 // 00633 // 00634 // Now we can decay the nuclear fragments if present. The secondaries are 00635 // collected and boosted as well. This is performed first for the projectile... 00636 // 00637 if (fragmentP !=NULL) 00638 { 00639 G4ReactionProductVector *products = NULL; 00640 if (fragmentP->GetZ() != fragmentP->GetA()) 00641 products = theExcitationHandler->BreakItUp(*fragmentP); 00642 else 00643 products = theExcitationHandlerx->BreakItUp(*fragmentP); 00644 delete fragmentP; 00645 fragmentP = NULL; 00646 00647 G4ReactionProductVector::iterator iter; 00648 for (iter = products->begin(); iter != products->end(); ++iter) 00649 { 00650 G4DynamicParticle *secondary = 00651 new G4DynamicParticle((*iter)->GetDefinition(), 00652 (*iter)->GetTotalEnergy(), (*iter)->GetMomentum()); 00653 theParticleChange.AddSecondary (secondary); // Added MHM 20050118 00654 G4String particleName = (*iter)->GetDefinition()->GetParticleName(); 00655 delete (*iter); // get rid of leftover particle def! // Added MHM 20050118 00656 if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size()) 00657 { 00658 G4cout <<"------------------------" <<G4endl; 00659 G4cout <<"The projectile fragment:" <<G4endl; 00660 G4cout <<"------------------------" <<G4endl; 00661 G4cout <<" fragmentP = " <<particleName 00662 <<" Energy = " <<secondary->GetKineticEnergy() 00663 <<G4endl; 00664 } 00665 } 00666 delete products; // Added MHM 20050118 00667 } 00668 // 00669 // 00670 // Now decay the target nucleus - no boost is applied since in this 00671 // approximation it is assumed that there is negligible momentum transfer from 00672 // the projectile. 00673 // 00674 if (fragmentT != NULL) 00675 { 00676 G4ReactionProductVector *products = NULL; 00677 if (fragmentT->GetZ() != fragmentT->GetA()) 00678 products = theExcitationHandler->BreakItUp(*fragmentT); 00679 else 00680 products = theExcitationHandlerx->BreakItUp(*fragmentT); 00681 delete fragmentT; 00682 fragmentT = NULL; 00683 00684 G4ReactionProductVector::iterator iter; 00685 for (iter = products->begin(); iter != products->end(); ++iter) 00686 { 00687 G4DynamicParticle *secondary = 00688 new G4DynamicParticle((*iter)->GetDefinition(), 00689 (*iter)->GetTotalEnergy(), (*iter)->GetMomentum()); 00690 theParticleChange.AddSecondary (secondary); 00691 G4String particleName = (*iter)->GetDefinition()->GetParticleName(); 00692 delete (*iter); // get rid of leftover particle def! // Added MHM 20050118 00693 if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size()) 00694 { 00695 G4cout <<"--------------------" <<G4endl; 00696 G4cout <<"The target fragment:" <<G4endl; 00697 G4cout <<"--------------------" <<G4endl; 00698 G4cout <<" fragmentT = " <<particleName 00699 <<" Energy = " <<secondary->GetKineticEnergy() 00700 <<G4endl; 00701 } 00702 } 00703 delete products; // Added MHM 20050118 00704 } 00705 00706 if (verboseLevel >= 2) 00707 G4cout <<"########################################" 00708 <<"########################################" 00709 <<G4endl; 00710 00711 delete theAbrasionGeometry; 00712 00713 return &theParticleChange; 00714 }
G4bool G4WilsonAbrasionModel::GetConserveMomentum | ( | ) | [inline] |
G4ExcitationHandler * G4WilsonAbrasionModel::GetExcitationHandler | ( | ) | [inline] |
G4bool G4WilsonAbrasionModel::GetUseAblation | ( | ) | [inline] |
void G4WilsonAbrasionModel::ModelDescription | ( | std::ostream & | ) | const [virtual] |
Reimplemented from G4HadronicInteraction.
Definition at line 167 of file G4WilsonAbrasionModel.cc.
00168 { 00169 outFile << "G4WilsonAbrasionModel is a macroscopic treatment of\n" 00170 << "nucleus-nucleus collisions using simple geometric arguments.\n" 00171 << "The smaller projectile nucleus gouges out a part of the larger\n" 00172 << "target nucleus, leaving a residual nucleus and a fireball\n" 00173 << "region where the projectile and target intersect. The fireball" 00174 << "is then treated as a highly excited nuclear fragment. This\n" 00175 << "model is based on the NUCFRG2 model and is valid for all\n" 00176 << "projectile energies between 70 MeV/n and 10.1 GeV/n. \n"; 00177 }
const G4WilsonAbrasionModel& G4WilsonAbrasionModel::operator= | ( | G4WilsonAbrasionModel & | right | ) |
void G4WilsonAbrasionModel::SetConserveMomentum | ( | G4bool | ) | [inline] |
void G4WilsonAbrasionModel::SetExcitationHandler | ( | G4ExcitationHandler * | ) | [inline] |
void G4WilsonAbrasionModel::SetUseAblation | ( | G4bool | ) |
Definition at line 875 of file G4WilsonAbrasionModel.cc.
References G4ExcitationHandler::SetEvaporation(), G4ExcitationHandler::SetFermiModel(), G4ExcitationHandler::SetMaxAandZForFermiBreakUp(), G4ExcitationHandler::SetMinEForMultiFrag(), G4ExcitationHandler::SetMultiFragmentation(), G4WilsonAblationModel::SetVerboseLevel(), and G4HadronicInteraction::verboseLevel.
00876 { 00877 if (useAblation != useAblation1) 00878 { 00879 useAblation = useAblation1; 00880 delete theExcitationHandler; 00881 delete theExcitationHandlerx; 00882 theExcitationHandler = new G4ExcitationHandler; 00883 theExcitationHandlerx = new G4ExcitationHandler; 00884 if (useAblation) 00885 { 00886 theAblation = new G4WilsonAblationModel; 00887 theAblation->SetVerboseLevel(verboseLevel); 00888 theExcitationHandler->SetEvaporation(theAblation); 00889 theExcitationHandlerx->SetEvaporation(theAblation); 00890 } 00891 else 00892 { 00893 theAblation = NULL; 00894 G4Evaporation * theEvaporation = new G4Evaporation; 00895 G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp; 00896 G4StatMF * theMF = new G4StatMF; 00897 theExcitationHandler->SetEvaporation(theEvaporation); 00898 theExcitationHandler->SetFermiModel(theFermiBreakUp); 00899 theExcitationHandler->SetMultiFragmentation(theMF); 00900 theExcitationHandler->SetMaxAandZForFermiBreakUp(12, 6); 00901 theExcitationHandler->SetMinEForMultiFrag(5.0*MeV); 00902 00903 theEvaporation = new G4Evaporation; 00904 theFermiBreakUp = new G4FermiBreakUp; 00905 theExcitationHandlerx->SetEvaporation(theEvaporation); 00906 theExcitationHandlerx->SetFermiModel(theFermiBreakUp); 00907 theExcitationHandlerx->SetMaxAandZForFermiBreakUp(12, 6); 00908 } 00909 } 00910 return; 00911 }
void G4WilsonAbrasionModel::SetVerboseLevel | ( | G4int | ) | [inline] |
Reimplemented from G4HadronicInteraction.
Definition at line 142 of file G4WilsonAbrasionModel.hh.
References G4WilsonAblationModel::SetVerboseLevel(), and G4HadronicInteraction::verboseLevel.
00143 { 00144 verboseLevel = verboseLevel1; 00145 if (useAblation) theAblation->SetVerboseLevel(verboseLevel); 00146 }