#include <DsG4NNDCCaptureGammas.h>
Public Member Functions | |
DsG4NNDCCaptureGammas () | |
~DsG4NNDCCaptureGammas () | |
void | Init (int A, int Z) |
G4ReactionProductVector * | GetGammas () |
std::vector< double > | GetEnergy () |
bool | hasData () |
Private Attributes | |
bool | m_hasData |
std::vector< double > | m_probabilities |
std::vector< std::vector< double > > | m_gammaEnergies |
Definition at line 15 of file DsG4NNDCCaptureGammas.h.
DsG4NNDCCaptureGammas::DsG4NNDCCaptureGammas | ( | ) |
Definition at line 23 of file DsG4NNDCCaptureGammas.cc.
00023 : 00024 m_hasData(false) 00025 { 00026 }
DsG4NNDCCaptureGammas::~DsG4NNDCCaptureGammas | ( | ) |
void DsG4NNDCCaptureGammas::Init | ( | int | A, | |
int | Z | |||
) |
Definition at line 32 of file DsG4NNDCCaptureGammas.cc.
00033 { 00034 // Gamma files should be placed in $DAYA_DATA_DIR directory. 00035 00036 // define a map for Z 00037 map<int, string> elementName; 00038 elementName[64] = "Gd"; 00039 elementName[26] = "Fe"; 00040 elementName[28] = "Ni"; 00041 elementName[14] = "Si"; 00042 elementName[15] = "P"; 00043 elementName[25] = "Mn"; 00044 elementName[16] = "S"; 00045 elementName[24] = "Cr"; 00046 elementName[8] = "O"; 00047 elementName[7] = "N"; 00048 elementName[6] = "C"; 00049 00050 string fileName; 00051 fileName.append(elementName[Z]); 00052 ostringstream strA; 00053 strA << A; 00054 fileName.append(strA.str()); 00055 fileName.append("NeutronCaptureGammas.txt"); 00056 00057 const char* prefix = getenv("DAYA_DATA_DIR"); 00058 if(prefix) { 00059 if(strlen(prefix)>0) { 00060 std::string newname = prefix; 00061 newname += "/"; 00062 newname += fileName; 00063 fileName = newname; 00064 } 00065 } 00066 00067 ifstream file(fileName.c_str()); 00068 00069 //read data file 00070 string line; 00071 double totalProbability = 0; 00072 while(getline(file,line)) 00073 { 00074 stringstream streamLine(line); 00075 double probability = 0.0; 00076 int numGammas = 0; 00077 if(streamLine.peek() != '#') 00078 { 00079 streamLine >> probability >> numGammas; 00080 m_probabilities.push_back(probability); 00081 totalProbability += probability; 00082 vector<double> gammas; 00083 for(int i = 1; i <= numGammas; i++) 00084 { 00085 double gammaEnergy = 0.0; 00086 streamLine >> gammaEnergy; 00087 gammas.push_back(gammaEnergy/1000.0); // Convert to MeV 00088 } 00089 m_gammaEnergies.push_back(gammas); 00090 } 00091 } 00092 file.close(); 00093 00094 if( m_probabilities.size() > 0){ 00095 std::cout << "Loaded " << m_probabilities.size() 00096 << " gamma decay combinations from file " << fileName << endl; 00097 std::cout << "Rescaling all gamma lines by total probability: " 00098 << totalProbability << endl; 00099 for(unsigned int i=0; i<m_probabilities.size(); i++){ 00100 m_probabilities[i] /= totalProbability; 00101 } 00102 m_hasData = true; 00103 }else{ 00104 std::cout << "DsG4NNDCCaptureGammas WARNING: No data loaded for isotope " 00105 << "A=" << A << " Z=" << Z 00106 << ". Cross section will be set to zero." 00107 << endl; 00108 m_hasData = false; 00109 } 00110 }
G4ReactionProductVector * DsG4NNDCCaptureGammas::GetGammas | ( | ) |
Definition at line 112 of file DsG4NNDCCaptureGammas.cc.
00113 { 00114 G4ReactionProductVector* theGammas = new G4ReactionProductVector; 00115 vector<double> energy = GetEnergy(); 00116 for(unsigned int i=0; i<energy.size(); i++) { 00117 00118 G4ReactionProduct* theOne = new G4ReactionProduct; 00119 theOne->SetDefinition( G4Gamma::Gamma() ); 00120 00121 // Get the gammas direction. 00122 // Isotropic emission. 00123 G4double costheta = 2.*G4UniformRand()-1; 00124 G4double theta = acos(costheta); 00125 G4double phi = twopi*G4UniformRand(); 00126 G4double sinth = sin(theta); 00127 theOne->SetTotalEnergy( energy[i] ); 00128 G4ThreeVector temp(energy[i]*sinth*cos(phi), 00129 energy[i]*sinth*sin(phi), 00130 energy[i]*cos(theta) ); 00131 theOne->SetMomentum( temp ) ; 00132 theGammas->push_back(theOne); 00133 } 00134 return theGammas; 00135 }
vector< double > DsG4NNDCCaptureGammas::GetEnergy | ( | ) |
Definition at line 137 of file DsG4NNDCCaptureGammas.cc.
00138 { 00139 //pick a decay scheme 00140 G4double random = G4UniformRand(); 00141 double sum = 0.0; 00142 unsigned int i = 0; 00143 while(i < m_probabilities.size() ) 00144 { 00145 sum += m_probabilities[i]; 00146 if(sum >= random) 00147 { 00148 // cout << "Returning gamma branch " << i << endl; 00149 return m_gammaEnergies[i]; 00150 } 00151 i++; 00152 } 00153 std::cerr << "Failed to generate gammas!" << endl; 00154 return vector<double>(0); 00155 }
bool DsG4NNDCCaptureGammas::hasData | ( | ) | [inline] |
bool DsG4NNDCCaptureGammas::m_hasData [private] |
Definition at line 28 of file DsG4NNDCCaptureGammas.h.
std::vector<double> DsG4NNDCCaptureGammas::m_probabilities [private] |
Definition at line 29 of file DsG4NNDCCaptureGammas.h.
std::vector< std::vector<double> > DsG4NNDCCaptureGammas::m_gammaEnergies [private] |
Definition at line 30 of file DsG4NNDCCaptureGammas.h.