Ignore:
Timestamp:
Apr 6, 2009, 12:30:29 PM (15 years ago)
Author:
garnier
Message:

update processes

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPInelasticBaseFS.cc

    r819 r962  
    2828// A prototype of the low energy neutron transport model.
    2929//
     30// 080801 Give a warning message for irregular mass value in data file by T. Koi
     31//        Introduce theNDLDataA,Z which has A and Z of NDL data by T. Koi
     32// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
     33//
    3034#include "G4NeutronHPInelasticBaseFS.hh"
    3135#include "G4Nucleus.hh"
    32 #include "G4NucleiPropertiesTable.hh"
     36#include "G4NucleiProperties.hh"
    3337#include "G4He3.hh"
    3438#include "G4Alpha.hh"
    3539#include "G4Electron.hh"
    3640#include "G4NeutronHPDataUsed.hh"
     41
     42#include "G4ParticleTable.hh"
    3743
    3844void G4NeutronHPInelasticBaseFS::InitGammas(G4double AR, G4double ZR)
     
    5561   G4double eps = 0.001;
    5662   theNuclearMassDifference =
    57        G4NucleiPropertiesTable::GetBindingEnergy(static_cast<G4int>(ZR+eps),static_cast<G4int>(AR+eps)) -
    58        G4NucleiPropertiesTable::GetBindingEnergy(static_cast<G4int>(theBaseZ+eps), static_cast<G4int>(theBaseA+eps));
     63       G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) -
     64       G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps));
    5965   theGammas.Init(theGammaData);
    6066   //   delete aName;
     
    7480  theBaseA = aFile.GetA();
    7581  theBaseZ = aFile.GetZ();
     82   theNDLDataA = (int)aFile.GetA();
     83   theNDLDataZ = aFile.GetZ();
    7684  if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
    7785  {
     
    170178  G4double targetMass;
    171179  G4double eps = 0.0001;
    172   targetMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theBaseZ+eps), static_cast<G4int>(theBaseA+eps))) /
     180  targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
    173181               G4Neutron::Neutron()->GetPDGMass();
    174182  if(theEnergyAngData!=0)
     
    176184  if(theAngularDistribution!=0)
    177185     { targetMass = theAngularDistribution->GetTargetMass(); }
     186//080731a
     187if ( 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;
    178188  G4Nucleus aNucleus;
    179189  G4ReactionProduct theTarget;
     
    255265    }
    256266    G4ReactionProduct * aHadron;
    257     G4double localMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theBaseZ+eps), static_cast<G4int>(theBaseA+eps)));
     267    G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)));
    258268    G4ThreeVector bufferedDirection(0,0,0);
    259269    for(i0=0; i0<nDef; i0++)
     
    291301            G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge());
    292302            G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber();
    293             G4double concreteMass = G4NucleiPropertiesTable::GetNuclearMass(z1, a1);
     303            G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1);
    294304            G4double availableEnergy = eKinetic+mn+localMass-m1-m2-concreteMass;
    295305            // available kinetic energy in CMS (non relativistic)
     
    357367    G4double eBindN = 0;
    358368    G4double eBindP = 0;
    359     G4double eBindD = G4NucleiPropertiesTable::GetBindingEnergy(1,2);
    360     G4double eBindT = G4NucleiPropertiesTable::GetBindingEnergy(1,3);
    361     G4double eBindHe3 = G4NucleiPropertiesTable::GetBindingEnergy(2,3);
    362     G4double eBindA = G4NucleiPropertiesTable::GetBindingEnergy(2,4);
     369    G4double eBindD = G4NucleiProperties::GetBindingEnergy(2,1);
     370    G4double eBindT = G4NucleiProperties::GetBindingEnergy(3,1);
     371    G4double eBindHe3 = G4NucleiProperties::GetBindingEnergy(3,2);
     372    G4double eBindA = G4NucleiProperties::GetBindingEnergy(4,2);
    363373    for(i=0; i<tmpHadrons->size(); i++)
    364374    {
     
    455465  delete tmpHadrons;
    456466
     467//080721
     468   G4ParticleDefinition* targ_pd = G4ParticleTable::GetParticleTable()->GetIon ( (G4int)theBaseZ , (G4int)theBaseA , 0.0 );
     469   G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
     470   G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
     471   G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
     472   adjust_final_state ( init_4p_lab );
     473
    457474// clean up the primary neutron
    458475  theResult.SetStatusChange(stopAndKill);
Note: See TracChangeset for help on using the changeset viewer.