G4NeutronHPInelastic.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 // this code implementation is the intellectual property of
00027 // neutron_hp -- source file
00028 // J.P. Wellisch, Nov-1996
00029 // A prototype of the low energy neutron transport model.
00030 //
00031 // By copying, distributing or modifying the Program (or any work
00032 // based on the Program) you indicate your acceptance of this statement,
00033 // and all its terms.
00034 //
00035 // $Id$
00036 //
00037 // 070523 bug fix for G4FPE_DEBUG on by A. Howard (and T. Koi)
00038 // 081203 limit maximum trial for creating final states add protection for 1H isotope case by T. Koi
00039 //
00040 #include "G4NeutronHPInelastic.hh"
00041 #include "G4SystemOfUnits.hh"
00042 
00043 #include "G4NeutronHPManager.hh"
00044 
00045   G4NeutronHPInelastic::G4NeutronHPInelastic()
00046     :G4HadronicInteraction("NeutronHPInelastic")
00047   {
00048     SetMinEnergy( 0.0 );
00049     SetMaxEnergy( 20.*MeV );
00050 
00051     G4int istatus = system("echo $G4NEUTRONHPDATA");
00052     if ( istatus < 0 )
00053     {
00054       G4cout << "Warning! system(\"echo $G4NEUTRONHPDATA\") returns error value at G4NeutronHPInelastic" << G4endl;
00055     } 
00056 
00057 //    G4cout << " entering G4NeutronHPInelastic constructor"<<G4endl;
00058     if(!getenv("G4NEUTRONHPDATA")) 
00059        throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
00060     dirName = getenv("G4NEUTRONHPDATA");
00061     G4String tString = "/Inelastic";
00062     dirName = dirName + tString;
00063     numEle = G4Element::GetNumberOfElements();
00064 /*
00065     theInelastic = new G4NeutronHPChannelList[numEle];
00066     for (G4int i=0; i<numEle; i++)
00067     { 
00068       theInelastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
00069       G4int itry = 0;
00070       do
00071       {
00072         theInelastic[i].Register(&theNFS, "F01"); // has
00073         theInelastic[i].Register(&theNXFS, "F02");
00074         theInelastic[i].Register(&the2NDFS, "F03");
00075         theInelastic[i].Register(&the2NFS, "F04"); // has, E Done
00076         theInelastic[i].Register(&the3NFS, "F05"); // has, E Done
00077         theInelastic[i].Register(&theNAFS, "F06");
00078         theInelastic[i].Register(&theN3AFS, "F07");
00079         theInelastic[i].Register(&the2NAFS, "F08");
00080         theInelastic[i].Register(&the3NAFS, "F09");
00081         theInelastic[i].Register(&theNPFS, "F10");
00082         theInelastic[i].Register(&theN2AFS, "F11");
00083         theInelastic[i].Register(&the2N2AFS, "F12");
00084         theInelastic[i].Register(&theNDFS, "F13");
00085         theInelastic[i].Register(&theNTFS, "F14");
00086         theInelastic[i].Register(&theNHe3FS, "F15");
00087         theInelastic[i].Register(&theND2AFS, "F16");
00088         theInelastic[i].Register(&theNT2AFS, "F17");
00089         theInelastic[i].Register(&the4NFS, "F18"); // has, E Done
00090         theInelastic[i].Register(&the2NPFS, "F19");
00091         theInelastic[i].Register(&the3NPFS, "F20");
00092         theInelastic[i].Register(&theN2PFS, "F21");
00093         theInelastic[i].Register(&theNPAFS, "F22");
00094         theInelastic[i].Register(&thePFS, "F23");
00095         theInelastic[i].Register(&theDFS, "F24");
00096         theInelastic[i].Register(&theTFS, "F25");
00097         theInelastic[i].Register(&theHe3FS, "F26");
00098         theInelastic[i].Register(&theAFS, "F27");
00099         theInelastic[i].Register(&the2AFS, "F28");
00100         theInelastic[i].Register(&the3AFS, "F29");
00101         theInelastic[i].Register(&the2PFS, "F30");
00102         theInelastic[i].Register(&thePAFS, "F31");
00103         theInelastic[i].Register(&theD2AFS, "F32");
00104         theInelastic[i].Register(&theT2AFS, "F33");
00105         theInelastic[i].Register(&thePDFS, "F34");
00106         theInelastic[i].Register(&thePTFS, "F35");
00107         theInelastic[i].Register(&theDAFS, "F36");
00108         theInelastic[i].RestartRegistration();
00109         itry++;
00110       }
00111       //while(!theInelastic[i].HasDataInAnyFinalState());
00112       while( !theInelastic[i].HasDataInAnyFinalState() && itry < 6 );
00113                                                               // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK  
00114 
00115       if ( itry == 6 ) 
00116       {
00117          // No Final State at all.
00118          G4bool exceptional = false;
00119          if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
00120          {
00121             if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true;  //1H
00122          } 
00123          if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
00124       }
00125     }
00126 */
00127 
00128     for (G4int i=0; i<numEle; i++)
00129     { 
00130       theInelastic.push_back( new G4NeutronHPChannelList );
00131       (*theInelastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
00132       G4int itry = 0;
00133       do
00134       {
00135         (*theInelastic[i]).Register(&theNFS, "F01"); // has
00136         (*theInelastic[i]).Register(&theNXFS, "F02");
00137         (*theInelastic[i]).Register(&the2NDFS, "F03");
00138         (*theInelastic[i]).Register(&the2NFS, "F04"); // has, E Done
00139         (*theInelastic[i]).Register(&the3NFS, "F05"); // has, E Done
00140         (*theInelastic[i]).Register(&theNAFS, "F06");
00141         (*theInelastic[i]).Register(&theN3AFS, "F07");
00142         (*theInelastic[i]).Register(&the2NAFS, "F08");
00143         (*theInelastic[i]).Register(&the3NAFS, "F09");
00144         (*theInelastic[i]).Register(&theNPFS, "F10");
00145         (*theInelastic[i]).Register(&theN2AFS, "F11");
00146         (*theInelastic[i]).Register(&the2N2AFS, "F12");
00147         (*theInelastic[i]).Register(&theNDFS, "F13");
00148         (*theInelastic[i]).Register(&theNTFS, "F14");
00149         (*theInelastic[i]).Register(&theNHe3FS, "F15");
00150         (*theInelastic[i]).Register(&theND2AFS, "F16");
00151         (*theInelastic[i]).Register(&theNT2AFS, "F17");
00152         (*theInelastic[i]).Register(&the4NFS, "F18"); // has, E Done
00153         (*theInelastic[i]).Register(&the2NPFS, "F19");
00154         (*theInelastic[i]).Register(&the3NPFS, "F20");
00155         (*theInelastic[i]).Register(&theN2PFS, "F21");
00156         (*theInelastic[i]).Register(&theNPAFS, "F22");
00157         (*theInelastic[i]).Register(&thePFS, "F23");
00158         (*theInelastic[i]).Register(&theDFS, "F24");
00159         (*theInelastic[i]).Register(&theTFS, "F25");
00160         (*theInelastic[i]).Register(&theHe3FS, "F26");
00161         (*theInelastic[i]).Register(&theAFS, "F27");
00162         (*theInelastic[i]).Register(&the2AFS, "F28");
00163         (*theInelastic[i]).Register(&the3AFS, "F29");
00164         (*theInelastic[i]).Register(&the2PFS, "F30");
00165         (*theInelastic[i]).Register(&thePAFS, "F31");
00166         (*theInelastic[i]).Register(&theD2AFS, "F32");
00167         (*theInelastic[i]).Register(&theT2AFS, "F33");
00168         (*theInelastic[i]).Register(&thePDFS, "F34");
00169         (*theInelastic[i]).Register(&thePTFS, "F35");
00170         (*theInelastic[i]).Register(&theDAFS, "F36");
00171         (*theInelastic[i]).RestartRegistration();
00172         itry++;
00173       }
00174       while( !(*theInelastic[i]).HasDataInAnyFinalState() && itry < 6 );
00175                                                               // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK  
00176 
00177       if ( itry == 6 ) 
00178       {
00179          // No Final State at all.
00180          G4bool exceptional = false;
00181          if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
00182          {
00183             if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true;  //1H
00184          } 
00185          if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
00186       }
00187 
00188     }
00189   }
00190 
00191   G4NeutronHPInelastic::~G4NeutronHPInelastic()
00192   {
00193 //    delete [] theInelastic;
00194      for ( std::vector<G4NeutronHPChannelList*>::iterator 
00195            it = theInelastic.begin() ; it != theInelastic.end() ; it++ )
00196      {
00197         delete *it;
00198      }
00199      theInelastic.clear();
00200   }
00201   
00202   #include "G4NeutronHPThermalBoost.hh"
00203   
00204   G4HadFinalState * G4NeutronHPInelastic::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aNucleus )
00205   {
00206     if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
00207     G4NeutronHPManager::GetInstance()->OpenReactionWhiteBoard();
00208     const G4Material * theMaterial = aTrack.GetMaterial();
00209     G4int n = theMaterial->GetNumberOfElements();
00210     G4int index = theMaterial->GetElement(0)->GetIndex();
00211     G4int it=0;
00212     if(n!=1)
00213     {
00214       xSec = new G4double[n];
00215       G4double sum=0;
00216       G4int i;
00217       const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
00218       G4double rWeight;    
00219       G4NeutronHPThermalBoost aThermalE;
00220       for (i=0; i<n; i++)
00221       {
00222         index = theMaterial->GetElement(i)->GetIndex();
00223         rWeight = NumAtomsPerVolume[i];
00224         //xSec[i] = theInelastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
00225         xSec[i] = (*theInelastic[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
00226                                                                          theMaterial->GetElement(i),
00227                                                                          theMaterial->GetTemperature()));
00228         xSec[i] *= rWeight;
00229         sum+=xSec[i];
00230       }
00231       G4double random = G4UniformRand();
00232       G4double running = 0;
00233       for (i=0; i<n; i++)
00234       {
00235         running += xSec[i];
00236         index = theMaterial->GetElement(i)->GetIndex();
00237         it = i;
00238         //if(random<=running/sum) break;
00239         if( sum == 0 || random<=running/sum) break;
00240       }
00241       delete [] xSec;
00242     }
00243 
00244     //return theInelastic[index].ApplyYourself(theMaterial->GetElement(it), aTrack);
00245     G4HadFinalState* result = (*theInelastic[index]).ApplyYourself(theMaterial->GetElement(it), aTrack);
00246     //
00247     aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
00248     G4NeutronHPManager::GetInstance()->CloseReactionWhiteBoard();
00249     return result;
00250   }
00251 
00252 const std::pair<G4double, G4double> G4NeutronHPInelastic::GetFatalEnergyCheckLevels() const
00253 {
00254       // max energy non-conservation is mass of heavy nucleus
00255 //      if ( getenv("G4NEUTRONHP_DO_NOT_ADJUST_FINAL_STATE") ) return std::pair<G4double, G4double>(5*perCent,250*GeV);
00256       // This should be same to the hadron default value
00257 //      return std::pair<G4double, G4double>(10*perCent,10*GeV);
00258       return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
00259 }
00260 
00261 void G4NeutronHPInelastic::addChannelForNewElement()
00262 {
00263    for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) 
00264    {
00265       G4cout << "G4NeutronHPInelastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
00266 
00267       theInelastic.push_back( new G4NeutronHPChannelList );
00268       (*theInelastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
00269       G4int itry = 0;
00270       do
00271       {
00272          (*theInelastic[i]).Register(&theNFS, "F01"); // has
00273          (*theInelastic[i]).Register(&theNXFS, "F02");
00274          (*theInelastic[i]).Register(&the2NDFS, "F03");
00275          (*theInelastic[i]).Register(&the2NFS, "F04"); // has, E Done
00276          (*theInelastic[i]).Register(&the3NFS, "F05"); // has, E Done
00277          (*theInelastic[i]).Register(&theNAFS, "F06");
00278          (*theInelastic[i]).Register(&theN3AFS, "F07");
00279          (*theInelastic[i]).Register(&the2NAFS, "F08");
00280          (*theInelastic[i]).Register(&the3NAFS, "F09");
00281          (*theInelastic[i]).Register(&theNPFS, "F10");
00282          (*theInelastic[i]).Register(&theN2AFS, "F11");
00283          (*theInelastic[i]).Register(&the2N2AFS, "F12");
00284          (*theInelastic[i]).Register(&theNDFS, "F13");
00285          (*theInelastic[i]).Register(&theNTFS, "F14");
00286          (*theInelastic[i]).Register(&theNHe3FS, "F15");
00287          (*theInelastic[i]).Register(&theND2AFS, "F16");
00288          (*theInelastic[i]).Register(&theNT2AFS, "F17");
00289          (*theInelastic[i]).Register(&the4NFS, "F18"); // has, E Done
00290          (*theInelastic[i]).Register(&the2NPFS, "F19");
00291          (*theInelastic[i]).Register(&the3NPFS, "F20");
00292          (*theInelastic[i]).Register(&theN2PFS, "F21");
00293          (*theInelastic[i]).Register(&theNPAFS, "F22");
00294          (*theInelastic[i]).Register(&thePFS, "F23");
00295          (*theInelastic[i]).Register(&theDFS, "F24");
00296          (*theInelastic[i]).Register(&theTFS, "F25");
00297          (*theInelastic[i]).Register(&theHe3FS, "F26");
00298          (*theInelastic[i]).Register(&theAFS, "F27");
00299          (*theInelastic[i]).Register(&the2AFS, "F28");
00300          (*theInelastic[i]).Register(&the3AFS, "F29");
00301          (*theInelastic[i]).Register(&the2PFS, "F30");
00302          (*theInelastic[i]).Register(&thePAFS, "F31");
00303          (*theInelastic[i]).Register(&theD2AFS, "F32");
00304          (*theInelastic[i]).Register(&theT2AFS, "F33");
00305          (*theInelastic[i]).Register(&thePDFS, "F34");
00306          (*theInelastic[i]).Register(&thePTFS, "F35");
00307          (*theInelastic[i]).Register(&theDAFS, "F36");
00308          (*theInelastic[i]).RestartRegistration();
00309          itry++;
00310       }
00311       while( !(*theInelastic[i]).HasDataInAnyFinalState() && itry < 6 );
00312                                                                   // 6 is corresponding to the value(5) of G4NeutronHPChannel. TK  
00313 
00314       if ( itry == 6 ) 
00315       {
00316          // No Final State at all.
00317          G4bool exceptional = false;
00318          if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
00319          {
00320             if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true;  //1H
00321          } 
00322          if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
00323       }
00324    }
00325 
00326    numEle = (G4int)G4Element::GetNumberOfElements();
00327 }

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