[819] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | // |
---|
| 27 | // 05-11-21 NeutronHP or Low Energy Parameterization Models |
---|
| 28 | // Implemented by T. Koi (SLAC/SCCS) |
---|
| 29 | // If NeutronHP data do not available for an element, then Low Energy |
---|
| 30 | // Parameterization models handle the interactions of the element. |
---|
| 31 | // |
---|
| 32 | |
---|
| 33 | // neutron_hp -- source file |
---|
| 34 | // J.P. Wellisch, Nov-1996 |
---|
| 35 | // A prototype of the low energy neutron transport model. |
---|
| 36 | // |
---|
| 37 | #include "G4NeutronHPorLEInelastic.hh" |
---|
| 38 | //#include "G4NeutronHPInelasticFS.hh" |
---|
| 39 | |
---|
| 40 | G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic() |
---|
| 41 | :G4HadronicInteraction("NeutronHPorLEInelastic") |
---|
| 42 | { |
---|
| 43 | SetMinEnergy(0.*eV); |
---|
| 44 | SetMaxEnergy(20.*MeV); |
---|
| 45 | |
---|
| 46 | // G4NeutronHPInelasticFS * theFS = new G4NeutronHPInelasticFS; |
---|
| 47 | if(!getenv("G4NEUTRONHPDATA")) |
---|
| 48 | throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files."); |
---|
| 49 | dirName = getenv("G4NEUTRONHPDATA"); |
---|
| 50 | G4String tString = "/Inelastic/"; |
---|
| 51 | dirName = dirName + tString; |
---|
| 52 | // G4cout <<"G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic testit "<<dirName<<G4endl; |
---|
| 53 | numEle = G4Element::GetNumberOfElements(); |
---|
| 54 | theInelastic = new G4NeutronHPChannelList[numEle]; |
---|
| 55 | unavailable_elements.clear(); |
---|
| 56 | for (G4int i=0; i<numEle; i++) |
---|
| 57 | { |
---|
| 58 | theInelastic[i].Init( (*(G4Element::GetElementTable()))[i] , dirName ); |
---|
| 59 | do |
---|
| 60 | { |
---|
| 61 | |
---|
| 62 | try |
---|
| 63 | { |
---|
| 64 | theInelastic[i].Register(&theNFS, "F01"); // has |
---|
| 65 | theInelastic[i].Register(&theNXFS, "F02"); |
---|
| 66 | theInelastic[i].Register(&the2NDFS, "F03"); |
---|
| 67 | theInelastic[i].Register(&the2NFS, "F04"); // has, E Done |
---|
| 68 | theInelastic[i].Register(&the3NFS, "F05"); // has, E Done |
---|
| 69 | theInelastic[i].Register(&theNAFS, "F06"); |
---|
| 70 | theInelastic[i].Register(&theN3AFS, "F07"); |
---|
| 71 | theInelastic[i].Register(&the2NAFS, "F08"); |
---|
| 72 | theInelastic[i].Register(&the3NAFS, "F09"); |
---|
| 73 | theInelastic[i].Register(&theNPFS, "F10"); |
---|
| 74 | theInelastic[i].Register(&theN2AFS, "F11"); |
---|
| 75 | theInelastic[i].Register(&the2N2AFS, "F12"); |
---|
| 76 | theInelastic[i].Register(&theNDFS, "F13"); |
---|
| 77 | theInelastic[i].Register(&theNTFS, "F14"); |
---|
| 78 | theInelastic[i].Register(&theNHe3FS, "F15"); |
---|
| 79 | theInelastic[i].Register(&theND2AFS, "F16"); |
---|
| 80 | theInelastic[i].Register(&theNT2AFS, "F17"); |
---|
| 81 | theInelastic[i].Register(&the4NFS, "F18"); // has, E Done |
---|
| 82 | theInelastic[i].Register(&the2NPFS, "F19"); |
---|
| 83 | theInelastic[i].Register(&the3NPFS, "F20"); |
---|
| 84 | theInelastic[i].Register(&theN2PFS, "F21"); |
---|
| 85 | theInelastic[i].Register(&theNPAFS, "F22"); |
---|
| 86 | theInelastic[i].Register(&thePFS, "F23"); |
---|
| 87 | theInelastic[i].Register(&theDFS, "F24"); |
---|
| 88 | theInelastic[i].Register(&theTFS, "F25"); |
---|
| 89 | theInelastic[i].Register(&theHe3FS, "F26"); |
---|
| 90 | theInelastic[i].Register(&theAFS, "F27"); |
---|
| 91 | theInelastic[i].Register(&the2AFS, "F28"); |
---|
| 92 | theInelastic[i].Register(&the3AFS, "F29"); |
---|
| 93 | theInelastic[i].Register(&the2PFS, "F30"); |
---|
| 94 | theInelastic[i].Register(&thePAFS, "F31"); |
---|
| 95 | theInelastic[i].Register(&theD2AFS, "F32"); |
---|
| 96 | theInelastic[i].Register(&theT2AFS, "F33"); |
---|
| 97 | theInelastic[i].Register(&thePDFS, "F34"); |
---|
| 98 | theInelastic[i].Register(&thePTFS, "F35"); |
---|
| 99 | theInelastic[i].Register(&theDAFS, "F36"); |
---|
| 100 | } |
---|
| 101 | |
---|
| 102 | catch ( G4HadronicException ) |
---|
| 103 | { |
---|
| 104 | unavailable_elements.insert ( (*(G4Element::GetElementTable()))[i]->GetName() ); |
---|
| 105 | } |
---|
| 106 | theInelastic[i].RestartRegistration(); |
---|
| 107 | } |
---|
| 108 | while( !theInelastic[i].HasDataInAnyFinalState()); |
---|
| 109 | |
---|
| 110 | } |
---|
| 111 | |
---|
| 112 | // delete theFS; |
---|
| 113 | if ( unavailable_elements.size() > 0 ) |
---|
| 114 | { |
---|
| 115 | std::set< G4String>::iterator it; |
---|
| 116 | G4cout << "HP Inelastic data are not available for thess elements "<< G4endl; |
---|
| 117 | for ( it = unavailable_elements.begin() ; it != unavailable_elements.end() ; it++ ) |
---|
| 118 | G4cout << *it << G4endl; |
---|
| 119 | G4cout << "Low Energy Parameterization Models will be used."<< G4endl; |
---|
| 120 | } |
---|
| 121 | |
---|
| 122 | createXSectionDataSet(); |
---|
| 123 | } |
---|
| 124 | |
---|
| 125 | G4NeutronHPorLEInelastic::~G4NeutronHPorLEInelastic() |
---|
| 126 | { |
---|
| 127 | delete [] theInelastic; |
---|
| 128 | delete theDataSet; |
---|
| 129 | } |
---|
| 130 | |
---|
| 131 | #include "G4NeutronHPThermalBoost.hh" |
---|
| 132 | |
---|
| 133 | G4HadFinalState * G4NeutronHPorLEInelastic::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& ) |
---|
| 134 | { |
---|
| 135 | G4int it=0; |
---|
| 136 | const G4Material * theMaterial = aTrack.GetMaterial(); |
---|
| 137 | G4int n = theMaterial->GetNumberOfElements(); |
---|
| 138 | G4int index = theMaterial->GetElement(0)->GetIndex(); |
---|
| 139 | if(n!=1) |
---|
| 140 | { |
---|
| 141 | G4int i; |
---|
| 142 | xSec = new G4double[n]; |
---|
| 143 | G4double sum=0; |
---|
| 144 | const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume(); |
---|
| 145 | G4double rWeight; |
---|
| 146 | G4NeutronHPThermalBoost aThermalE; |
---|
| 147 | for (i=0; i<n; i++) |
---|
| 148 | { |
---|
| 149 | index = theMaterial->GetElement(i)->GetIndex(); |
---|
| 150 | rWeight = NumAtomsPerVolume[i]; |
---|
| 151 | G4double x = aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i), theMaterial->GetTemperature()); |
---|
| 152 | |
---|
| 153 | //xSec[i] = theInelastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack, |
---|
| 154 | // theMaterial->GetElement(i), |
---|
| 155 | // theMaterial->GetTemperature())); |
---|
| 156 | xSec[i] = theInelastic[index].GetXsec(x); |
---|
| 157 | |
---|
| 158 | xSec[i] *= rWeight; |
---|
| 159 | sum+=xSec[i]; |
---|
| 160 | } |
---|
| 161 | G4double random = G4UniformRand(); |
---|
| 162 | G4double running = 0; |
---|
| 163 | for (i=0; i<n; i++) |
---|
| 164 | { |
---|
| 165 | running += xSec[i]; |
---|
| 166 | index = theMaterial->GetElement(i)->GetIndex(); |
---|
| 167 | it = i; |
---|
| 168 | if(random<=running/sum) break; |
---|
| 169 | } |
---|
| 170 | delete [] xSec; |
---|
| 171 | // it is element-wise initialised. |
---|
| 172 | } |
---|
| 173 | //return theInelastic[index].ApplyYourself(aTrack); |
---|
| 174 | return theInelastic[index].ApplyYourself( theMaterial->GetElement(it) , aTrack ); |
---|
| 175 | } |
---|
| 176 | |
---|
| 177 | |
---|
| 178 | |
---|
| 179 | G4bool G4NeutronHPorLEInelastic::IsThisElementOK( G4String name ) |
---|
| 180 | { |
---|
| 181 | if ( unavailable_elements.find( name ) == unavailable_elements.end() ) |
---|
| 182 | return TRUE; |
---|
| 183 | else |
---|
| 184 | return FALSE; |
---|
| 185 | } |
---|
| 186 | |
---|
| 187 | |
---|
| 188 | |
---|
| 189 | void G4NeutronHPorLEInelastic::createXSectionDataSet() |
---|
| 190 | { |
---|
| 191 | theDataSet = new G4NeutronHPorLEInelasticData ( theInelastic , &unavailable_elements ); |
---|
| 192 | } |
---|