Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 // 26 // 27 // ------------------------------------------- 27 // ------------------------------------------------------------------- 28 // 28 // 29 // GEANT4 Class file 29 // GEANT4 Class file 30 // 30 // 31 // 31 // 32 // File name: G4PhotoElectricAngularGenera 32 // File name: G4PhotoElectricAngularGeneratorPolarized 33 // 33 // 34 // Author: A. C. Farinha, L. Peralta, P. Rodri 34 // Author: A. C. Farinha, L. Peralta, P. Rodrigues and A. Trindade 35 // << 35 // 36 // Creation date: 36 // Creation date: 37 // 37 // 38 // Modifications: << 38 // Modifications: 39 // 10 January 2006 << 39 // 10 January 2006 40 // 06 May 2011, Replace FILE with std::ifstrea << 41 // 40 // 42 // Class Description: << 41 // Class Description: 43 // 42 // 44 // Concrete class for PhotoElectric Electron A << 43 // Concrete class for PhotoElectric Electron Angular Polarized Distribution Generation 45 // 44 // 46 // Class Description: << 45 // Class Description: 47 // PhotoElectric Electron Angular Generator ba 46 // PhotoElectric Electron Angular Generator based on the general Gavrila photoelectron angular distribution. 48 // Includes polarization effects for K and L1 47 // Includes polarization effects for K and L1 atomic shells, according to Gavrila (1959, 1961). 49 // For higher shells the L1 cross-section is u << 48 // For higher shells the L1 cross-section is used. 50 // 49 // 51 // The Gavrila photoelectron angular distribut 50 // The Gavrila photoelectron angular distribution is a complex function which can not be sampled using 52 // the inverse-transform method (James 1980). << 51 // the inverse-transform method (James 1980). Instead a more general approach based on the one already 53 // used to sample bremsstrahlung 2BN cross sec 52 // used to sample bremsstrahlung 2BN cross section (G4Generator2BN, Peralta, 2005) was used. 54 // 53 // 55 // M. Gavrila, "Relativistic K-Shell Photoeffe 54 // M. Gavrila, "Relativistic K-Shell Photoeffect", Phys. Rev. 113, 514-526 (1959) 56 // M. Gavrila, "Relativistic L-Shell Photoeffe 55 // M. Gavrila, "Relativistic L-Shell Photoeffect", Phys. Rev. 124, 1132-1141 (1961) 57 // F. James, Rept. on Prog. in Phys. 43, 1145 56 // F. James, Rept. on Prog. in Phys. 43, 1145 (1980) 58 // L. Peralta et al., "A new low-energy bremss 57 // L. Peralta et al., "A new low-energy bremsstrahlung generator for GEANT4", Radiat. Prot. Dosimetry. 116, 59-64 (2005) 59 // 58 // 60 // 59 // 61 // ------------------------------------------- 60 // ------------------------------------------------------------------- 62 // 61 // 63 // << 62 // 64 63 65 #include "G4PhotoElectricAngularGeneratorPolar 64 #include "G4PhotoElectricAngularGeneratorPolarized.hh" 66 #include "G4PhysicalConstants.hh" << 67 #include "G4RotationMatrix.hh" 65 #include "G4RotationMatrix.hh" 68 #include "Randomize.hh" 66 #include "Randomize.hh" 69 #include "G4Exp.hh" << 70 67 71 G4PhotoElectricAngularGeneratorPolarized::G4Ph << 68 // 72 :G4VEmAngularDistribution("AngularGenSauterG << 69 >> 70 G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized(const G4String& name):G4VPhotoElectricAngularDistribution(name) 73 { 71 { 74 const G4int arrayDim = 980; 72 const G4int arrayDim = 980; 75 73 76 //minimum electron beta parameter allowed 74 //minimum electron beta parameter allowed 77 betaArray[0] = 0.02; 75 betaArray[0] = 0.02; 78 //beta step 76 //beta step 79 betaArray[1] = 0.001; 77 betaArray[1] = 0.001; 80 //maximum index array for a and c tables 78 //maximum index array for a and c tables 81 betaArray[2] = arrayDim - 1; 79 betaArray[2] = arrayDim - 1; 82 80 83 // read Majorant Surface Parameters. This ar << 81 // read Majorant Surface Parameters. This are required in order to generate Gavrila angular photoelectron distribution 84 //Gavrila angular photoelectron distribution << 85 for(G4int level = 0; level < 2; level++){ 82 for(G4int level = 0; level < 2; level++){ >> 83 86 char nameChar0[100] = "ftab0.dat"; // K-sh 84 char nameChar0[100] = "ftab0.dat"; // K-shell Majorant Surface Parameters 87 char nameChar1[100] = "ftab1.dat"; // L-sh 85 char nameChar1[100] = "ftab1.dat"; // L-shell Majorant Surface Parameters 88 86 89 G4String filename; 87 G4String filename; 90 if(level == 0) filename = nameChar0; 88 if(level == 0) filename = nameChar0; 91 if(level == 1) filename = nameChar1; 89 if(level == 1) filename = nameChar1; 92 90 93 const char* path = G4FindDataDir("G4LEDATA << 91 char* path = getenv("G4LEDATA"); 94 if (!path) 92 if (!path) 95 { 93 { 96 G4String excep = "G4EMDataSet - G4LEDA 94 G4String excep = "G4EMDataSet - G4LEDATA environment variable not set"; 97 G4Exception("G4PhotoElectricAngularGen << 95 G4Exception(excep); 98 "em0006",FatalException,"G4LEDATA envi << 99 return; << 100 } 96 } 101 97 102 G4String pathString(path); 98 G4String pathString(path); 103 G4String dirFile = pathString + "/photoele 99 G4String dirFile = pathString + "/photoelectric_angular/" + filename; 104 std::ifstream infile(dirFile); << 100 FILE *infile; 105 if (!infile.is_open()) << 101 infile = fopen(dirFile,"r"); >> 102 if (infile == 0) 106 { 103 { 107 G4String excep = "data file: " + dirFile + " << 104 G4String excep = "G4PhotoElectricAngularGeneratorPolarized - data file: " + dirFile + " not found"; 108 G4Exception("G4PhotoElectricAngularGen << 105 G4Exception(excep); 109 "em0003",FatalException,excep); << 110 return; << 111 } 106 } 112 107 113 // Read parameters into tables. The parame << 108 // Read parameters into tables. The parameters are function of incident electron energy and shell level 114 // energy and shell level << 109 G4float aRead,cRead, beta; 115 G4float aRead=0,cRead=0, beta=0; << 110 for(G4int i=0 ; i<arrayDim ;i++){ 116 for(G4int i=0 ; i<arrayDim ;++i){ << 111 fscanf(infile,"%f\t %e\t %e",&beta,&aRead,&cRead); 117 infile >> beta >> aRead >> cRead; << 112 aMajorantSurfaceParameterTable[i][level] = aRead; 118 aMajorantSurfaceParameterTable[i][level] << 119 cMajorantSurfaceParameterTable[i][level] 113 cMajorantSurfaceParameterTable[i][level] = cRead; 120 } 114 } 121 infile.close(); << 115 fclose(infile); >> 116 122 } 117 } 123 } 118 } 124 119 125 //....oooOO0OOooo........oooOO0OOooo........oo << 120 // 126 121 127 G4PhotoElectricAngularGeneratorPolarized::~G4P << 122 G4PhotoElectricAngularGeneratorPolarized::~G4PhotoElectricAngularGeneratorPolarized() 128 {} << 123 {;} 129 124 130 //....oooOO0OOooo........oooOO0OOooo........oo << 125 // 131 126 132 G4ThreeVector& << 127 G4ThreeVector G4PhotoElectricAngularGeneratorPolarized::GetPhotoElectronDirection(const G4ThreeVector& direction, const G4double eKineticEnergy, 133 G4PhotoElectricAngularGeneratorPolarized::Samp << 128 const G4ThreeVector& polarization, const G4int shellId) const 134 const G4Dynam << 135 G4double eKinEnergy, << 136 G4int shellId, << 137 const G4Material*) << 138 { 129 { 139 // (shellId == 0) - K-shell - Polarized mod << 140 // (shellId > 0) - L1-shell - Polarized mod << 141 << 142 // Calculate Lorentz term (gamma) and beta p 130 // Calculate Lorentz term (gamma) and beta parameters 143 G4double gamma = 1 + eKinEnergy/electron_mas << 131 G4double gamma = 1. + eKineticEnergy/electron_mass_c2; 144 G4double beta = std::sqrt((gamma - 1)*(gamm << 132 G4double beta = std::sqrt(gamma*gamma-1.)/gamma; 145 << 146 const G4ThreeVector& direction = dp->GetMome << 147 const G4ThreeVector& polarization = dp->GetP << 148 133 149 G4double theta, phi = 0; 134 G4double theta, phi = 0; 150 // Majorant surface parameters << 135 G4double aBeta = 0; // Majorant surface parameter (function of the outgoing electron kinetic energy) 151 // function of the outgoing electron kinetic << 136 G4double cBeta = 0; // Majorant surface parameter (function of the outgoing electron kinetic energy) 152 G4double aBeta = 0; << 137 153 G4double cBeta = 0; << 138 G4int shellLevel = 0; 154 << 139 if(shellId < 2) shellLevel = 0; // K-shell // Polarized model for K-shell 155 // For the outgoing kinetic energy << 140 if(shellId >= 2) shellLevel = 1; // L1-shell // Polarized model for L1 and higher shells 156 // find the current majorant surface paramet << 141 157 PhotoElectronGetMajorantSurfaceAandCParamete << 142 // For the outgoing kinetic energy find the current majorant surface parameters 158 << 143 PhotoElectronGetMajorantSurfaceAandCParameters( shellLevel, beta, &aBeta, &cBeta); 159 // Generate pho and theta according to the s << 144 160 // and beta parameter of the electron << 145 // Generate pho and theta according to the shell level and beta parameter of the electron 161 PhotoElectronGeneratePhiAndTheta(shellId, be << 146 PhotoElectronGeneratePhiAndTheta(shellLevel, beta, aBeta, cBeta, &phi, &theta); 162 147 163 // Determine the rotation matrix 148 // Determine the rotation matrix 164 const G4RotationMatrix rotation = << 149 G4RotationMatrix rotation = PhotoElectronRotationMatrix(direction, polarization); 165 PhotoElectronRotationMatrix(direction, pol << 166 150 167 // Compute final direction of the outgoing e 151 // Compute final direction of the outgoing electron 168 fLocalDirection = PhotoElectronComputeFinalD << 152 G4ThreeVector final_direction = PhotoElectronComputeFinalDirection(rotation, theta, phi); 169 153 170 return fLocalDirection; << 154 return final_direction; 171 } 155 } 172 156 173 //....oooOO0OOooo........oooOO0OOooo........oo << 157 // 174 158 175 void << 159 void G4PhotoElectricAngularGeneratorPolarized::PhotoElectronGeneratePhiAndTheta(const G4int shellLevel, const G4double beta, 176 G4PhotoElectricAngularGeneratorPolarized::Phot << 160 const G4double aBeta, const G4double cBeta, 177 G4int shellLevel, G4double beta, G4doubl << 161 G4double *pphi, G4double *ptheta) const 178 G4double *pphi, G4double *ptheta) const << 179 { 162 { 180 G4double rand1, rand2, rand3 = 0; 163 G4double rand1, rand2, rand3 = 0; 181 G4double phi = 0; 164 G4double phi = 0; 182 G4double theta = 0; 165 G4double theta = 0; 183 G4double crossSectionValue = 0; 166 G4double crossSectionValue = 0; 184 G4double crossSectionMajorantFunctionValue = 167 G4double crossSectionMajorantFunctionValue = 0; 185 G4double maxBeta = 0; 168 G4double maxBeta = 0; 186 169 187 do { 170 do { 188 171 189 rand1 = G4UniformRand(); 172 rand1 = G4UniformRand(); 190 rand2 = G4UniformRand(); 173 rand2 = G4UniformRand(); 191 rand3 = G4UniformRand(); 174 rand3 = G4UniformRand(); 192 << 175 193 phi=2*pi*rand1; 176 phi=2*pi*rand1; 194 177 195 if(shellLevel == 0){ 178 if(shellLevel == 0){ >> 179 196 // Polarized Gavrila Cross-Section for K 180 // Polarized Gavrila Cross-Section for K-shell (1959) 197 theta=std::sqrt(((G4Exp(rand2*std::log(1 << 181 theta=std::sqrt(((std::exp(rand2*std::log(1+cBeta*pi*pi)))-1)/cBeta); 198 crossSectionMajorantFunctionValue = << 182 crossSectionMajorantFunctionValue = CrossSectionMajorantFunction(theta, cBeta); 199 CrossSectionMajorantFunction(theta, cBeta); << 200 crossSectionValue = DSigmaKshellGavrila1 183 crossSectionValue = DSigmaKshellGavrila1959(beta, theta, phi); >> 184 201 } else { 185 } else { >> 186 202 // Polarized Gavrila Cross-Section for 187 // Polarized Gavrila Cross-Section for other shells (L1-shell) (1961) 203 theta = std::sqrt(((G4Exp(rand2*std::log << 188 theta = std::sqrt(((std::exp(rand2*std::log(1+cBeta*pi*pi)))-1)/cBeta); 204 crossSectionMajorantFunctionValue = << 189 crossSectionMajorantFunctionValue = CrossSectionMajorantFunction(theta, cBeta); 205 CrossSectionMajorantFunction(theta, cBeta); << 206 crossSectionValue = DSigmaL1shellGavrila 190 crossSectionValue = DSigmaL1shellGavrila(beta, theta, phi); >> 191 207 } 192 } 208 193 209 maxBeta=rand3*aBeta*crossSectionMajorantFu 194 maxBeta=rand3*aBeta*crossSectionMajorantFunctionValue; 210 if(crossSectionValue < 0.0) { crossSection << 211 195 212 } while(maxBeta > crossSectionValue || theta << 196 }while(maxBeta > crossSectionValue); 213 197 214 *pphi = phi; 198 *pphi = phi; 215 *ptheta = theta; 199 *ptheta = theta; 216 } 200 } 217 201 218 //....oooOO0OOooo........oooOO0OOooo........oo << 202 // 219 203 220 G4double << 204 G4double G4PhotoElectricAngularGeneratorPolarized::CrossSectionMajorantFunction(const G4double theta, const G4double cBeta) const 221 G4PhotoElectricAngularGeneratorPolarized::Cros << 222 G4double theta, G4double cBeta) const << 223 { 205 { 224 // Compute Majorant Function 206 // Compute Majorant Function 225 G4double crossSectionMajorantFunctionValue = << 207 G4double crossSectionMajorantFunctionValue = 0; 226 crossSectionMajorantFunctionValue = theta/(1 208 crossSectionMajorantFunctionValue = theta/(1+cBeta*theta*theta); 227 return crossSectionMajorantFunctionValue; 209 return crossSectionMajorantFunctionValue; 228 } 210 } 229 211 230 //....oooOO0OOooo........oooOO0OOooo........oo << 212 // 231 213 232 G4double << 214 G4double G4PhotoElectricAngularGeneratorPolarized::DSigmaKshellGavrila1959(const G4double beta, const G4double theta, const G4double phi) const 233 G4PhotoElectricAngularGeneratorPolarized::DSig << 234 G4double beta, G4double theta, G4doub << 235 { 215 { >> 216 236 //Double differential K shell cross-section 217 //Double differential K shell cross-section (Gavrila 1959) 237 218 238 G4double beta2 = beta*beta; 219 G4double beta2 = beta*beta; 239 G4double oneBeta2 = 1 - beta2; 220 G4double oneBeta2 = 1 - beta2; 240 G4double sqrtOneBeta2 = std::sqrt(oneBeta2); 221 G4double sqrtOneBeta2 = std::sqrt(oneBeta2); 241 G4double oneBeta2_to_3_2 = std::pow(oneBeta2 222 G4double oneBeta2_to_3_2 = std::pow(oneBeta2,1.5); 242 G4double cosTheta = std::cos(theta); 223 G4double cosTheta = std::cos(theta); 243 G4double sinTheta2 = std::sin(theta)*std::si 224 G4double sinTheta2 = std::sin(theta)*std::sin(theta); 244 G4double cosPhi2 = std::cos(phi)*std::cos(ph 225 G4double cosPhi2 = std::cos(phi)*std::cos(phi); 245 G4double oneBetaCosTheta = 1-beta*cosTheta; 226 G4double oneBetaCosTheta = 1-beta*cosTheta; 246 G4double dsigma = 0; 227 G4double dsigma = 0; 247 G4double firstTerm = 0; 228 G4double firstTerm = 0; 248 G4double secondTerm = 0; 229 G4double secondTerm = 0; 249 230 250 firstTerm = sinTheta2*cosPhi2/std::pow(oneBe << 231 firstTerm = sinTheta2*cosPhi2/std::pow(oneBetaCosTheta,4)-(1 - sqrtOneBeta2)/(2*oneBeta2) * 251 (sinTheta2 * cosPhi2)/std::pow(o 232 (sinTheta2 * cosPhi2)/std::pow(oneBetaCosTheta,3) + (1-sqrtOneBeta2)* 252 (1-sqrtOneBeta2)/(4*oneBeta2_to_ 233 (1-sqrtOneBeta2)/(4*oneBeta2_to_3_2) * sinTheta2/std::pow(oneBetaCosTheta,3); 253 234 254 secondTerm = std::sqrt(1 - sqrtOneBeta2)/(st 235 secondTerm = std::sqrt(1 - sqrtOneBeta2)/(std::pow(2.,3.5)*beta2*std::pow(oneBetaCosTheta,2.5)) * 255 (4*beta2/sqrtOneBeta2 * sinThet 236 (4*beta2/sqrtOneBeta2 * sinTheta2*cosPhi2/oneBetaCosTheta + 4*beta/oneBeta2 * cosTheta * cosPhi2 256 - 4*(1-sqrtOneBeta2)/oneBeta2 * 237 - 4*(1-sqrtOneBeta2)/oneBeta2 *(1+cosPhi2) - beta2 * (1-sqrtOneBeta2)/oneBeta2 * sinTheta2/oneBetaCosTheta 257 + 4*beta2*(1-sqrtOneBeta2)/oneB 238 + 4*beta2*(1-sqrtOneBeta2)/oneBeta2_to_3_2 - 4*beta*(1-sqrtOneBeta2)*(1-sqrtOneBeta2)/oneBeta2_to_3_2 * cosTheta) 258 + (1-sqrtOneBeta2)/(4*beta2*one << 239 + (1-sqrtOneBeta2)/(4*beta2*oneBetaCosTheta*oneBetaCosTheta) * (beta/oneBeta2 - 2/oneBeta2 * cosTheta * cosPhi2 + 259 (1-sqrtOneBeta2)/oneBeta2_to_3_ 240 (1-sqrtOneBeta2)/oneBeta2_to_3_2 * cosTheta - beta * (1-sqrtOneBeta2)/oneBeta2_to_3_2); 260 241 261 dsigma = ( firstTerm*(1-pi*fine_structure_co << 242 dsigma = ( firstTerm*(1-pi*fine_structure_const/beta) + secondTerm*(pi*fine_structure_const) ); 262 243 263 return dsigma; 244 return dsigma; 264 } 245 } 265 246 266 //....oooOO0OOooo........oooOO0OOooo........oo << 247 // 267 248 268 G4double << 249 G4double G4PhotoElectricAngularGeneratorPolarized::DSigmaL1shellGavrila(const G4double beta, const G4double theta, const G4double phi) const 269 G4PhotoElectricAngularGeneratorPolarized::DSig << 270 G4double beta, G4double theta, G4doubl << 271 { 250 { >> 251 272 //Double differential L1 shell cross-section 252 //Double differential L1 shell cross-section (Gavrila 1961) 273 // 26Oct2022: included factor (1/8) in dsigm << 253 274 G4double beta2 = beta*beta; 254 G4double beta2 = beta*beta; 275 G4double oneBeta2 = 1-beta2; 255 G4double oneBeta2 = 1-beta2; 276 G4double sqrtOneBeta2 = std::sqrt(oneBeta2); 256 G4double sqrtOneBeta2 = std::sqrt(oneBeta2); 277 G4double oneBeta2_to_3_2=std::pow(oneBeta2,1 257 G4double oneBeta2_to_3_2=std::pow(oneBeta2,1.5); 278 G4double cosTheta = std::cos(theta); 258 G4double cosTheta = std::cos(theta); 279 G4double sinTheta2 =std::sin(theta)*std::sin 259 G4double sinTheta2 =std::sin(theta)*std::sin(theta); 280 G4double cosPhi2 = std::cos(phi)*std::cos(ph 260 G4double cosPhi2 = std::cos(phi)*std::cos(phi); 281 G4double oneBetaCosTheta = 1-beta*cosTheta; 261 G4double oneBetaCosTheta = 1-beta*cosTheta; 282 << 262 283 G4double dsigma = 0; 263 G4double dsigma = 0; 284 G4double firstTerm = 0; 264 G4double firstTerm = 0; 285 G4double secondTerm = 0; 265 G4double secondTerm = 0; 286 266 287 firstTerm = sinTheta2*cosPhi2/std::pow(oneBe 267 firstTerm = sinTheta2*cosPhi2/std::pow(oneBetaCosTheta,4)-(1 - sqrtOneBeta2)/(2*oneBeta2) 288 * (sinTheta2 * cosPhi2)/std::po 268 * (sinTheta2 * cosPhi2)/std::pow(oneBetaCosTheta,3) + (1-sqrtOneBeta2)* 289 (1-sqrtOneBeta2)/(4*oneBeta2_to_ 269 (1-sqrtOneBeta2)/(4*oneBeta2_to_3_2) * sinTheta2/std::pow(oneBetaCosTheta,3); 290 270 291 secondTerm = std::sqrt(1 - sqrtOneBeta2)/(st 271 secondTerm = std::sqrt(1 - sqrtOneBeta2)/(std::pow(2.,3.5)*beta2*std::pow(oneBetaCosTheta,2.5)) * 292 (4*beta2/sqrtOneBeta2 * sinThet 272 (4*beta2/sqrtOneBeta2 * sinTheta2*cosPhi2/oneBetaCosTheta + 4*beta/oneBeta2 * cosTheta * cosPhi2 293 - 4*(1-sqrtOneBeta2)/oneBeta2 * 273 - 4*(1-sqrtOneBeta2)/oneBeta2 *(1+cosPhi2) - beta2 * (1-sqrtOneBeta2)/oneBeta2 * sinTheta2/oneBetaCosTheta 294 + 4*beta2*(1-sqrtOneBeta2)/oneB 274 + 4*beta2*(1-sqrtOneBeta2)/oneBeta2_to_3_2 - 4*beta*(1-sqrtOneBeta2)*(1-sqrtOneBeta2)/oneBeta2_to_3_2 * cosTheta) 295 + (1-sqrtOneBeta2)/(4*beta2*one << 275 + (1-sqrtOneBeta2)/(4*beta2*oneBetaCosTheta*oneBetaCosTheta) * (beta/oneBeta2 - 2/oneBeta2 * cosTheta * cosPhi2 + 296 (1-sqrtOneBeta2)/oneBeta2_to_3_ 276 (1-sqrtOneBeta2)/oneBeta2_to_3_2*cosTheta - beta*(1-sqrtOneBeta2)/oneBeta2_to_3_2); 297 277 298 dsigma = ( firstTerm*(1-pi*fine_structure_co << 278 dsigma = ( firstTerm*(1-pi*fine_structure_const/beta) + secondTerm*(pi*fine_structure_const) ); 299 279 300 return dsigma; 280 return dsigma; 301 } 281 } 302 282 303 //....oooOO0OOooo........oooOO0OOooo........oo << 283 G4double G4PhotoElectricAngularGeneratorPolarized::GetMax(const G4double arg1, const G4double arg2) const >> 284 { >> 285 if (arg1 > arg2) >> 286 return arg1; >> 287 else >> 288 return arg2; >> 289 } 304 290 305 G4RotationMatrix << 291 // 306 G4PhotoElectricAngularGeneratorPolarized::Phot << 292 307 const G4ThreeVector& direction, << 293 G4RotationMatrix G4PhotoElectricAngularGeneratorPolarized::PhotoElectronRotationMatrix(const G4ThreeVector& direction, 308 const G4ThreeVector& polarization) << 294 const G4ThreeVector& polarization) const 309 { 295 { 310 G4double mK = direction.mag(); 296 G4double mK = direction.mag(); 311 G4double mS = polarization.mag(); 297 G4double mS = polarization.mag(); 312 G4ThreeVector polarization2 = polarization; 298 G4ThreeVector polarization2 = polarization; 313 const G4double kTolerance = 1e-6; 299 const G4double kTolerance = 1e-6; 314 300 315 if(!(polarization.isOrthogonal(direction,kTo 301 if(!(polarization.isOrthogonal(direction,kTolerance)) || mS == 0){ 316 G4ThreeVector d0 = direction.unit(); 302 G4ThreeVector d0 = direction.unit(); 317 G4ThreeVector a1 = PerpendicularVector(d0) << 303 G4ThreeVector a1 = SetPerpendicularVector(d0); 318 G4ThreeVector a0 = a1.unit(); << 304 G4ThreeVector a0 = a1.unit(); 319 G4double rand1 = G4UniformRand(); 305 G4double rand1 = G4UniformRand(); 320 G4double angle = twopi*rand1; << 306 G4double angle = twopi*rand1; 321 G4ThreeVector b0 = d0.cross(a0); << 307 G4ThreeVector b0 = d0.cross(a0); 322 G4ThreeVector c; 308 G4ThreeVector c; 323 c.setX(std::cos(angle)*(a0.x())+std::sin(a 309 c.setX(std::cos(angle)*(a0.x())+std::sin(angle)*b0.x()); 324 c.setY(std::cos(angle)*(a0.y())+std::sin(a 310 c.setY(std::cos(angle)*(a0.y())+std::sin(angle)*b0.y()); 325 c.setZ(std::cos(angle)*(a0.z())+std::sin(a 311 c.setZ(std::cos(angle)*(a0.z())+std::sin(angle)*b0.z()); 326 polarization2 = c.unit(); 312 polarization2 = c.unit(); 327 mS = polarization2.mag(); 313 mS = polarization2.mag(); 328 }else 314 }else 329 { 315 { 330 if ( polarization.howOrthogonal(directio 316 if ( polarization.howOrthogonal(direction) != 0) 331 { 317 { 332 polarization2 = polarization << 318 polarization2 = polarization - polarization.dot(direction)/direction.dot(direction) * direction; 333 - polarization.dot(direction)/direction. << 334 } 319 } 335 } 320 } 336 321 337 G4ThreeVector direction2 = direction/mK; 322 G4ThreeVector direction2 = direction/mK; 338 polarization2 = polarization2/mS; 323 polarization2 = polarization2/mS; 339 324 340 G4ThreeVector y = direction2.cross(polarizat 325 G4ThreeVector y = direction2.cross(polarization2); 341 << 326 342 G4RotationMatrix R(polarization2,y,direction 327 G4RotationMatrix R(polarization2,y,direction2); 343 return R; 328 return R; 344 } 329 } 345 330 346 //....oooOO0OOooo........oooOO0OOooo........oo << 331 void G4PhotoElectricAngularGeneratorPolarized::PhotoElectronGetMajorantSurfaceAandCParameters(const G4int shellLevel, const G4double beta,G4double *majorantSurfaceParameterA, G4double *majorantSurfaceParameterC) const 347 << 348 void << 349 G4PhotoElectricAngularGeneratorPolarized::Phot << 350 { 332 { 351 // This member function finds for a given sh << 333 // This member function finds for a given shell and beta value of the outgoing electron the correct Majorant Surface parameters 352 // of the outgoing electron the correct Majo << 334 353 G4double aBeta,cBeta; 335 G4double aBeta,cBeta; 354 G4double bMin,bStep; 336 G4double bMin,bStep; 355 G4int indexMax; 337 G4int indexMax; 356 G4int level = 0; << 338 G4int level = shellLevel; 357 if(shellId > 0) { level = 1; } << 339 if(shellLevel > 1) level = 1; // protection since only K and L1 polarized double differential cross-sections were implemented 358 << 340 359 bMin = betaArray[0]; 341 bMin = betaArray[0]; 360 bStep = betaArray[1]; 342 bStep = betaArray[1]; 361 indexMax = (G4int)betaArray[2]; 343 indexMax = (G4int)betaArray[2]; 362 const G4double kBias = 1e-9; 344 const G4double kBias = 1e-9; 363 345 364 G4int k = (G4int)((beta-bMin+kBias)/bStep); << 346 G4int k = (G4int)((beta-bMin+kBias)/bStep); 365 << 347 366 if(k < 0) 348 if(k < 0) 367 k = 0; 349 k = 0; 368 if(k > indexMax) 350 if(k > indexMax) 369 k = indexMax; << 351 k = indexMax; 370 << 352 371 if(k == 0) << 353 if(k == 0) 372 aBeta = std::max(aMajorantSurfaceParameter << 354 aBeta = GetMax(aMajorantSurfaceParameterTable[k][level],aMajorantSurfaceParameterTable[k+1][level]); 373 else if(k==indexMax) 355 else if(k==indexMax) 374 aBeta = std::max(aMajorantSurfaceParameter << 356 aBeta = GetMax(aMajorantSurfaceParameterTable[k-1][level],aMajorantSurfaceParameterTable[k][level]); 375 else{ 357 else{ 376 aBeta = std::max(aMajorantSurfaceParameter << 358 aBeta = GetMax(aMajorantSurfaceParameterTable[k-1][level],aMajorantSurfaceParameterTable[k][level]); 377 aBeta = std::max(aBeta,aMajorantSurfacePar << 359 aBeta = GetMax(aBeta,aMajorantSurfaceParameterTable[k+1][level]); 378 } << 360 } 379 << 361 380 if(k == 0) 362 if(k == 0) 381 cBeta = std::max(cMajorantSurfaceParameter << 363 cBeta = GetMax(cMajorantSurfaceParameterTable[k][level],cMajorantSurfaceParameterTable[k+1][level]); 382 else if(k == indexMax) 364 else if(k == indexMax) 383 cBeta = std::max(cMajorantSurfaceParameter << 365 cBeta = GetMax(cMajorantSurfaceParameterTable[k-1][level],cMajorantSurfaceParameterTable[k][level]); 384 else{ 366 else{ 385 cBeta = std::max(cMajorantSurfaceParameter << 367 cBeta = GetMax(cMajorantSurfaceParameterTable[k-1][level],cMajorantSurfaceParameterTable[k][level]); 386 cBeta = std::max(cBeta,cMajorantSurfacePar << 368 cBeta = GetMax(cBeta,cMajorantSurfaceParameterTable[k+1][level]); 387 } 369 } 388 370 389 *majorantSurfaceParameterA = aBeta; 371 *majorantSurfaceParameterA = aBeta; 390 *majorantSurfaceParameterC = cBeta; 372 *majorantSurfaceParameterC = cBeta; >> 373 391 } 374 } 392 375 393 //....oooOO0OOooo........oooOO0OOooo........oo << 394 376 395 G4ThreeVector G4PhotoElectricAngularGeneratorP << 377 // >> 378 G4ThreeVector G4PhotoElectricAngularGeneratorPolarized::PhotoElectronComputeFinalDirection(const G4RotationMatrix& rotation, const G4double theta, const G4double phi) const 396 { 379 { 397 //computes the photoelectron momentum unitar << 380 398 G4double sint = std::sin(theta); << 381 //computes the photoelectron momentum unitary vector 399 G4double px = std::cos(phi)*sint; << 382 G4double px = std::cos(phi)*std::sin(theta); 400 G4double py = std::sin(phi)*sint; << 383 G4double py = std::sin(phi)*std::sin(theta); 401 G4double pz = std::cos(theta); 384 G4double pz = std::cos(theta); 402 385 403 G4ThreeVector samplingDirection(px,py,pz); 386 G4ThreeVector samplingDirection(px,py,pz); 404 387 405 G4ThreeVector outgoingDirection = rotation*s 388 G4ThreeVector outgoingDirection = rotation*samplingDirection; 406 return outgoingDirection; 389 return outgoingDirection; 407 } 390 } 408 391 409 //....oooOO0OOooo........oooOO0OOooo........oo << 392 // 410 393 411 void G4PhotoElectricAngularGeneratorPolarized: 394 void G4PhotoElectricAngularGeneratorPolarized::PrintGeneratorInformation() const 412 { 395 { 413 G4cout << "\n" << G4endl; 396 G4cout << "\n" << G4endl; 414 G4cout << "Polarized Photoelectric Angular G 397 G4cout << "Polarized Photoelectric Angular Generator" << G4endl; 415 G4cout << "PhotoElectric Electron Angular Ge 398 G4cout << "PhotoElectric Electron Angular Generator based on the general Gavrila photoelectron angular distribution" << G4endl; 416 G4cout << "Includes polarization effects for 399 G4cout << "Includes polarization effects for K and L1 atomic shells, according to Gavrilla (1959, 1961)." << G4endl; 417 G4cout << "For higher shells the L1 cross-se 400 G4cout << "For higher shells the L1 cross-section is used." << G4endl; 418 G4cout << "(see Physics Reference Manual) \n 401 G4cout << "(see Physics Reference Manual) \n" << G4endl; 419 } << 402 } 420 << 421 //....oooOO0OOooo........oooOO0OOooo........oo << 422 403 423 G4ThreeVector << 404 G4ThreeVector G4PhotoElectricAngularGeneratorPolarized::SetPerpendicularVector(const G4ThreeVector& a) const 424 G4PhotoElectricAngularGeneratorPolarized::Perp << 425 const G4ThreeVector& a) const << 426 { 405 { 427 G4double dx = a.x(); 406 G4double dx = a.x(); 428 G4double dy = a.y(); 407 G4double dy = a.y(); 429 G4double dz = a.z(); 408 G4double dz = a.z(); 430 G4double x = dx < 0.0 ? -dx : dx; 409 G4double x = dx < 0.0 ? -dx : dx; 431 G4double y = dy < 0.0 ? -dy : dy; 410 G4double y = dy < 0.0 ? -dy : dy; 432 G4double z = dz < 0.0 ? -dz : dz; 411 G4double z = dz < 0.0 ? -dz : dz; 433 if (x < y) { 412 if (x < y) { 434 return x < z ? G4ThreeVector(-dy,dx,0) : G 413 return x < z ? G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy); 435 }else{ 414 }else{ 436 return y < z ? G4ThreeVector(dz,0,-dx) : G 415 return y < z ? G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0); 437 } 416 } 438 } 417 } 439 418