G4NeutronHPFissionFS.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // neutron_hp -- source file
00027 // J.P. Wellisch, Nov-1996
00028 // A prototype of the low energy neutron transport model.
00029 //
00030 // 12-Apr-06 fix in delayed neutron and photon emission without FS data by T. Koi
00031 // 07-Sep-11 M. Kelsey -- Follow change to G4HadFinalState interface
00032 
00033 #include "G4NeutronHPFissionFS.hh"
00034 #include "G4PhysicalConstants.hh"
00035 #include "G4Nucleus.hh"
00036 #include "G4DynamicParticleVector.hh"
00037 #include "G4NeutronHPFissionERelease.hh"
00038 #include "G4ParticleTable.hh"
00039 
00040  void G4NeutronHPFissionFS::Init (G4double A, G4double Z, G4int M, G4String & dirName, G4String & aFSType)
00041  {
00042     //G4cout << "G4NeutronHPFissionFS::Init " << A << " " << Z << " " << M << G4endl;
00043     theFS.Init(A, Z, M, dirName, aFSType);
00044     theFC.Init(A, Z, M, dirName, aFSType);
00045     theSC.Init(A, Z, M, dirName, aFSType);
00046     theTC.Init(A, Z, M, dirName, aFSType);
00047     theLC.Init(A, Z, M, dirName, aFSType);
00048 
00049     theFF.Init(A, Z, M, dirName, aFSType);
00050     if ( getenv("G4NEUTRONHP_PRODUCE_FISSION_FRAGMENTS") && theFF.HasFSData() ) 
00051     {
00052        G4cout << "Activate Fission Fragments Prodcution for the target isotope of " 
00053        << "Z = " << (G4int)Z
00054        << ", A = " << (G4int)A
00055        //<< "M = " << M
00056        << G4endl;
00057        G4cout << "As the result, delayed neutrons are omitted and they should be taken care by RadioaActiveDecay." 
00058        << G4endl;
00059        produceFissionFragments = true; 
00060     }
00061  }
00062  G4HadFinalState * G4NeutronHPFissionFS::ApplyYourself(const G4HadProjectile & theTrack)
00063  {  
00064  //G4cout << "G4NeutronHPFissionFS::ApplyYourself " << G4endl;
00065 // prepare neutron
00066    theResult.Clear();
00067    G4double eKinetic = theTrack.GetKineticEnergy();
00068    const G4HadProjectile *incidentParticle = &theTrack;
00069    G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
00070    theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
00071    theNeutron.SetKineticEnergy( eKinetic );
00072 
00073 // prepare target
00074    G4Nucleus aNucleus;
00075    G4ReactionProduct theTarget; 
00076    G4double targetMass = theFS.GetMass();
00077    G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
00078    theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
00079 
00080 // set neutron and target in the FS classes 
00081   theFS.SetNeutron(theNeutron);
00082   theFS.SetTarget(theTarget);
00083   theFC.SetNeutron(theNeutron);
00084   theFC.SetTarget(theTarget);
00085   theSC.SetNeutron(theNeutron);
00086   theSC.SetTarget(theTarget);
00087   theTC.SetNeutron(theNeutron);
00088   theTC.SetTarget(theTarget);
00089   theLC.SetNeutron(theNeutron);
00090   theLC.SetTarget(theTarget);
00091 
00092 
00093   theFF.SetNeutron(theNeutron);
00094   theFF.SetTarget(theTarget);
00095 
00096 //TKWORK 120531
00097 //G4cout << theTarget.GetDefinition() << G4endl; this should be NULL
00098 //G4cout << "Z = " << theBaseZ << ", A = " << theBaseA << ", M = " << theBaseM << G4endl;
00099 // theNDLDataZ,A,M should be filled in each FS (theFS, theFC, theSC, theTC, theLC and theFF)
00101 
00102 // boost to target rest system and decide on channel.
00103    theNeutron.Lorentz(theNeutron, -1*theTarget);
00104 
00105 // dice the photons
00106 
00107    G4DynamicParticleVector * thePhotons;    
00108    thePhotons = theFS.GetPhotons();
00109 
00110 // select the FS in charge
00111 
00112    eKinetic = theNeutron.GetKineticEnergy();    
00113    G4double xSec[4];
00114    xSec[0] = theFC.GetXsec(eKinetic);
00115    xSec[1] = xSec[0]+theSC.GetXsec(eKinetic);
00116    xSec[2] = xSec[1]+theTC.GetXsec(eKinetic);
00117    xSec[3] = xSec[2]+theLC.GetXsec(eKinetic);
00118    G4int it;
00119    unsigned int i=0;
00120    G4double random = G4UniformRand();
00121    if(xSec[3]==0) 
00122    {
00123      it=-1;
00124    }
00125    else
00126    {
00127      for(i=0; i<4; i++)
00128      {
00129        it =i;
00130        if(random<xSec[i]/xSec[3]) break;
00131      }
00132    }
00133 
00134 // dice neutron multiplicities, energies and momenta in Lab. @@
00135 // no energy conservation on an event-to-event basis. we rely on the data to be ok. @@
00136 // also for mean, we rely on the consistancy of the data. @@
00137 
00138    G4int Prompt=0, delayed=0, all=0;
00139    G4DynamicParticleVector * theNeutrons = 0;
00140    switch(it) // check logic, and ask, if partials can be assumed to correspond to individual particles @@@
00141    {
00142      case 0:
00143        theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0);
00144        if(Prompt==0&&delayed==0) Prompt=all;
00145        theNeutrons = theFC.ApplyYourself(Prompt); // delayed always in FS 
00146        // take 'U' into account explicitely (see 5.4) in the sampling of energy @@@@
00147        break;
00148      case 1:
00149        theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 1);
00150        if(Prompt==0&&delayed==0) Prompt=all;
00151        theNeutrons = theSC.ApplyYourself(Prompt); // delayed always in FS, off done in FSFissionFS
00152        break;
00153      case 2:
00154        theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 2);
00155        if(Prompt==0&&delayed==0) Prompt=all;
00156        theNeutrons = theTC.ApplyYourself(Prompt); // delayed always in FS
00157        break;
00158      case 3:
00159        theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 3);
00160        if(Prompt==0&&delayed==0) Prompt=all;
00161        theNeutrons = theLC.ApplyYourself(Prompt); // delayed always in FS
00162        break;
00163      default:
00164        break;
00165    }
00166 
00167 // dice delayed neutrons and photons, and fallback 
00168 // for Prompt in case channel had no FS data; add all paricles to FS.
00169 
00170    //TKWORK120531
00171    if ( produceFissionFragments ) delayed=0;
00172 
00173    G4double * theDecayConstants;
00174 
00175    if( theNeutrons != 0)
00176    {
00177      theDecayConstants = new G4double[delayed];
00178      //
00179      //110527TKDB  Unused codes, Detected by gcc4.6 compiler 
00180      //G4int nPhotons = 0;
00181      //if(thePhotons!=0) nPhotons = thePhotons->size();
00182      for(i=0; i<theNeutrons->size(); i++)
00183      {
00184        theResult.AddSecondary(theNeutrons->operator[](i));
00185      }
00186      delete theNeutrons;  
00187 
00188      G4DynamicParticleVector * theDelayed = 0;
00189 //   G4cout << "delayed" << G4endl;
00190      theDelayed = theFS.ApplyYourself(0, delayed, theDecayConstants);
00191      for(i=0; i<theDelayed->size(); i++)
00192      {
00193        G4double time = -std::log(G4UniformRand())/theDecayConstants[i];
00194        time += theTrack.GetGlobalTime();
00195        theResult.AddSecondary(theDelayed->operator[](i));
00196        theResult.GetSecondary(theResult.GetNumberOfSecondaries()-1)->SetTime(time);
00197      }
00198      delete theDelayed;                  
00199    }
00200    else
00201    {
00202 //    cout << " all = "<<all<<G4endl;
00203      theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0);
00204      theDecayConstants = new G4double[delayed];
00205      if(Prompt==0&&delayed==0) Prompt=all;
00206      theNeutrons = theFS.ApplyYourself(Prompt, delayed, theDecayConstants);
00207      //110527TKDB  Unused codes, Detected by gcc4.6 compiler 
00208      //G4int nPhotons = 0;
00209      //if(thePhotons!=0) nPhotons = thePhotons->size();
00210      G4int i0;
00211      for(i0=0; i0<Prompt; i0++)
00212      {
00213        theResult.AddSecondary(theNeutrons->operator[](i0));
00214      }
00215 
00216 //G4cout << "delayed" << G4endl;
00217      for(i0=Prompt; i0<Prompt+delayed; i0++)
00218      {
00219        G4double time = -std::log(G4UniformRand())/theDecayConstants[i0-Prompt];
00220        time += theTrack.GetGlobalTime();        
00221        theResult.AddSecondary(theNeutrons->operator[](i0));
00222        theResult.GetSecondary(theResult.GetNumberOfSecondaries()-1)->SetTime(time);
00223      }
00224      delete theNeutrons;   
00225    }
00226    delete [] theDecayConstants;
00227 //    cout << "all delayed "<<delayed<<G4endl; 
00228    unsigned int nPhotons = 0;
00229    if(thePhotons!=0)
00230    {
00231      nPhotons = thePhotons->size();
00232      for(i=0; i<nPhotons; i++)
00233      {
00234        theResult.AddSecondary(thePhotons->operator[](i));
00235      }
00236      delete thePhotons; 
00237    }
00238 
00239 // finally deal with local energy depositions.
00240 //    G4cout <<"Number of secondaries = "<<theResult.GetNumberOfSecondaries()<< G4endl;
00241 //    G4cout <<"Number of photons = "<<nPhotons<<G4endl;
00242 //    G4cout <<"Number of Prompt = "<<Prompt<<G4endl;
00243 //    G4cout <<"Number of delayed = "<<delayed<<G4endl;
00244 
00245    G4NeutronHPFissionERelease * theERelease = theFS.GetEnergyRelease();
00246    G4double eDepByFragments = theERelease->GetFragmentKinetic();
00247    //theResult.SetLocalEnergyDeposit(eDepByFragments);
00248    if ( !produceFissionFragments ) theResult.SetLocalEnergyDeposit(eDepByFragments);
00249 //    cout << "local energy deposit" << eDepByFragments<<G4endl;
00250 // clean up the primary neutron
00251    theResult.SetStatusChange(stopAndKill);
00252    //G4cout << "Prompt = " << Prompt << ", Delayed = " << delayed << ", All= " << all << G4endl;
00253    //G4cout << "local energy deposit " << eDepByFragments/MeV << "MeV " << G4endl;
00254 
00255    //TKWORK120531
00256    if ( produceFissionFragments )
00257    {
00258       G4int fragA_Z=0;
00259       G4int fragA_A=0;
00260       G4int fragA_M=0;
00261       // System is traget rest!
00262       theFF.GetAFissionFragment(eKinetic,fragA_Z,fragA_A,fragA_M);
00263       G4int fragB_Z=(G4int)theBaseZ-fragA_Z;
00264       G4int fragB_A=(G4int)theBaseA-fragA_A-Prompt;
00265       //fragA_M ignored 
00266       //G4int fragB_M=theBaseM-fragA_M; 
00267       //G4cout << fragA_Z << " " << fragA_A << " " << fragA_M << G4endl;
00268       //G4cout << fragB_Z << " " << fragB_A << G4endl;
00269 
00270       G4ParticleTable* pt = G4ParticleTable::GetParticleTable();
00271       //Excitation energy is not taken into account 
00272       G4ParticleDefinition* pdA = pt->GetIon( fragA_Z , fragA_A , 0.0 );
00273       G4ParticleDefinition* pdB = pt->GetIon( fragB_Z , fragB_A , 0.0 );
00274 
00275       //Isotropic Distribution 
00276       G4double phi = twopi*G4UniformRand();
00277       G4double theta = pi*G4UniformRand();
00278       G4double sinth = std::sin(theta);
00279       G4ThreeVector direction (sinth*std::cos(phi) , sinth*std::sin(phi), std::cos(theta) );
00280 
00281       // Just use ENDF value for this 
00282       G4double ER = eDepByFragments; 
00283       G4double ma = pdA->GetPDGMass();
00284       G4double mb = pdB->GetPDGMass();
00285       G4double EA = ER / ( 1 + ma/mb);
00286       G4double EB = ER - EA;
00287       G4DynamicParticle* dpA = new G4DynamicParticle( pdA , direction , EA);
00288       G4DynamicParticle* dpB = new G4DynamicParticle( pdB , -direction , EB);
00289       theResult.AddSecondary(dpA);
00290       theResult.AddSecondary(dpB);
00291    }
00292    //TKWORK 120531 END
00293 
00294    return &theResult;
00295  }

Generated on Mon May 27 17:49:01 2013 for Geant4 by  doxygen 1.4.7