G4NeutronHPFission.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 // 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
00031 // 08-08-06 delete unnecessary and harmed declaration; Bug Report[857]
00032 //
00033 #include "G4NeutronHPFission.hh"
00034 #include "G4SystemOfUnits.hh"
00035 
00036 #include "G4NeutronHPManager.hh"
00037 
00038   G4NeutronHPFission::G4NeutronHPFission()
00039     :G4HadronicInteraction("NeutronHPFission")
00040   {
00041     SetMinEnergy( 0.0 );
00042     SetMaxEnergy( 20.*MeV );
00043     if(!getenv("G4NEUTRONHPDATA")) 
00044        throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
00045     dirName = getenv("G4NEUTRONHPDATA");
00046     G4String tString = "/Fission";
00047     dirName = dirName + tString;
00048     numEle = G4Element::GetNumberOfElements();
00049     //theFission = new G4NeutronHPChannel[numEle];
00050 
00051     //for (G4int i=0; i<numEle; i++)
00052     //{ 
00053       //if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
00054     //  if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
00055     //  {
00056     //    theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
00057     //    theFission[i].Register(&theFS);
00058     //  }
00059     //}
00060 
00061     for ( G4int i = 0 ; i < numEle ; i++ ) 
00062     {
00063       theFission.push_back( new G4NeutronHPChannel );
00064       if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
00065       {
00066        (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
00067        (*theFission[i]).Register(&theFS);
00068       }
00069     }
00070   }
00071   
00072   G4NeutronHPFission::~G4NeutronHPFission()
00073   {
00074     //delete [] theFission;
00075      for ( std::vector<G4NeutronHPChannel*>::iterator 
00076            it = theFission.begin() ; it != theFission.end() ; it++ )
00077      {
00078         delete *it;
00079      }
00080      theFission.clear();
00081   }
00082   
00083   #include "G4NeutronHPThermalBoost.hh"
00084   G4HadFinalState * G4NeutronHPFission::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aNucleus )
00085   {
00086 
00087     if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
00088 
00089     G4NeutronHPManager::GetInstance()->OpenReactionWhiteBoard();
00090     const G4Material * theMaterial = aTrack.GetMaterial();
00091     G4int n = theMaterial->GetNumberOfElements();
00092     G4int index = theMaterial->GetElement(0)->GetIndex();
00093     if(n!=1)
00094     {
00095       xSec = new G4double[n];
00096       G4double sum=0;
00097       G4int i;
00098       const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
00099       G4double rWeight;    
00100       G4NeutronHPThermalBoost aThermalE;
00101       for (i=0; i<n; i++)
00102       {
00103         index = theMaterial->GetElement(i)->GetIndex();
00104         rWeight = NumAtomsPerVolume[i];
00105         xSec[i] = (*theFission[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
00106                                                                       theMaterial->GetElement(i),
00107                                                                       theMaterial->GetTemperature()));
00108         xSec[i] *= rWeight;
00109         sum+=xSec[i];
00110       }
00111       G4double random = G4UniformRand();
00112       G4double running = 0;
00113       for (i=0; i<n; i++)
00114       {
00115         running += xSec[i];
00116         index = theMaterial->GetElement(i)->GetIndex();
00117         //if(random<=running/sum) break;
00118         if( sum == 0 ||  random <= running/sum ) break;
00119       }
00120       delete [] xSec;
00121     }
00122     //return theFission[index].ApplyYourself(aTrack);
00123     G4HadFinalState* result = (*theFission[index]).ApplyYourself(aTrack);
00124     aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
00125     G4NeutronHPManager::GetInstance()->CloseReactionWhiteBoard();
00126     return result; 
00127   }
00128 
00129 const std::pair<G4double, G4double> G4NeutronHPFission::GetFatalEnergyCheckLevels() const
00130 {
00131         // max energy non-conservation is mass of heavy nucleus
00132         //return std::pair<G4double, G4double>(5*perCent,250*GeV);
00133         return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
00134 }
00135 
00136 
00137 
00138 void G4NeutronHPFission::addChannelForNewElement()
00139 {
00140    for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) 
00141    {
00142       theFission.push_back( new G4NeutronHPChannel );
00143       if ( (*(G4Element::GetElementTable()))[i]->GetZ() > 87 ) //TK modified for ENDF-VII
00144       {
00145          G4cout << "G4NeutronHPFission Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
00146          (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
00147          (*theFission[i]).Register(&theFS);
00148       }
00149    }
00150    numEle = (G4int)G4Element::GetNumberOfElements();
00151 }

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