#include <G4NeutronHPInelasticBaseFS.hh>
Inheritance diagram for G4NeutronHPInelasticBaseFS:
Public Member Functions | |
G4NeutronHPInelasticBaseFS () | |
virtual | ~G4NeutronHPInelasticBaseFS () |
void | Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &bit) |
void | BaseApply (const G4HadProjectile &theTrack, G4ParticleDefinition **theDefs, G4int nDef) |
void | InitGammas (G4double AR, G4double ZR) |
virtual G4HadFinalState * | ApplyYourself (const G4HadProjectile &theTrack)=0 |
virtual G4NeutronHPFinalState * | New ()=0 |
virtual G4double | GetXsec (G4double anEnergy) |
virtual G4NeutronHPVector * | GetXsec () |
Protected Attributes | |
G4NeutronHPVector * | theXsection |
G4NeutronHPEnergyDistribution * | theEnergyDistribution |
G4NeutronHPAngular * | theAngularDistribution |
G4NeutronHPEnAngCorrelation * | theEnergyAngData |
G4NeutronHPPhotonDist * | theFinalStatePhotons |
G4double | theNuclearMassDifference |
G4NeutronHPDeExGammas | theGammas |
G4String | gammaPath |
Definition at line 42 of file G4NeutronHPInelasticBaseFS.hh.
G4NeutronHPInelasticBaseFS::G4NeutronHPInelasticBaseFS | ( | ) | [inline] |
Definition at line 46 of file G4NeutronHPInelasticBaseFS.hh.
References G4NeutronHPFinalState::hasXsec, theAngularDistribution, theEnergyAngData, theEnergyDistribution, theFinalStatePhotons, and theXsection.
00047 { 00048 hasXsec = true; 00049 theXsection = new G4NeutronHPVector; 00050 00051 theEnergyDistribution = 0; 00052 theFinalStatePhotons = 0; 00053 theEnergyAngData = 0; 00054 theAngularDistribution = 0; 00055 00056 }
virtual G4NeutronHPInelasticBaseFS::~G4NeutronHPInelasticBaseFS | ( | ) | [inline, virtual] |
Definition at line 57 of file G4NeutronHPInelasticBaseFS.hh.
References theAngularDistribution, theEnergyAngData, theEnergyDistribution, theFinalStatePhotons, and theXsection.
00058 { 00059 delete theXsection; 00060 if(theEnergyDistribution!=0) delete theEnergyDistribution; 00061 if(theFinalStatePhotons!=0) delete theFinalStatePhotons; 00062 if(theEnergyAngData!=0) delete theEnergyAngData; 00063 if(theAngularDistribution!=0) delete theAngularDistribution; 00064 }
virtual G4HadFinalState* G4NeutronHPInelasticBaseFS::ApplyYourself | ( | const G4HadProjectile & | theTrack | ) | [pure virtual] |
Reimplemented from G4NeutronHPFinalState.
Implemented in G4NeutronHP2AInelasticFS, G4NeutronHP2N2AInelasticFS, G4NeutronHP2NAInelasticFS, G4NeutronHP2NDInelasticFS, G4NeutronHP2NInelasticFS, G4NeutronHP2NPInelasticFS, G4NeutronHP2PInelasticFS, G4NeutronHP3AInelasticFS, G4NeutronHP3NAInelasticFS, G4NeutronHP3NInelasticFS, G4NeutronHP3NPInelasticFS, G4NeutronHP4NInelasticFS, G4NeutronHPD2AInelasticFS, G4NeutronHPDAInelasticFS, G4NeutronHPN2AInelasticFS, G4NeutronHPN2PInelasticFS, G4NeutronHPN3AInelasticFS, G4NeutronHPNAInelasticFS, G4NeutronHPND2AInelasticFS, G4NeutronHPNDInelasticFS, G4NeutronHPNHe3InelasticFS, G4NeutronHPNPAInelasticFS, G4NeutronHPNPInelasticFS, G4NeutronHPNT2AInelasticFS, G4NeutronHPNTInelasticFS, G4NeutronHPNXInelasticFS, G4NeutronHPPAInelasticFS, G4NeutronHPPDInelasticFS, G4NeutronHPPTInelasticFS, and G4NeutronHPT2AInelasticFS.
void G4NeutronHPInelasticBaseFS::BaseApply | ( | const G4HadProjectile & | theTrack, | |
G4ParticleDefinition ** | theDefs, | |||
G4int | nDef | |||
) |
Definition at line 168 of file G4NeutronHPInelasticBaseFS.cc.
References G4HadFinalState::AddSecondary(), G4NeutronHPFinalState::adjust_final_state(), G4Alpha::Alpha(), G4HadFinalState::Clear(), G4Deuteron::Deuteron(), G4cout, G4endl, G4UniformRand, G4HadProjectile::Get4Momentum(), G4ParticleDefinition::GetBaryonNumber(), G4Nucleus::GetBiasedThermalNucleus(), G4NucleiProperties::GetBindingEnergy(), G4NeutronHPDeExGammas::GetDecayGammas(), G4ReactionProduct::GetDefinition(), G4HadProjectile::GetDefinition(), G4ParticleTable::GetIon(), G4ReactionProduct::GetKineticEnergy(), G4HadProjectile::GetKineticEnergy(), G4NeutronHPDeExGammas::GetLevelEnergy(), G4ReactionProduct::GetMass(), G4HadProjectile::GetMaterial(), G4ReactionProduct::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4NeutronHPDeExGammas::GetNumberOfLevels(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4NeutronHPPhotonDist::GetPhotons(), G4NeutronHPAngular::GetTargetMass(), G4NeutronHPEnAngCorrelation::GetTargetMass(), G4Material::GetTemperature(), G4ReactionProduct::GetTotalEnergy(), G4NeutronHPEnAngCorrelation::GetTotalMeanEnergy(), G4NeutronHPFinalState::HasFSData(), G4He3::He3(), G4NeutronHPNBodyPhaseSpace::Init(), G4ReactionProduct::Lorentz(), G4Neutron::Neutron(), G4Proton::Proton(), G4NeutronHPEnergyDistribution::Sample(), G4NeutronHPEnAngCorrelation::Sample(), G4NeutronHPNBodyPhaseSpace::Sample(), G4NeutronHPAngular::SampleAndUpdate(), G4ReactionProduct::SetDefinition(), G4DynamicParticle::SetDefinition(), G4ReactionProduct::SetKineticEnergy(), G4DynamicParticle::SetMomentum(), G4ReactionProduct::SetMomentum(), G4NeutronHPEnAngCorrelation::SetNeutron(), G4NeutronHPAngular::SetNeutron(), G4VNeutronHPEnergyAngular::SetNeutron(), G4HadFinalState::SetStatusChange(), G4NeutronHPEnAngCorrelation::SetTarget(), G4NeutronHPAngular::SetTarget(), G4VNeutronHPEnergyAngular::SetTarget(), G4ReactionProduct::SetTotalEnergy(), stopAndKill, theAngularDistribution, G4NeutronHPFinalState::theBaseA, G4NeutronHPFinalState::theBaseZ, theEnergyAngData, theEnergyDistribution, theFinalStatePhotons, theGammas, theNuclearMassDifference, G4NeutronHPFinalState::theResult, and G4Triton::Triton().
Referenced by G4NeutronHPT2AInelasticFS::ApplyYourself(), G4NeutronHPPTInelasticFS::ApplyYourself(), G4NeutronHPPDInelasticFS::ApplyYourself(), G4NeutronHPPAInelasticFS::ApplyYourself(), G4NeutronHPNXInelasticFS::ApplyYourself(), G4NeutronHPNTInelasticFS::ApplyYourself(), G4NeutronHPNT2AInelasticFS::ApplyYourself(), G4NeutronHPNPInelasticFS::ApplyYourself(), G4NeutronHPNPAInelasticFS::ApplyYourself(), G4NeutronHPNHe3InelasticFS::ApplyYourself(), G4NeutronHPNDInelasticFS::ApplyYourself(), G4NeutronHPND2AInelasticFS::ApplyYourself(), G4NeutronHPNAInelasticFS::ApplyYourself(), G4NeutronHPN3AInelasticFS::ApplyYourself(), G4NeutronHPN2PInelasticFS::ApplyYourself(), G4NeutronHPN2AInelasticFS::ApplyYourself(), G4NeutronHPDAInelasticFS::ApplyYourself(), G4NeutronHPD2AInelasticFS::ApplyYourself(), G4NeutronHP4NInelasticFS::ApplyYourself(), G4NeutronHP3NPInelasticFS::ApplyYourself(), G4NeutronHP3NInelasticFS::ApplyYourself(), G4NeutronHP3NAInelasticFS::ApplyYourself(), G4NeutronHP3AInelasticFS::ApplyYourself(), G4NeutronHP2PInelasticFS::ApplyYourself(), G4NeutronHP2NPInelasticFS::ApplyYourself(), G4NeutronHP2NInelasticFS::ApplyYourself(), G4NeutronHP2NDInelasticFS::ApplyYourself(), G4NeutronHP2NAInelasticFS::ApplyYourself(), G4NeutronHP2N2AInelasticFS::ApplyYourself(), and G4NeutronHP2AInelasticFS::ApplyYourself().
00171 { 00172 00173 // prepare neutron 00174 theResult.Clear(); 00175 G4double eKinetic = theTrack.GetKineticEnergy(); 00176 const G4HadProjectile *incidentParticle = &theTrack; 00177 G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) ); 00178 theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() ); 00179 theNeutron.SetKineticEnergy( eKinetic ); 00180 00181 // prepare target 00182 G4double targetMass; 00183 G4double eps = 0.0001; 00184 targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / 00185 G4Neutron::Neutron()->GetPDGMass(); 00186 00187 if(theEnergyAngData!=0) 00188 { targetMass = theEnergyAngData->GetTargetMass(); } 00189 if(theAngularDistribution!=0) 00190 { targetMass = theAngularDistribution->GetTargetMass(); } 00191 00192 //080731a 00193 //110512 TKDB ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded. 00194 //if ( targetMass == 0 ) G4cout << "080731a It looks like something wrong value in G4NDL, please update the latest version. If you use the latest, then please report this problem to Geant4 Hyper news." << G4endl; 00195 if ( targetMass == 0 ) 00196 { 00197 //G4cout << "TKDB targetMass = 0; ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded. This could be a similar situation." << G4endl; 00198 targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / G4Neutron::Neutron()->GetPDGMass(); 00199 } 00200 00201 G4Nucleus aNucleus; 00202 G4ReactionProduct theTarget; 00203 G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum(); 00204 theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature()); 00205 00206 // prepare energy in target rest frame 00207 G4ReactionProduct boosted; 00208 boosted.Lorentz(theNeutron, theTarget); 00209 eKinetic = boosted.GetKineticEnergy(); 00210 G4double orgMomentum = boosted.GetMomentum().mag(); 00211 00212 // Take N-body phase-space distribution, if no other data present. 00213 if(!HasFSData()) // adding the residual is trivial here @@@ 00214 { 00215 G4NeutronHPNBodyPhaseSpace thePhaseSpaceDistribution; 00216 G4double aPhaseMass=0; 00217 G4int ii; 00218 for(ii=0; ii<nDef; ii++) 00219 { 00220 aPhaseMass+=theDefs[ii]->GetPDGMass(); 00221 } 00222 thePhaseSpaceDistribution.Init(aPhaseMass, nDef); 00223 thePhaseSpaceDistribution.SetNeutron(&theNeutron); 00224 thePhaseSpaceDistribution.SetTarget(&theTarget); 00225 for(ii=0; ii<nDef; ii++) 00226 { 00227 G4double massCode = 1000.*std::abs(theDefs[ii]->GetPDGCharge()); 00228 massCode += theDefs[ii]->GetBaryonNumber(); 00229 G4double dummy = 0; 00230 G4ReactionProduct * aSec = thePhaseSpaceDistribution.Sample(eKinetic, massCode, dummy); 00231 aSec->Lorentz(*aSec, -1.*theTarget); 00232 G4DynamicParticle * aPart = new G4DynamicParticle(); 00233 aPart->SetDefinition(aSec->GetDefinition()); 00234 aPart->SetMomentum(aSec->GetMomentum()); 00235 delete aSec; 00236 theResult.AddSecondary(aPart); 00237 } 00238 theResult.SetStatusChange(stopAndKill); 00239 00240 //TK120607 00241 //Final momentum check should be done before return 00242 G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , 0.0 ); 00243 G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) ); 00244 G4LorentzVector proj_4p_lab = theTrack.Get4Momentum(); 00245 G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab; 00246 adjust_final_state ( init_4p_lab ); 00247 00248 return; 00249 } 00250 00251 // set target and neutron in the relevant exit channel 00252 if(theAngularDistribution!=0) 00253 { 00254 theAngularDistribution->SetTarget(theTarget); 00255 theAngularDistribution->SetNeutron(theNeutron); 00256 } 00257 else if(theEnergyAngData!=0) 00258 { 00259 theEnergyAngData->SetTarget(theTarget); 00260 theEnergyAngData->SetNeutron(theNeutron); 00261 } 00262 00263 G4ReactionProductVector * tmpHadrons = 0; 00264 G4int ii, dummy; 00265 unsigned int i; 00266 if(theEnergyAngData != 0) 00267 { 00268 tmpHadrons = theEnergyAngData->Sample(eKinetic); 00269 } 00270 else if(theAngularDistribution!= 0) 00271 { 00272 G4bool * Done = new G4bool[nDef]; 00273 G4int i0; 00274 for(i0=0; i0<nDef; i0++) Done[i0] = false; 00275 if(tmpHadrons == 0) 00276 { 00277 tmpHadrons = new G4ReactionProductVector; 00278 } 00279 else 00280 { 00281 for(i=0; i<tmpHadrons->size(); i++) 00282 { 00283 for(ii=0; ii<nDef; ii++) 00284 if(!Done[ii] && tmpHadrons->operator[](i)->GetDefinition() == theDefs[ii]) 00285 Done[ii] = true; 00286 } 00287 } 00288 G4ReactionProduct * aHadron; 00289 G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))); 00290 G4ThreeVector bufferedDirection(0,0,0); 00291 for(i0=0; i0<nDef; i0++) 00292 { 00293 if(!Done[i0]) 00294 { 00295 aHadron = new G4ReactionProduct; 00296 if(theEnergyDistribution!=0) 00297 { 00298 aHadron->SetDefinition(theDefs[i0]); 00299 aHadron->SetKineticEnergy(theEnergyDistribution->Sample(eKinetic, dummy)); 00300 } 00301 else if(nDef == 1) 00302 { 00303 aHadron->SetDefinition(theDefs[i0]); 00304 aHadron->SetKineticEnergy(eKinetic); 00305 } 00306 else if(nDef == 2) 00307 { 00308 aHadron->SetDefinition(theDefs[i0]); 00309 aHadron->SetKineticEnergy(50*MeV); 00310 } 00311 else 00312 { 00313 throw G4HadronicException(__FILE__, __LINE__, "No energy distribution to sample from in InelasticBaseFS::BaseApply"); 00314 } 00315 theAngularDistribution->SampleAndUpdate(*aHadron); 00316 if(theEnergyDistribution==0 && nDef == 2) 00317 { 00318 if(i0==0) 00319 { 00320 G4double mass1 = theDefs[0]->GetPDGMass(); 00321 G4double mass2 = theDefs[1]->GetPDGMass(); 00322 G4double massn = G4Neutron::Neutron()->GetPDGMass(); 00323 G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge()); 00324 G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber(); 00325 G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1); 00326 G4double availableEnergy = eKinetic+massn+localMass-mass1-mass2-concreteMass; 00327 // available kinetic energy in CMS (non relativistic) 00328 G4double emin = availableEnergy+mass1+mass2 - std::sqrt((mass1+mass2)*(mass1+mass2)+orgMomentum*orgMomentum); 00329 G4double p1=std::sqrt(2.*mass2*emin); 00330 bufferedDirection = p1*aHadron->GetMomentum().unit(); 00331 if(getenv("HTOKEN")) // @@@@@ verify the nucleon counting... 00332 { 00333 G4cout << "HTOKEN "<<z1<<" "<<theBaseZ<<" "<<a1<<" "<<theBaseA<<" "<<availableEnergy<<" " 00334 << emin<<G4endl; 00335 } 00336 } 00337 else 00338 { 00339 bufferedDirection = -bufferedDirection; 00340 } 00341 // boost from cms to lab 00342 if(getenv("HTOKEN")) 00343 { 00344 G4cout << " HTOKEN "<<bufferedDirection.mag2()<<G4endl; 00345 } 00346 aHadron->SetTotalEnergy( std::sqrt(aHadron->GetMass()*aHadron->GetMass() 00347 +bufferedDirection.mag2()) ); 00348 aHadron->SetMomentum(bufferedDirection); 00349 aHadron->Lorentz(*aHadron, -1.*(theTarget+theNeutron)); 00350 if(getenv("HTOKEN")) 00351 { 00352 G4cout << " HTOKEN "<<aHadron->GetTotalEnergy()<<" "<<aHadron->GetMomentum()<<G4endl; 00353 } 00354 } 00355 tmpHadrons->push_back(aHadron); 00356 } 00357 } 00358 delete [] Done; 00359 } 00360 else 00361 { 00362 throw G4HadronicException(__FILE__, __LINE__, "No data to create the neutrons in NInelasticFS"); 00363 } 00364 00365 G4ReactionProductVector * thePhotons = 0; 00366 if(theFinalStatePhotons!=0) 00367 { 00368 // the photon distributions are in the Nucleus rest frame. 00369 G4ReactionProduct boosted_tmp; 00370 boosted_tmp.Lorentz(theNeutron, theTarget); 00371 G4double anEnergy = boosted_tmp.GetKineticEnergy(); 00372 thePhotons = theFinalStatePhotons->GetPhotons(anEnergy); 00373 if(thePhotons!=0) 00374 { 00375 for(i=0; i<thePhotons->size(); i++) 00376 { 00377 // back to lab 00378 thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), -1.*theTarget); 00379 } 00380 } 00381 } 00382 else if(theEnergyAngData!=0) 00383 { 00384 00385 G4double theGammaEnergy = theEnergyAngData->GetTotalMeanEnergy(); 00386 G4double anEnergy = boosted.GetKineticEnergy(); 00387 theGammaEnergy = anEnergy-theGammaEnergy; 00388 theGammaEnergy += theNuclearMassDifference; 00389 G4double eBindProducts = 0; 00390 G4double eBindN = 0; 00391 G4double eBindP = 0; 00392 G4double eBindD = G4NucleiProperties::GetBindingEnergy(2,1); 00393 G4double eBindT = G4NucleiProperties::GetBindingEnergy(3,1); 00394 G4double eBindHe3 = G4NucleiProperties::GetBindingEnergy(3,2); 00395 G4double eBindA = G4NucleiProperties::GetBindingEnergy(4,2); 00396 G4int ia=0; 00397 for(i=0; i<tmpHadrons->size(); i++) 00398 { 00399 if(tmpHadrons->operator[](i)->GetDefinition() == G4Neutron::Neutron()) 00400 { 00401 eBindProducts+=eBindN; 00402 } 00403 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Proton::Proton()) 00404 { 00405 eBindProducts+=eBindP; 00406 } 00407 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Deuteron::Deuteron()) 00408 { 00409 eBindProducts+=eBindD; 00410 } 00411 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Triton::Triton()) 00412 { 00413 eBindProducts+=eBindT; 00414 } 00415 else if(tmpHadrons->operator[](i)->GetDefinition() == G4He3::He3()) 00416 { 00417 eBindProducts+=eBindHe3; 00418 } 00419 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Alpha::Alpha()) 00420 { 00421 eBindProducts+=eBindA; 00422 ia++; 00423 } 00424 } 00425 00426 theGammaEnergy += eBindProducts; 00427 00428 //101111 00429 //Special treatment for Be9 + n -> 2n + Be8 -> 2n + a + a 00430 if ( (G4int)(theBaseZ+eps) == 4 && (G4int)(theBaseA+eps) == 9 ) 00431 { 00432 // This only valid for G4NDL3.13,,, 00433 if ( std::abs( theNuclearMassDifference - 00434 ( G4NucleiProperties::GetBindingEnergy( 8 , 4 ) - 00435 G4NucleiProperties::GetBindingEnergy( 9 , 4 ) ) ) < 1*keV 00436 && ia == 2 ) 00437 { 00438 theGammaEnergy -= (2*eBindA); 00439 } 00440 } 00441 00442 G4ReactionProductVector * theOtherPhotons = 0; 00443 G4int iLevel; 00444 while(theGammaEnergy>=theGammas.GetLevelEnergy(0)) 00445 { 00446 for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--) 00447 { 00448 if(theGammas.GetLevelEnergy(iLevel)<theGammaEnergy) break; 00449 } 00450 if(iLevel==0||iLevel==theGammas.GetNumberOfLevels()-1) 00451 { 00452 theOtherPhotons = theGammas.GetDecayGammas(iLevel); 00453 } 00454 else 00455 { 00456 G4double random = G4UniformRand(); 00457 G4double eLow = theGammas.GetLevelEnergy(iLevel); 00458 G4double eHigh = theGammas.GetLevelEnergy(iLevel+1); 00459 if(random > (eHigh-eLow)/(theGammaEnergy-eLow)) iLevel++; 00460 theOtherPhotons = theGammas.GetDecayGammas(iLevel); 00461 } 00462 if(thePhotons==0) thePhotons = new G4ReactionProductVector; 00463 if(theOtherPhotons != 0) 00464 { 00465 for(unsigned int iii=0; iii<theOtherPhotons->size(); iii++) 00466 { 00467 thePhotons->push_back(theOtherPhotons->operator[](iii)); 00468 } 00469 delete theOtherPhotons; 00470 } 00471 theGammaEnergy -= theGammas.GetLevelEnergy(iLevel); 00472 if(iLevel == -1) break; 00473 } 00474 } 00475 00476 // fill the result 00477 unsigned int nSecondaries = tmpHadrons->size(); 00478 unsigned int nPhotons = 0; 00479 if(thePhotons!=0) { nPhotons = thePhotons->size(); } 00480 nSecondaries += nPhotons; 00481 G4DynamicParticle * theSec; 00482 00483 for(i=0; i<nSecondaries-nPhotons; i++) 00484 { 00485 theSec = new G4DynamicParticle; 00486 theSec->SetDefinition(tmpHadrons->operator[](i)->GetDefinition()); 00487 theSec->SetMomentum(tmpHadrons->operator[](i)->GetMomentum()); 00488 theResult.AddSecondary(theSec); 00489 delete tmpHadrons->operator[](i); 00490 } 00491 if(thePhotons != 0) 00492 { 00493 for(i=0; i<nPhotons; i++) 00494 { 00495 theSec = new G4DynamicParticle; 00496 theSec->SetDefinition(thePhotons->operator[](i)->GetDefinition()); 00497 theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum()); 00498 theResult.AddSecondary(theSec); 00499 delete thePhotons->operator[](i); 00500 } 00501 } 00502 00503 // some garbage collection 00504 delete thePhotons; 00505 delete tmpHadrons; 00506 00507 //080721 00508 G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , 0.0 ); 00509 G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) ); 00510 G4LorentzVector proj_4p_lab = theTrack.Get4Momentum(); 00511 G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab; 00512 adjust_final_state ( init_4p_lab ); 00513 00514 // clean up the primary neutron 00515 theResult.SetStatusChange(stopAndKill); 00516 }
virtual G4NeutronHPVector* G4NeutronHPInelasticBaseFS::GetXsec | ( | ) | [inline, virtual] |
Reimplemented from G4NeutronHPFinalState.
Definition at line 76 of file G4NeutronHPInelasticBaseFS.hh.
References theXsection.
00076 { return theXsection; }
Reimplemented from G4NeutronHPFinalState.
Definition at line 72 of file G4NeutronHPInelasticBaseFS.hh.
References G4NeutronHPVector::GetY(), and theXsection.
00073 { 00074 return std::max(0., theXsection->GetY(anEnergy)); 00075 }
void G4NeutronHPInelasticBaseFS::Init | ( | G4double | A, | |
G4double | Z, | |||
G4int | M, | |||
G4String & | dirName, | |||
G4String & | bit | |||
) | [virtual] |
Implements G4NeutronHPFinalState.
Reimplemented in G4NeutronHP2AInelasticFS, G4NeutronHP2N2AInelasticFS, G4NeutronHP2NAInelasticFS, G4NeutronHP2NDInelasticFS, G4NeutronHP2NInelasticFS, G4NeutronHP2NPInelasticFS, G4NeutronHP2PInelasticFS, G4NeutronHP3AInelasticFS, G4NeutronHP3NAInelasticFS, G4NeutronHP3NInelasticFS, G4NeutronHP3NPInelasticFS, G4NeutronHP4NInelasticFS, G4NeutronHPD2AInelasticFS, G4NeutronHPDAInelasticFS, G4NeutronHPN2AInelasticFS, G4NeutronHPN2PInelasticFS, G4NeutronHPN3AInelasticFS, G4NeutronHPNAInelasticFS, G4NeutronHPND2AInelasticFS, G4NeutronHPNDInelasticFS, G4NeutronHPNHe3InelasticFS, G4NeutronHPNPAInelasticFS, G4NeutronHPNPInelasticFS, G4NeutronHPNT2AInelasticFS, G4NeutronHPNTInelasticFS, G4NeutronHPNXInelasticFS, G4NeutronHPPAInelasticFS, G4NeutronHPPDInelasticFS, G4NeutronHPPTInelasticFS, and G4NeutronHPT2AInelasticFS.
Definition at line 71 of file G4NeutronHPInelasticBaseFS.cc.
References G4cout, G4endl, gammaPath, G4NeutronHPDataUsed::GetName(), G4NeutronHPNames::GetName(), G4NeutronHPFinalState::hasAnyData, G4NeutronHPFinalState::hasFSData, G4NeutronHPFinalState::hasXsec, G4NeutronHPEnAngCorrelation::Init(), G4NeutronHPEnergyDistribution::Init(), G4NeutronHPAngular::Init(), G4NeutronHPVector::Init(), G4NeutronHPPhotonDist::InitAngular(), G4NeutronHPPhotonDist::InitEnergies(), G4NeutronHPPhotonDist::InitMean(), G4NeutronHPPhotonDist::InitPartials(), INT_MAX, G4NeutronHPFinalState::SetAZMs(), theAngularDistribution, theEnergyAngData, theEnergyDistribution, theFinalStatePhotons, G4NeutronHPFinalState::theNames, G4NeutronHPFinalState::theNDLDataA, G4NeutronHPFinalState::theNDLDataZ, and theXsection.
Referenced by G4NeutronHPT2AInelasticFS::Init(), G4NeutronHPPTInelasticFS::Init(), G4NeutronHPPDInelasticFS::Init(), G4NeutronHPPAInelasticFS::Init(), G4NeutronHPNXInelasticFS::Init(), G4NeutronHPNTInelasticFS::Init(), G4NeutronHPNT2AInelasticFS::Init(), G4NeutronHPNPInelasticFS::Init(), G4NeutronHPNPAInelasticFS::Init(), G4NeutronHPNHe3InelasticFS::Init(), G4NeutronHPNDInelasticFS::Init(), G4NeutronHPND2AInelasticFS::Init(), G4NeutronHPNAInelasticFS::Init(), G4NeutronHPN3AInelasticFS::Init(), G4NeutronHPN2PInelasticFS::Init(), G4NeutronHPN2AInelasticFS::Init(), G4NeutronHPDAInelasticFS::Init(), G4NeutronHPD2AInelasticFS::Init(), G4NeutronHP4NInelasticFS::Init(), G4NeutronHP3NPInelasticFS::Init(), G4NeutronHP3NInelasticFS::Init(), G4NeutronHP3NAInelasticFS::Init(), G4NeutronHP3AInelasticFS::Init(), G4NeutronHP2PInelasticFS::Init(), G4NeutronHP2NPInelasticFS::Init(), G4NeutronHP2NInelasticFS::Init(), G4NeutronHP2NDInelasticFS::Init(), G4NeutronHP2NAInelasticFS::Init(), G4NeutronHP2N2AInelasticFS::Init(), and G4NeutronHP2AInelasticFS::Init().
00072 { 00073 gammaPath = "/Inelastic/Gammas/"; 00074 if(!getenv("G4NEUTRONHPDATA")) 00075 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files."); 00076 G4String tBase = getenv("G4NEUTRONHPDATA"); 00077 gammaPath = tBase+gammaPath; 00078 G4String tString = dirName; 00079 G4bool dbool; 00080 G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M,tString, bit, dbool); 00081 G4String filename = aFile.GetName(); 00082 SetAZMs( A, Z, M, aFile); 00083 //theBaseA = aFile.GetA(); 00084 //theBaseZ = aFile.GetZ(); 00085 // theNDLDataA = (int)aFile.GetA(); 00086 // theNDLDataZ = aFile.GetZ(); 00087 //if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001))) 00088 if ( !dbool || ( Z<2.5 && ( std::abs(theNDLDataZ - Z)>0.0001 || std::abs(theNDLDataA - A)>0.0001)) ) 00089 { 00090 if(getenv("NeutronHPNamesLogging")) G4cout << "Skipped = "<< filename <<" "<<A<<" "<<Z<<G4endl; 00091 hasAnyData = false; 00092 hasFSData = false; 00093 hasXsec = false; 00094 return; 00095 } 00096 //theBaseA = A; 00097 //theBaseZ = G4int(Z+.5); 00098 std::ifstream theData(filename, std::ios::in); 00099 if(!(theData)) 00100 { 00101 hasAnyData = false; 00102 hasFSData = false; 00103 hasXsec = false; 00104 theData.close(); 00105 return; // no data for exactly this isotope and FS 00106 } 00107 // here we go 00108 G4int infoType, dataType, dummy=INT_MAX; 00109 hasFSData = false; 00110 while (theData >> infoType) 00111 { 00112 theData >> dataType; 00113 if(dummy==INT_MAX) theData >> dummy >> dummy; 00114 if(dataType==3) 00115 { 00116 G4int total; 00117 theData >> total; 00118 theXsection->Init(theData, total, eV); 00119 } 00120 else if(dataType==4) 00121 { 00122 theAngularDistribution = new G4NeutronHPAngular; 00123 theAngularDistribution->Init(theData); 00124 hasFSData = true; 00125 } 00126 else if(dataType==5) 00127 { 00128 theEnergyDistribution = new G4NeutronHPEnergyDistribution; 00129 theEnergyDistribution->Init(theData); 00130 hasFSData = true; 00131 } 00132 else if(dataType==6) 00133 { 00134 theEnergyAngData = new G4NeutronHPEnAngCorrelation; 00135 theEnergyAngData->Init(theData); 00136 hasFSData = true; 00137 } 00138 else if(dataType==12) 00139 { 00140 theFinalStatePhotons = new G4NeutronHPPhotonDist; 00141 theFinalStatePhotons->InitMean(theData); 00142 hasFSData = true; 00143 } 00144 else if(dataType==13) 00145 { 00146 theFinalStatePhotons = new G4NeutronHPPhotonDist; 00147 theFinalStatePhotons->InitPartials(theData); 00148 hasFSData = true; 00149 } 00150 else if(dataType==14) 00151 { 00152 theFinalStatePhotons->InitAngular(theData); 00153 hasFSData = true; 00154 } 00155 else if(dataType==15) 00156 { 00157 theFinalStatePhotons->InitEnergies(theData); 00158 hasFSData = true; 00159 } 00160 else 00161 { 00162 throw G4HadronicException(__FILE__, __LINE__, "Data-type unknown to G4NeutronHPInelasticBaseFS"); 00163 } 00164 } 00165 theData.close(); 00166 }
Definition at line 46 of file G4NeutronHPInelasticBaseFS.cc.
References gammaPath, G4NucleiProperties::GetBindingEnergy(), G4NeutronHPDeExGammas::Init(), G4NeutronHPFinalState::theBaseA, G4NeutronHPFinalState::theBaseZ, theGammas, and theNuclearMassDifference.
Referenced by G4NeutronHPT2AInelasticFS::Init(), G4NeutronHPPTInelasticFS::Init(), G4NeutronHPPDInelasticFS::Init(), G4NeutronHPPAInelasticFS::Init(), G4NeutronHPNXInelasticFS::Init(), G4NeutronHPNTInelasticFS::Init(), G4NeutronHPNT2AInelasticFS::Init(), G4NeutronHPNPInelasticFS::Init(), G4NeutronHPNPAInelasticFS::Init(), G4NeutronHPNHe3InelasticFS::Init(), G4NeutronHPNDInelasticFS::Init(), G4NeutronHPND2AInelasticFS::Init(), G4NeutronHPNAInelasticFS::Init(), G4NeutronHPN3AInelasticFS::Init(), G4NeutronHPN2PInelasticFS::Init(), G4NeutronHPN2AInelasticFS::Init(), G4NeutronHPDAInelasticFS::Init(), G4NeutronHPD2AInelasticFS::Init(), G4NeutronHP4NInelasticFS::Init(), G4NeutronHP3NPInelasticFS::Init(), G4NeutronHP3NInelasticFS::Init(), G4NeutronHP3NAInelasticFS::Init(), G4NeutronHP3AInelasticFS::Init(), G4NeutronHP2PInelasticFS::Init(), G4NeutronHP2NPInelasticFS::Init(), G4NeutronHP2NInelasticFS::Init(), G4NeutronHP2NDInelasticFS::Init(), G4NeutronHP2NAInelasticFS::Init(), G4NeutronHP2N2AInelasticFS::Init(), and G4NeutronHP2AInelasticFS::Init().
00047 { 00048 // char the[100] = {""}; 00049 // std::ostrstream ost(the, 100, std::ios::out); 00050 // ost <<gammaPath<<"z"<<ZR<<".a"<<AR; 00051 // G4String * aName = new G4String(the); 00052 // std::ifstream from(*aName, std::ios::in); 00053 00054 std::ostringstream ost; 00055 ost <<gammaPath<<"z"<<ZR<<".a"<<AR; 00056 G4String aName = ost.str(); 00057 std::ifstream from(aName, std::ios::in); 00058 00059 if(!from) return; // no data found for this isotope 00060 // std::ifstream theGammaData(*aName, std::ios::in); 00061 std::ifstream theGammaData(aName, std::ios::in); 00062 00063 G4double eps = 0.001; 00064 theNuclearMassDifference = 00065 G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) - 00066 G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)); 00067 theGammas.Init(theGammaData); 00068 // delete aName; 00069 }
virtual G4NeutronHPFinalState* G4NeutronHPInelasticBaseFS::New | ( | ) | [pure virtual] |
Implements G4NeutronHPFinalState.
Implemented in G4NeutronHP2AInelasticFS, G4NeutronHP2N2AInelasticFS, G4NeutronHP2NAInelasticFS, G4NeutronHP2NDInelasticFS, G4NeutronHP2NInelasticFS, G4NeutronHP2NPInelasticFS, G4NeutronHP2PInelasticFS, G4NeutronHP3AInelasticFS, G4NeutronHP3NAInelasticFS, G4NeutronHP3NInelasticFS, G4NeutronHP3NPInelasticFS, G4NeutronHP4NInelasticFS, G4NeutronHPD2AInelasticFS, G4NeutronHPDAInelasticFS, G4NeutronHPN2AInelasticFS, G4NeutronHPN2PInelasticFS, G4NeutronHPN3AInelasticFS, G4NeutronHPNAInelasticFS, G4NeutronHPND2AInelasticFS, G4NeutronHPNDInelasticFS, G4NeutronHPNHe3InelasticFS, G4NeutronHPNPAInelasticFS, G4NeutronHPNPInelasticFS, G4NeutronHPNT2AInelasticFS, G4NeutronHPNTInelasticFS, G4NeutronHPNXInelasticFS, G4NeutronHPPAInelasticFS, G4NeutronHPPDInelasticFS, G4NeutronHPPTInelasticFS, and G4NeutronHPT2AInelasticFS.
G4String G4NeutronHPInelasticBaseFS::gammaPath [protected] |
Definition at line 88 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by Init(), and InitGammas().
Definition at line 82 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by BaseApply(), G4NeutronHPInelasticBaseFS(), Init(), and ~G4NeutronHPInelasticBaseFS().
Definition at line 83 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by BaseApply(), G4NeutronHPInelasticBaseFS(), Init(), and ~G4NeutronHPInelasticBaseFS().
Definition at line 81 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by BaseApply(), G4NeutronHPInelasticBaseFS(), Init(), and ~G4NeutronHPInelasticBaseFS().
Definition at line 85 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by BaseApply(), G4NeutronHPInelasticBaseFS(), Init(), and ~G4NeutronHPInelasticBaseFS().
Definition at line 87 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by BaseApply(), and InitGammas().
Definition at line 86 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by BaseApply(), and InitGammas().
G4NeutronHPVector* G4NeutronHPInelasticBaseFS::theXsection [protected] |
Definition at line 80 of file G4NeutronHPInelasticBaseFS.hh.
Referenced by G4NeutronHPInelasticBaseFS(), GetXsec(), Init(), and ~G4NeutronHPInelasticBaseFS().