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 // | G4LowEPPolarizedComptonModel-- Geant 28 // | G4LowEPPolarizedComptonModel-- Geant4 Monash University | 29 // | polarised low energy Compton scat 29 // | polarised low energy Compton scattering model. | 30 // | J. M. C. Brown, Monash Universit 30 // | J. M. C. Brown, Monash University, Australia | 31 // | 31 // | | 32 // | 32 // | | 33 // ******************************************* 33 // ********************************************************************* 34 // | 34 // | | 35 // | The following is a Geant4 class to simula 35 // | The following is a Geant4 class to simulate the process of | 36 // | bound electron Compton scattering. Genera 36 // | bound electron Compton scattering. General code structure is | 37 // | based on G4LowEnergyCompton.cc and 37 // | based on G4LowEnergyCompton.cc and | 38 // | G4LivermorePolarizedComptonModel.cc. 38 // | G4LivermorePolarizedComptonModel.cc. | 39 // | Algorithms for photon energy, and ejected 39 // | Algorithms for photon energy, and ejected Compton electron | 40 // | direction taken from: 40 // | direction taken from: | 41 // | 41 // | | 42 // | J. M. C. Brown, M. R. Dimmock, J. E. Gill 42 // | J. M. C. Brown, M. R. Dimmock, J. E. Gillam and D. M. Paganin, | 43 // | "A low energy bound atomic electron Compt 43 // | "A low energy bound atomic electron Compton scattering model | 44 // | for Geant4", NIMB, Vol. 338, 77-88, 2014 44 // | for Geant4", NIMB, Vol. 338, 77-88, 2014. | 45 // | 45 // | | 46 // | The author acknowledges the work of the G 46 // | The author acknowledges the work of the Geant4 collaboration | 47 // | in developing the following algorithms th 47 // | in developing the following algorithms that have been employed | 48 // | or adapeted for the present software: 48 // | or adapeted for the present software: | 49 // | 49 // | | 50 // | # sampling of photon scattering angle, 50 // | # sampling of photon scattering angle, | 51 // | # target element selection in composite 51 // | # target element selection in composite materials, | 52 // | # target shell selection in element, 52 // | # target shell selection in element, | 53 // | # and sampling of bound electron momentu 53 // | # and sampling of bound electron momentum from Compton profiles. | 54 // | 54 // | | 55 // ******************************************* 55 // ********************************************************************* 56 // | 56 // | | 57 // | History: 57 // | History: | 58 // | -------- 58 // | -------- | 59 // | 59 // | | 60 // | Jan. 2015 JMCB - 1st Version based 60 // | Jan. 2015 JMCB - 1st Version based on G4LowEPPComptonModel | 61 // | Feb. 2016 JMCB - Geant4 10.2 FPE fi 61 // | Feb. 2016 JMCB - Geant4 10.2 FPE fix for bug 1676 | 62 // | Nov. 2016 JMCB - Polarisation track 62 // | Nov. 2016 JMCB - Polarisation tracking fix in collaboration | 63 // | of Dr. Merlin Reyn 63 // | of Dr. Merlin Reynaard Kole, | 64 // | University of Gene 64 // | University of Geneva | 65 // | 65 // | | 66 // ******************************************* 66 // ********************************************************************* 67 67 68 #include "G4LowEPPolarizedComptonModel.hh" 68 #include "G4LowEPPolarizedComptonModel.hh" 69 #include "G4AutoLock.hh" 69 #include "G4AutoLock.hh" 70 #include "G4PhysicalConstants.hh" 70 #include "G4PhysicalConstants.hh" 71 #include "G4SystemOfUnits.hh" 71 #include "G4SystemOfUnits.hh" 72 72 73 //******************************************** 73 //**************************************************************************** 74 74 75 using namespace std; 75 using namespace std; 76 namespace { G4Mutex LowEPPolarizedComptonModel 76 namespace { G4Mutex LowEPPolarizedComptonModelMutex = G4MUTEX_INITIALIZER; } 77 77 78 78 79 G4PhysicsFreeVector* G4LowEPPolarizedComptonMo 79 G4PhysicsFreeVector* G4LowEPPolarizedComptonModel::data[] = {nullptr}; 80 G4ShellData* G4LowEPPolarizedComptonMode 80 G4ShellData* G4LowEPPolarizedComptonModel::shellData = nullptr; 81 G4DopplerProfile* G4LowEPPolarizedComptonMode 81 G4DopplerProfile* G4LowEPPolarizedComptonModel::profileData = nullptr; 82 82 83 static const G4double ln10 = G4Log(10.); 83 static const G4double ln10 = G4Log(10.); 84 84 85 G4LowEPPolarizedComptonModel::G4LowEPPolarized 85 G4LowEPPolarizedComptonModel::G4LowEPPolarizedComptonModel(const G4ParticleDefinition*, 86 86 const G4String& nam) 87 : G4VEmModel(nam),isInitialised(false) 87 : G4VEmModel(nam),isInitialised(false) 88 { 88 { 89 verboseLevel=1 ; 89 verboseLevel=1 ; 90 // Verbosity scale: 90 // Verbosity scale: 91 // 0 = nothing 91 // 0 = nothing 92 // 1 = warning for energy non-conservation 92 // 1 = warning for energy non-conservation 93 // 2 = details of energy budget 93 // 2 = details of energy budget 94 // 3 = calculation of cross sections, file o 94 // 3 = calculation of cross sections, file openings, sampling of atoms 95 // 4 = entering in methods 95 // 4 = entering in methods 96 96 97 if( verboseLevel>1 ) { 97 if( verboseLevel>1 ) { 98 G4cout << "Low energy photon Compton model 98 G4cout << "Low energy photon Compton model is constructed " << G4endl; 99 } 99 } 100 100 101 //Mark this model as "applicable" for atomic 101 //Mark this model as "applicable" for atomic deexcitation 102 SetDeexcitationFlag(true); 102 SetDeexcitationFlag(true); 103 103 104 fParticleChange = nullptr; 104 fParticleChange = nullptr; 105 fAtomDeexcitation = nullptr; 105 fAtomDeexcitation = nullptr; 106 } 106 } 107 107 108 //******************************************** 108 //**************************************************************************** 109 109 110 G4LowEPPolarizedComptonModel::~G4LowEPPolarize 110 G4LowEPPolarizedComptonModel::~G4LowEPPolarizedComptonModel() 111 { 111 { 112 if(IsMaster()) { 112 if(IsMaster()) { 113 delete shellData; 113 delete shellData; 114 shellData = nullptr; 114 shellData = nullptr; 115 delete profileData; 115 delete profileData; 116 profileData = nullptr; 116 profileData = nullptr; 117 } 117 } 118 } 118 } 119 119 120 //******************************************** 120 //**************************************************************************** 121 121 122 void G4LowEPPolarizedComptonModel::Initialise( 122 void G4LowEPPolarizedComptonModel::Initialise(const G4ParticleDefinition* particle, 123 const 123 const G4DataVector& cuts) 124 { 124 { 125 if (verboseLevel > 1) { 125 if (verboseLevel > 1) { 126 G4cout << "Calling G4LowEPPolarizedCompton 126 G4cout << "Calling G4LowEPPolarizedComptonModel::Initialise()" << G4endl; 127 } 127 } 128 128 129 // Initialise element selector 129 // Initialise element selector 130 if(IsMaster()) { 130 if(IsMaster()) { 131 // Access to elements 131 // Access to elements 132 const char* path = G4FindDataDir("G4LEDATA << 132 char* path = std::getenv("G4LEDATA"); 133 133 134 G4ProductionCutsTable* theCoupleTable = 134 G4ProductionCutsTable* theCoupleTable = 135 G4ProductionCutsTable::GetProductionCuts 135 G4ProductionCutsTable::GetProductionCutsTable(); 136 G4int numOfCouples = (G4int)theCoupleTable << 136 G4int numOfCouples = theCoupleTable->GetTableSize(); 137 137 138 for(G4int i=0; i<numOfCouples; ++i) { 138 for(G4int i=0; i<numOfCouples; ++i) { 139 const G4Material* material = 139 const G4Material* material = 140 theCoupleTable->GetMaterialCutsCouple( 140 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial(); 141 const G4ElementVector* theElementVector 141 const G4ElementVector* theElementVector = material->GetElementVector(); 142 std::size_t nelm = material->GetNumberOf << 142 G4int nelm = material->GetNumberOfElements(); 143 143 144 for (std::size_t j=0; j<nelm; ++j) { << 144 for (G4int j=0; j<nelm; ++j) { 145 G4int Z = G4lrint((*theElementVector)[ 145 G4int Z = G4lrint((*theElementVector)[j]->GetZ()); 146 if(Z < 1) { Z = 1; } 146 if(Z < 1) { Z = 1; } 147 else if(Z > maxZ){ Z = maxZ; } 147 else if(Z > maxZ){ Z = maxZ; } 148 148 149 if( (!data[Z]) ) { ReadData(Z, path); 149 if( (!data[Z]) ) { ReadData(Z, path); } 150 } 150 } 151 } 151 } 152 152 153 // For Doppler broadening 153 // For Doppler broadening 154 if(!shellData) { 154 if(!shellData) { 155 shellData = new G4ShellData(); 155 shellData = new G4ShellData(); 156 shellData->SetOccupancyData(); 156 shellData->SetOccupancyData(); 157 G4String file = "/doppler/shell-doppler" 157 G4String file = "/doppler/shell-doppler"; 158 shellData->LoadData(file); 158 shellData->LoadData(file); 159 } 159 } 160 if(!profileData) { profileData = new G4Dop 160 if(!profileData) { profileData = new G4DopplerProfile(); } 161 161 162 InitialiseElementSelectors(particle, cuts) 162 InitialiseElementSelectors(particle, cuts); 163 } 163 } 164 164 165 if (verboseLevel > 2) { 165 if (verboseLevel > 2) { 166 G4cout << "Loaded cross section files" << 166 G4cout << "Loaded cross section files" << G4endl; 167 } 167 } 168 168 169 if( verboseLevel>1 ) { 169 if( verboseLevel>1 ) { 170 G4cout << "G4LowEPPolarizedComptonModel is 170 G4cout << "G4LowEPPolarizedComptonModel is initialized " << G4endl 171 << "Energy range: " 171 << "Energy range: " 172 << LowEnergyLimit() / eV << " eV - 172 << LowEnergyLimit() / eV << " eV - " 173 << HighEnergyLimit() / GeV << " GeV 173 << HighEnergyLimit() / GeV << " GeV" 174 << G4endl; 174 << G4endl; 175 } 175 } 176 176 177 if(isInitialised) { return; } 177 if(isInitialised) { return; } 178 178 179 fParticleChange = GetParticleChangeForGamma( 179 fParticleChange = GetParticleChangeForGamma(); 180 fAtomDeexcitation = G4LossTableManager::Ins 180 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation(); 181 isInitialised = true; 181 isInitialised = true; 182 } 182 } 183 183 184 //******************************************** 184 //**************************************************************************** 185 185 186 void G4LowEPPolarizedComptonModel::InitialiseL 186 void G4LowEPPolarizedComptonModel::InitialiseLocal(const G4ParticleDefinition*, 187 187 G4VEmModel* masterModel) 188 { 188 { 189 SetElementSelectors(masterModel->GetElementS 189 SetElementSelectors(masterModel->GetElementSelectors()); 190 } 190 } 191 191 192 //******************************************** 192 //**************************************************************************** 193 193 194 void G4LowEPPolarizedComptonModel::ReadData(st << 194 void G4LowEPPolarizedComptonModel::ReadData(size_t Z, const char* path) 195 { 195 { 196 if (verboseLevel > 1) 196 if (verboseLevel > 1) 197 { 197 { 198 G4cout << "G4LowEPPolarizedComptonModel::R 198 G4cout << "G4LowEPPolarizedComptonModel::ReadData()" 199 << G4endl; 199 << G4endl; 200 } 200 } 201 if(data[Z]) { return; } 201 if(data[Z]) { return; } 202 const char* datadir = path; 202 const char* datadir = path; 203 if(!datadir) 203 if(!datadir) 204 { 204 { 205 datadir = G4FindDataDir("G4LEDATA"); << 205 datadir = std::getenv("G4LEDATA"); 206 if(!datadir) 206 if(!datadir) 207 { 207 { 208 G4Exception("G4LowEPPolarizedComptonMode 208 G4Exception("G4LowEPPolarizedComptonModel::ReadData()", 209 "em0006",FatalException, 209 "em0006",FatalException, 210 "Environment variable G4LEDA 210 "Environment variable G4LEDATA not defined"); 211 return; 211 return; 212 } 212 } 213 } 213 } 214 214 215 data[Z] = new G4PhysicsFreeVector(); 215 data[Z] = new G4PhysicsFreeVector(); 216 216 217 std::ostringstream ost; 217 std::ostringstream ost; 218 ost << datadir << "/livermore/comp/ce-cs-" < 218 ost << datadir << "/livermore/comp/ce-cs-" << Z <<".dat"; 219 std::ifstream fin(ost.str().c_str()); 219 std::ifstream fin(ost.str().c_str()); 220 220 221 if( !fin.is_open()) 221 if( !fin.is_open()) 222 { 222 { 223 G4ExceptionDescription ed; 223 G4ExceptionDescription ed; 224 ed << "G4LowEPPolarizedComptonModel data 224 ed << "G4LowEPPolarizedComptonModel data file <" << ost.str().c_str() 225 << "> is not opened!" << G4endl; 225 << "> is not opened!" << G4endl; 226 G4Exception("G4LowEPPolarizedComptonModel: 226 G4Exception("G4LowEPPolarizedComptonModel::ReadData()", 227 "em0003",FatalException, 227 "em0003",FatalException, 228 ed,"G4LEDATA version should be 228 ed,"G4LEDATA version should be G4EMLOW6.34 or later"); 229 return; 229 return; 230 } else { 230 } else { 231 if(verboseLevel > 3) { 231 if(verboseLevel > 3) { 232 G4cout << "File " << ost.str() 232 G4cout << "File " << ost.str() 233 << " is opened by G4LowEPPolari 233 << " is opened by G4LowEPPolarizedComptonModel" << G4endl; 234 } 234 } 235 data[Z]->Retrieve(fin, true); 235 data[Z]->Retrieve(fin, true); 236 data[Z]->ScaleVector(MeV, MeV*barn); 236 data[Z]->ScaleVector(MeV, MeV*barn); 237 } 237 } 238 fin.close(); 238 fin.close(); 239 } 239 } 240 240 241 //******************************************** 241 //**************************************************************************** 242 242 243 G4double 243 G4double 244 G4LowEPPolarizedComptonModel::ComputeCrossSect 244 G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*, 245 245 G4double GammaEnergy, 246 246 G4double Z, G4double, 247 247 G4double, G4double) 248 { 248 { 249 if (verboseLevel > 3) { 249 if (verboseLevel > 3) { 250 G4cout << "G4LowEPPolarizedComptonModel::C 250 G4cout << "G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom()" 251 << G4endl; 251 << G4endl; 252 } 252 } 253 G4double cs = 0.0; 253 G4double cs = 0.0; 254 254 255 if (GammaEnergy < LowEnergyLimit()) { return 255 if (GammaEnergy < LowEnergyLimit()) { return 0.0; } 256 256 257 G4int intZ = G4lrint(Z); 257 G4int intZ = G4lrint(Z); 258 if(intZ < 1 || intZ > maxZ) { return cs; } 258 if(intZ < 1 || intZ > maxZ) { return cs; } 259 259 260 G4PhysicsFreeVector* pv = data[intZ]; 260 G4PhysicsFreeVector* pv = data[intZ]; 261 261 262 // if element was not initialised 262 // if element was not initialised 263 // do initialisation safely for MT mode 263 // do initialisation safely for MT mode 264 if(!pv) 264 if(!pv) 265 { 265 { 266 InitialiseForElement(0, intZ); 266 InitialiseForElement(0, intZ); 267 pv = data[intZ]; 267 pv = data[intZ]; 268 if(!pv) { return cs; } 268 if(!pv) { return cs; } 269 } 269 } 270 270 271 G4int n = G4int(pv->GetVectorLength() - 1); << 271 G4int n = pv->GetVectorLength() - 1; 272 G4double e1 = pv->Energy(0); 272 G4double e1 = pv->Energy(0); 273 G4double e2 = pv->Energy(n); 273 G4double e2 = pv->Energy(n); 274 274 275 if(GammaEnergy <= e1) { cs = GammaEnerg 275 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); } 276 else if(GammaEnergy <= e2) { cs = pv->Value( 276 else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; } 277 else if(GammaEnergy > e2) { cs = pv->Value( 277 else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; } 278 278 279 return cs; 279 return cs; 280 } 280 } 281 281 282 //******************************************** 282 //**************************************************************************** 283 283 284 void G4LowEPPolarizedComptonModel::SampleSecon 284 void G4LowEPPolarizedComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, 285 const G4MaterialCutsCouple* c 285 const G4MaterialCutsCouple* couple, 286 const G4DynamicParticle* aDyn 286 const G4DynamicParticle* aDynamicGamma, 287 G4double, G4double) 287 G4double, G4double) 288 { 288 { 289 289 290 //Determine number of digits (in decimal bas 290 //Determine number of digits (in decimal base) that G4double can accurately represent 291 G4double g4d_order = G4double(numeric_limits 291 G4double g4d_order = G4double(numeric_limits<G4double>::digits10); 292 G4double g4d_limit = std::pow(10.,-g4d_order 292 G4double g4d_limit = std::pow(10.,-g4d_order); 293 293 294 // The scattered gamma energy is sampled acc 294 // The scattered gamma energy is sampled according to Klein - Nishina formula. 295 // then accepted or rejected depending on th 295 // then accepted or rejected depending on the Scattering Function multiplied 296 // by factor from Klein - Nishina formula. 296 // by factor from Klein - Nishina formula. 297 // Expression of the angular distribution as 297 // Expression of the angular distribution as Klein Nishina 298 // angular and energy distribution and Scatt 298 // angular and energy distribution and Scattering fuctions is taken from 299 // D. E. Cullen "A simple model of photon tr 299 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth. 300 // Phys. Res. B 101 (1995). Method of sampli 300 // Phys. Res. B 101 (1995). Method of sampling with form factors is different 301 // data are interpolated while in the articl 301 // data are interpolated while in the article they are fitted. 302 // Reference to the article is from J. Stepa 302 // Reference to the article is from J. Stepanek New Photon, Positron 303 // and Electron Interaction Data for GEANT i 303 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10 304 // TeV (draft). 304 // TeV (draft). 305 // The random number techniques of Butcher & 305 // The random number techniques of Butcher & Messel are used 306 // (Nucl Phys 20(1960),15). 306 // (Nucl Phys 20(1960),15). 307 307 308 G4double photonEnergy0 = aDynamicGamma->GetK 308 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy()/MeV; 309 309 310 if (verboseLevel > 3) { 310 if (verboseLevel > 3) { 311 G4cout << "G4LowEPPolarizedComptonModel::S 311 G4cout << "G4LowEPPolarizedComptonModel::SampleSecondaries() E(MeV)= " 312 << photonEnergy0/MeV << " in " << c 312 << photonEnergy0/MeV << " in " << couple->GetMaterial()->GetName() 313 << G4endl; 313 << G4endl; 314 } 314 } 315 // do nothing below the threshold 315 // do nothing below the threshold 316 // should never get here because the XS is z 316 // should never get here because the XS is zero below the limit 317 if (photonEnergy0 < LowEnergyLimit()) 317 if (photonEnergy0 < LowEnergyLimit()) 318 return ; 318 return ; 319 319 320 G4double e0m = photonEnergy0 / electron_mass 320 G4double e0m = photonEnergy0 / electron_mass_c2 ; 321 G4ParticleMomentum photonDirection0 = aDynam 321 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); 322 322 323 // Polarisation: check orientation of photon 323 // Polarisation: check orientation of photon propagation direction and polarisation 324 // Fix if needed 324 // Fix if needed 325 325 326 G4ThreeVector photonPolarization0 = aDynamic 326 G4ThreeVector photonPolarization0 = aDynamicGamma->GetPolarization(); 327 // Check if polarisation vector is perpendic 327 // Check if polarisation vector is perpendicular and fix if not 328 328 329 if (!(photonPolarization0.isOrthogonal(photo 329 if (!(photonPolarization0.isOrthogonal(photonDirection0, 1e-6))||(photonPolarization0.mag()==0)) 330 { 330 { 331 photonPolarization0 = GetRandomPolarizat 331 photonPolarization0 = GetRandomPolarization(photonDirection0); 332 } 332 } 333 else 333 else 334 { 334 { 335 if ((photonPolarization0.howOrthogonal(p 335 if ((photonPolarization0.howOrthogonal(photonDirection0) !=0) && (photonPolarization0.howOrthogonal(photonDirection0) > g4d_limit)) 336 { 336 { 337 photonPolarization0 = GetPerpendicul 337 photonPolarization0 = GetPerpendicularPolarization(photonDirection0,photonPolarization0); 338 } 338 } 339 } 339 } 340 340 341 // Select randomly one element in the curren 341 // Select randomly one element in the current material 342 const G4ParticleDefinition* particle = aDyn 342 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition(); 343 const G4Element* elm = SelectRandomAtom(coup 343 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0); 344 G4int Z = (G4int)elm->GetZ(); 344 G4int Z = (G4int)elm->GetZ(); 345 345 346 G4double LowEPPCepsilon0 = 1. / (1. + 2. * e 346 G4double LowEPPCepsilon0 = 1. / (1. + 2. * e0m); 347 G4double LowEPPCepsilon0Sq = LowEPPCepsilon0 347 G4double LowEPPCepsilon0Sq = LowEPPCepsilon0 * LowEPPCepsilon0; 348 G4double alpha1 = -std::log(LowEPPCepsilon0) 348 G4double alpha1 = -std::log(LowEPPCepsilon0); 349 G4double alpha2 = 0.5 * (1. - LowEPPCepsilon 349 G4double alpha2 = 0.5 * (1. - LowEPPCepsilon0Sq); 350 350 351 G4double wlPhoton = h_Planck*c_light/photonE 351 G4double wlPhoton = h_Planck*c_light/photonEnergy0; 352 352 353 // Sample the energy of the scattered photon 353 // Sample the energy of the scattered photon 354 G4double LowEPPCepsilon; 354 G4double LowEPPCepsilon; 355 G4double LowEPPCepsilonSq; 355 G4double LowEPPCepsilonSq; 356 G4double oneCosT; 356 G4double oneCosT; 357 G4double sinT2; 357 G4double sinT2; 358 G4double gReject; 358 G4double gReject; 359 359 360 if (verboseLevel > 3) { 360 if (verboseLevel > 3) { 361 G4cout << "Started loop to sample gamma en 361 G4cout << "Started loop to sample gamma energy" << G4endl; 362 } 362 } 363 363 364 do 364 do 365 { 365 { 366 if ( alpha1/(alpha1+alpha2) > G4UniformR 366 if ( alpha1/(alpha1+alpha2) > G4UniformRand()) 367 { 367 { 368 LowEPPCepsilon = G4Exp(-alpha1 * G4U 368 LowEPPCepsilon = G4Exp(-alpha1 * G4UniformRand()); 369 LowEPPCepsilonSq = LowEPPCepsilon * 369 LowEPPCepsilonSq = LowEPPCepsilon * LowEPPCepsilon; 370 } 370 } 371 else 371 else 372 { 372 { 373 LowEPPCepsilonSq = LowEPPCepsilon0Sq 373 LowEPPCepsilonSq = LowEPPCepsilon0Sq + (1. - LowEPPCepsilon0Sq) * G4UniformRand(); 374 LowEPPCepsilon = std::sqrt(LowEPPCep 374 LowEPPCepsilon = std::sqrt(LowEPPCepsilonSq); 375 } 375 } 376 376 377 oneCosT = (1. - LowEPPCepsilon) / ( LowE 377 oneCosT = (1. - LowEPPCepsilon) / ( LowEPPCepsilon * e0m); 378 sinT2 = oneCosT * (2. - oneCosT); 378 sinT2 = oneCosT * (2. - oneCosT); 379 G4double x = std::sqrt(oneCosT/2.) / (wl 379 G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm); 380 G4double scatteringFunction = ComputeSca 380 G4double scatteringFunction = ComputeScatteringFunction(x, Z); 381 gReject = (1. - LowEPPCepsilon * sinT2 / 381 gReject = (1. - LowEPPCepsilon * sinT2 / (1. + LowEPPCepsilonSq)) * scatteringFunction; 382 382 383 } while(gReject < G4UniformRand()*Z); 383 } while(gReject < G4UniformRand()*Z); 384 384 385 G4double cosTheta = 1. - oneCosT; 385 G4double cosTheta = 1. - oneCosT; 386 G4double sinTheta = std::sqrt(sinT2); 386 G4double sinTheta = std::sqrt(sinT2); 387 G4double phi = SetPhi(LowEPPCepsilon,sinT2); 387 G4double phi = SetPhi(LowEPPCepsilon,sinT2); 388 G4double dirx = sinTheta * std::cos(phi); 388 G4double dirx = sinTheta * std::cos(phi); 389 G4double diry = sinTheta * std::sin(phi); 389 G4double diry = sinTheta * std::sin(phi); 390 G4double dirz = cosTheta ; 390 G4double dirz = cosTheta ; 391 391 392 // Set outgoing photon polarization 392 // Set outgoing photon polarization 393 393 394 G4ThreeVector photonPolarization1 = SetNewPo 394 G4ThreeVector photonPolarization1 = SetNewPolarization(LowEPPCepsilon, 395 sinT2, 395 sinT2, 396 phi, 396 phi, 397 cosTheta); 397 cosTheta); 398 398 399 // Scatter photon energy and Compton electro 399 // Scatter photon energy and Compton electron direction - Method based on: 400 // J. M. C. Brown, M. R. Dimmock, J. E. Gill 400 // J. M. C. Brown, M. R. Dimmock, J. E. Gillam and D. M. Paganin' 401 // "A low energy bound atomic electron Compt 401 // "A low energy bound atomic electron Compton scattering model for Geant4" 402 // NIMB, Vol. 338, 77-88, 2014. 402 // NIMB, Vol. 338, 77-88, 2014. 403 403 404 // Set constants and initialize scattering p 404 // Set constants and initialize scattering parameters 405 const G4double vel_c = c_light / (m/s); 405 const G4double vel_c = c_light / (m/s); 406 const G4double momentum_au_to_nat = halfpi* 406 const G4double momentum_au_to_nat = halfpi* hbar_Planck / Bohr_radius / (kg*m/s); 407 const G4double e_mass_kg = electron_mass_c2 407 const G4double e_mass_kg = electron_mass_c2 / c_squared / kg ; 408 408 409 const G4int maxDopplerIterations = 1000; 409 const G4int maxDopplerIterations = 1000; 410 G4double bindingE = 0.; 410 G4double bindingE = 0.; 411 G4double pEIncident = photonEnergy0 ; 411 G4double pEIncident = photonEnergy0 ; 412 G4double pERecoil = -1.; 412 G4double pERecoil = -1.; 413 G4double eERecoil = -1.; 413 G4double eERecoil = -1.; 414 G4double e_alpha =0.; 414 G4double e_alpha =0.; 415 G4double e_beta = 0.; 415 G4double e_beta = 0.; 416 416 417 G4double CE_emission_flag = 0.; 417 G4double CE_emission_flag = 0.; 418 G4double ePAU = -1; 418 G4double ePAU = -1; 419 G4int shellIdx = 0; 419 G4int shellIdx = 0; 420 G4double u_temp = 0; 420 G4double u_temp = 0; 421 G4double cosPhiE =0; 421 G4double cosPhiE =0; 422 G4double sinThetaE =0; 422 G4double sinThetaE =0; 423 G4double cosThetaE =0; 423 G4double cosThetaE =0; 424 G4int iteration = 0; 424 G4int iteration = 0; 425 425 426 if (verboseLevel > 3) { 426 if (verboseLevel > 3) { 427 G4cout << "Started loop to sample photon e 427 G4cout << "Started loop to sample photon energy and electron direction" << G4endl; 428 } 428 } 429 429 430 do{ 430 do{ 431 // *************************************** 431 // ****************************************** 432 // | Determine scatter photon energy 432 // | Determine scatter photon energy | 433 // *************************************** 433 // ****************************************** 434 do 434 do 435 { 435 { 436 iteration++; 436 iteration++; 437 437 438 // ***************************************** 438 // ******************************************** 439 // | Sample bound electron information 439 // | Sample bound electron information | 440 // ***************************************** 440 // ******************************************** 441 441 442 // Select shell based on shell occupancy 442 // Select shell based on shell occupancy 443 shellIdx = shellData->SelectRandomShell(Z); 443 shellIdx = shellData->SelectRandomShell(Z); 444 bindingE = shellData->BindingEnergy(Z,shellI 444 bindingE = shellData->BindingEnergy(Z,shellIdx)/MeV; 445 445 446 // Randomly sample bound electron momentum ( 446 // Randomly sample bound electron momentum (memento: the data set is in Atomic Units) 447 ePAU = profileData->RandomSelectMomentum(Z,s 447 ePAU = profileData->RandomSelectMomentum(Z,shellIdx); 448 448 449 // Convert to SI units 449 // Convert to SI units 450 G4double ePSI = ePAU * momentum_au_to_nat; 450 G4double ePSI = ePAU * momentum_au_to_nat; 451 451 452 //Calculate bound electron velocity and norm 452 //Calculate bound electron velocity and normalise to natural units 453 u_temp = sqrt( ((ePSI*ePSI)*(vel_c*vel_c)) / 453 u_temp = sqrt( ((ePSI*ePSI)*(vel_c*vel_c)) / ((e_mass_kg*e_mass_kg)*(vel_c*vel_c)+(ePSI*ePSI)) )/vel_c; 454 454 455 // Sample incident electron direction, amorp 455 // Sample incident electron direction, amorphous material, to scattering photon scattering plane 456 e_alpha = pi*G4UniformRand(); 456 e_alpha = pi*G4UniformRand(); 457 e_beta = twopi*G4UniformRand(); 457 e_beta = twopi*G4UniformRand(); 458 458 459 // Total energy of system 459 // Total energy of system 460 G4double eEIncident = electron_mass_c2 / sqr 460 G4double eEIncident = electron_mass_c2 / sqrt( 1 - (u_temp*u_temp)); 461 G4double systemE = eEIncident + pEIncident; 461 G4double systemE = eEIncident + pEIncident; 462 462 463 G4double gamma_temp = 1.0 / sqrt( 1 - (u_tem 463 G4double gamma_temp = 1.0 / sqrt( 1 - (u_temp*u_temp)); 464 G4double numerator = gamma_temp*electron_mas 464 G4double numerator = gamma_temp*electron_mass_c2*(1 - u_temp * std::cos(e_alpha)); 465 G4double subdenom1 = u_temp*cosTheta*std::c 465 G4double subdenom1 = u_temp*cosTheta*std::cos(e_alpha); 466 G4double subdenom2 = u_temp*sinTheta*std::si 466 G4double subdenom2 = u_temp*sinTheta*std::sin(e_alpha)*std::cos(e_beta); 467 G4double denominator = (1.0 - cosTheta) + ( 467 G4double denominator = (1.0 - cosTheta) + (gamma_temp*electron_mass_c2*(1 - subdenom1 - subdenom2) / pEIncident); 468 pERecoil = (numerator/denominator); 468 pERecoil = (numerator/denominator); 469 eERecoil = systemE - pERecoil; 469 eERecoil = systemE - pERecoil; 470 CE_emission_flag = pEIncident - pERecoil; 470 CE_emission_flag = pEIncident - pERecoil; 471 } while ( (iteration <= maxDopplerIterat 471 } while ( (iteration <= maxDopplerIterations) && (CE_emission_flag < bindingE)); 472 472 473 // End of recalculation of photon energy w 473 // End of recalculation of photon energy with Doppler broadening 474 // *************************************** 474 // ******************************************************* 475 // | Determine ejected Compton electro 475 // | Determine ejected Compton electron direction | 476 // *************************************** 476 // ******************************************************* 477 477 478 // Calculate velocity of ejected Compton e 478 // Calculate velocity of ejected Compton electron 479 479 480 G4double a_temp = eERecoil / electron_mass 480 G4double a_temp = eERecoil / electron_mass_c2; 481 G4double u_p_temp = sqrt(1 - (1 / (a_temp* 481 G4double u_p_temp = sqrt(1 - (1 / (a_temp*a_temp))); 482 482 483 // Coefficients and terms from simulatenou 483 // Coefficients and terms from simulatenous equations 484 G4double sinAlpha = std::sin(e_alpha); 484 G4double sinAlpha = std::sin(e_alpha); 485 G4double cosAlpha = std::cos(e_alpha); 485 G4double cosAlpha = std::cos(e_alpha); 486 G4double sinBeta = std::sin(e_beta); 486 G4double sinBeta = std::sin(e_beta); 487 G4double cosBeta = std::cos(e_beta); 487 G4double cosBeta = std::cos(e_beta); 488 488 489 G4double gamma = 1.0 / sqrt(1 - (u_temp*u_ 489 G4double gamma = 1.0 / sqrt(1 - (u_temp*u_temp)); 490 G4double gamma_p = 1.0 / sqrt(1 - (u_p_tem 490 G4double gamma_p = 1.0 / sqrt(1 - (u_p_temp*u_p_temp)); 491 491 492 G4double var_A = pERecoil*u_p_temp*sinThet 492 G4double var_A = pERecoil*u_p_temp*sinTheta; 493 G4double var_B = u_p_temp* (pERecoil*cosTh 493 G4double var_B = u_p_temp* (pERecoil*cosTheta-pEIncident); 494 G4double var_C = (pERecoil-pEIncident) - ( 494 G4double var_C = (pERecoil-pEIncident) - ( (pERecoil*pEIncident) / (gamma_p*electron_mass_c2))*(1 - cosTheta); 495 495 496 G4double var_D1 = gamma*electron_mass_c2*p 496 G4double var_D1 = gamma*electron_mass_c2*pERecoil; 497 G4double var_D2 = (1 - (u_temp*cosTheta*co 497 G4double var_D2 = (1 - (u_temp*cosTheta*cosAlpha) - (u_temp*sinTheta*cosBeta*sinAlpha)); 498 G4double var_D3 = ((electron_mass_c2*elect 498 G4double var_D3 = ((electron_mass_c2*electron_mass_c2)*(gamma*gamma_p - 1)) - (gamma_p*electron_mass_c2*pERecoil); 499 G4double var_D = var_D1*var_D2 + var_D3; 499 G4double var_D = var_D1*var_D2 + var_D3; 500 500 501 G4double var_E1 = ((gamma*gamma_p)*(electr 501 G4double var_E1 = ((gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*cosAlpha); 502 G4double var_E2 = gamma_p*electron_mass_c2 502 G4double var_E2 = gamma_p*electron_mass_c2*pERecoil*u_p_temp*cosTheta; 503 G4double var_E = var_E1 - var_E2; 503 G4double var_E = var_E1 - var_E2; 504 504 505 G4double var_F1 = ((gamma*gamma_p)*(electr 505 G4double var_F1 = ((gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*cosBeta*sinAlpha); 506 G4double var_F2 = (gamma_p*electron_mass_c 506 G4double var_F2 = (gamma_p*electron_mass_c2*pERecoil*u_p_temp*sinTheta); 507 G4double var_F = var_F1 - var_F2; 507 G4double var_F = var_F1 - var_F2; 508 508 509 G4double var_G = (gamma*gamma_p)*(electron 509 G4double var_G = (gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*sinBeta*sinAlpha; 510 510 511 // Two equations form a quadratic form of 511 // Two equations form a quadratic form of Wx^2 + Yx + Z = 0 512 // Coefficents and solution to quadratic 512 // Coefficents and solution to quadratic 513 G4double var_W1 = (var_F*var_B - var_E*var 513 G4double var_W1 = (var_F*var_B - var_E*var_A)*(var_F*var_B - var_E*var_A); 514 G4double var_W2 = (var_G*var_G)*(var_A*var 514 G4double var_W2 = (var_G*var_G)*(var_A*var_A) + (var_G*var_G)*(var_B*var_B); 515 G4double var_W = var_W1 + var_W2; 515 G4double var_W = var_W1 + var_W2; 516 516 517 G4double var_Y = 2.0*(((var_A*var_D-var_F* 517 G4double var_Y = 2.0*(((var_A*var_D-var_F*var_C)*(var_F*var_B-var_E*var_A)) - ((var_G*var_G)*var_B*var_C)); 518 518 519 G4double var_Z1 = (var_A*var_D - var_F*var 519 G4double var_Z1 = (var_A*var_D - var_F*var_C)*(var_A*var_D - var_F*var_C); 520 G4double var_Z2 = (var_G*var_G)*(var_C*var 520 G4double var_Z2 = (var_G*var_G)*(var_C*var_C) - (var_G*var_G)*(var_A*var_A); 521 G4double var_Z = var_Z1 + var_Z2; 521 G4double var_Z = var_Z1 + var_Z2; 522 G4double diff1 = var_Y*var_Y; 522 G4double diff1 = var_Y*var_Y; 523 G4double diff2 = 4*var_W*var_Z; 523 G4double diff2 = 4*var_W*var_Z; 524 G4double diff = diff1 - diff2; 524 G4double diff = diff1 - diff2; 525 525 526 // Check if diff is less than zero, if so 526 // Check if diff is less than zero, if so ensure it is due to FPE 527 //Confirm that diff less than zero is due 527 //Confirm that diff less than zero is due FPE, i.e if abs of diff / diff1 and diff/ diff2 is less 528 //than 10^(-g4d_order), then set diff to z 528 //than 10^(-g4d_order), then set diff to zero 529 529 530 if ((diff < 0.0) && (abs(diff / diff1) < g 530 if ((diff < 0.0) && (abs(diff / diff1) < g4d_limit) && (abs(diff / diff2) < g4d_limit) ) 531 { 531 { 532 diff = 0.0; 532 diff = 0.0; 533 } 533 } 534 534 535 // Plus and minus of quadratic 535 // Plus and minus of quadratic 536 G4double X_p = (-var_Y + sqrt (diff))/(2*v 536 G4double X_p = (-var_Y + sqrt (diff))/(2*var_W); 537 G4double X_m = (-var_Y - sqrt (diff))/(2*v 537 G4double X_m = (-var_Y - sqrt (diff))/(2*var_W); 538 538 539 // Floating point precision protection 539 // Floating point precision protection 540 // Check if X_p and X_m are greater than o 540 // Check if X_p and X_m are greater than or less than 1 or -1, if so clean up FPE 541 // Issue due to propagation of FPE and onl 541 // Issue due to propagation of FPE and only impacts 8th sig fig onwards 542 if(X_p >1){X_p=1;} if(X_p<-1){X_p=-1;} 542 if(X_p >1){X_p=1;} if(X_p<-1){X_p=-1;} 543 if(X_m >1){X_m=1;} if(X_m<-1){X_m=-1;} 543 if(X_m >1){X_m=1;} if(X_m<-1){X_m=-1;} 544 // End of FP protection 544 // End of FP protection 545 545 546 G4double ThetaE = 0.; 546 G4double ThetaE = 0.; 547 547 548 // Randomly sample one of the two possible 548 // Randomly sample one of the two possible solutions and determin theta angle of ejected Compton electron 549 G4double sol_select = G4UniformRand(); 549 G4double sol_select = G4UniformRand(); 550 550 551 if (sol_select < 0.5) 551 if (sol_select < 0.5) 552 { 552 { 553 ThetaE = std::acos(X_p); 553 ThetaE = std::acos(X_p); 554 } 554 } 555 if (sol_select > 0.5) 555 if (sol_select > 0.5) 556 { 556 { 557 ThetaE = std::acos(X_m); 557 ThetaE = std::acos(X_m); 558 } 558 } 559 559 560 cosThetaE = std::cos(ThetaE); 560 cosThetaE = std::cos(ThetaE); 561 sinThetaE = std::sin(ThetaE); 561 sinThetaE = std::sin(ThetaE); 562 G4double Theta = std::acos(cosTheta); 562 G4double Theta = std::acos(cosTheta); 563 563 564 //Calculate electron Phi 564 //Calculate electron Phi 565 G4double iSinThetaE = std::sqrt(1+std::tan 565 G4double iSinThetaE = std::sqrt(1+std::tan((pi/2.0)-ThetaE)*std::tan((pi/2.0)-ThetaE)); 566 G4double iSinTheta = std::sqrt(1+std::tan( 566 G4double iSinTheta = std::sqrt(1+std::tan((pi/2.0)-Theta)*std::tan((pi/2.0)-Theta)); 567 G4double ivar_A = iSinTheta/ (pERecoil*u_p 567 G4double ivar_A = iSinTheta/ (pERecoil*u_p_temp); 568 // Trigs 568 // Trigs 569 cosPhiE = (var_C - var_B*cosThetaE)*(ivar_ 569 cosPhiE = (var_C - var_B*cosThetaE)*(ivar_A*iSinThetaE); 570 // End of calculation of ejection Compton 570 // End of calculation of ejection Compton electron direction 571 //Fix for floating point errors 571 //Fix for floating point errors 572 } while ( (iteration <= maxDopplerIterations 572 } while ( (iteration <= maxDopplerIterations) && (abs(cosPhiE) > 1)); 573 573 574 // Revert to original if maximum number of i 574 // Revert to original if maximum number of iterations threshold has been reached 575 if (iteration >= maxDopplerIterations) 575 if (iteration >= maxDopplerIterations) 576 { 576 { 577 pERecoil = photonEnergy0 ; 577 pERecoil = photonEnergy0 ; 578 bindingE = 0.; 578 bindingE = 0.; 579 dirx=0.0; 579 dirx=0.0; 580 diry=0.0; 580 diry=0.0; 581 dirz=1.0; 581 dirz=1.0; 582 } 582 } 583 583 584 // Set "scattered" photon direction and ener 584 // Set "scattered" photon direction and energy 585 G4ThreeVector photonDirection1(dirx,diry,dir 585 G4ThreeVector photonDirection1(dirx,diry,dirz); 586 SystemOfRefChange(photonDirection0,photonDir 586 SystemOfRefChange(photonDirection0,photonDirection1, 587 photonPolarization0,photonPolarization 587 photonPolarization0,photonPolarization1); 588 588 589 if (pERecoil > 0.) 589 if (pERecoil > 0.) 590 { 590 { 591 fParticleChange->SetProposedKineticEnerg 591 fParticleChange->SetProposedKineticEnergy(pERecoil) ; 592 fParticleChange->ProposeMomentumDirectio 592 fParticleChange->ProposeMomentumDirection(photonDirection1) ; 593 fParticleChange->ProposePolarization(pho 593 fParticleChange->ProposePolarization(photonPolarization1); 594 594 595 // Set ejected Compton electron directio 595 // Set ejected Compton electron direction and energy 596 G4double PhiE = std::acos(cosPhiE); 596 G4double PhiE = std::acos(cosPhiE); 597 G4double eDirX = sinThetaE * std::cos(ph 597 G4double eDirX = sinThetaE * std::cos(phi+PhiE); 598 G4double eDirY = sinThetaE * std::sin(ph 598 G4double eDirY = sinThetaE * std::sin(phi+PhiE); 599 G4double eDirZ = cosThetaE; 599 G4double eDirZ = cosThetaE; 600 600 601 G4double eKineticEnergy = pEIncident - p 601 G4double eKineticEnergy = pEIncident - pERecoil - bindingE; 602 602 603 G4ThreeVector eDirection(eDirX,eDirY,eDi 603 G4ThreeVector eDirection(eDirX,eDirY,eDirZ); 604 SystemOfRefChangeElect(photonDirection0, 604 SystemOfRefChangeElect(photonDirection0,eDirection, 605 photonPolarization0); 605 photonPolarization0); 606 606 607 G4DynamicParticle* dp = new G4DynamicPar 607 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(), 608 eDirection,eKineticEnergy) ; 608 eDirection,eKineticEnergy) ; 609 fvect->push_back(dp); 609 fvect->push_back(dp); 610 } 610 } 611 else 611 else 612 { 612 { 613 fParticleChange->SetProposedKineticEnerg 613 fParticleChange->SetProposedKineticEnergy(0.); 614 fParticleChange->ProposeTrackStatus(fSto 614 fParticleChange->ProposeTrackStatus(fStopAndKill); 615 } 615 } 616 616 617 // sample deexcitation 617 // sample deexcitation 618 // 618 // 619 if (verboseLevel > 3) { 619 if (verboseLevel > 3) { 620 G4cout << "Started atomic de-excitation " 620 G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl; 621 } 621 } 622 622 623 if(fAtomDeexcitation && iteration < maxDoppl 623 if(fAtomDeexcitation && iteration < maxDopplerIterations) { 624 G4int index = couple->GetIndex(); 624 G4int index = couple->GetIndex(); 625 if(fAtomDeexcitation->CheckDeexcitationAct 625 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) { 626 std::size_t nbefore = fvect->size(); << 626 size_t nbefore = fvect->size(); 627 G4AtomicShellEnumerator as = G4AtomicShe 627 G4AtomicShellEnumerator as = G4AtomicShellEnumerator(shellIdx); 628 const G4AtomicShell* shell = fAtomDeexci 628 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as); 629 fAtomDeexcitation->GenerateParticles(fve 629 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index); 630 std::size_t nafter = fvect->size(); << 630 size_t nafter = fvect->size(); 631 if(nafter > nbefore) { 631 if(nafter > nbefore) { 632 for (std::size_t i=nbefore; i<nafter; << 632 for (size_t i=nbefore; i<nafter; ++i) { 633 //Check if there is enough residual 633 //Check if there is enough residual energy 634 if (bindingE >= ((*fvect)[i])->GetKi 634 if (bindingE >= ((*fvect)[i])->GetKineticEnergy()) 635 { 635 { 636 //Ok, this is a valid secondary: keep 636 //Ok, this is a valid secondary: keep it 637 bindingE -= ((*fvect)[i])->GetKineticE 637 bindingE -= ((*fvect)[i])->GetKineticEnergy(); 638 } 638 } 639 else 639 else 640 { 640 { 641 //Invalid secondary: not enough energy 641 //Invalid secondary: not enough energy to create it! 642 //Keep its energy in the local deposit 642 //Keep its energy in the local deposit 643 delete (*fvect)[i]; 643 delete (*fvect)[i]; 644 (*fvect)[i]=nullptr; 644 (*fvect)[i]=nullptr; 645 } 645 } 646 } 646 } 647 } 647 } 648 } 648 } 649 } 649 } 650 650 651 //This should never happen 651 //This should never happen 652 if(bindingE < 0.0) 652 if(bindingE < 0.0) 653 G4Exception("G4LowEPPolarizedComptonModel: 653 G4Exception("G4LowEPPolarizedComptonModel::SampleSecondaries()", 654 "em2051",FatalException,"Negative local en 654 "em2051",FatalException,"Negative local energy deposit"); 655 655 656 fParticleChange->ProposeLocalEnergyDeposit(b 656 fParticleChange->ProposeLocalEnergyDeposit(bindingE); 657 } 657 } 658 658 659 //******************************************** 659 //**************************************************************************** 660 660 661 G4double 661 G4double 662 G4LowEPPolarizedComptonModel::ComputeScatterin 662 G4LowEPPolarizedComptonModel::ComputeScatteringFunction(G4double x, G4int Z) 663 { 663 { 664 G4double value = Z; 664 G4double value = Z; 665 if (x <= ScatFuncFitParam[Z][2]) { 665 if (x <= ScatFuncFitParam[Z][2]) { 666 666 667 G4double lgq = G4Log(x)/ln10; 667 G4double lgq = G4Log(x)/ln10; 668 668 669 if (lgq < ScatFuncFitParam[Z][1]) { 669 if (lgq < ScatFuncFitParam[Z][1]) { 670 value = ScatFuncFitParam[Z][3] + lgq*Sca 670 value = ScatFuncFitParam[Z][3] + lgq*ScatFuncFitParam[Z][4]; 671 } else { 671 } else { 672 value = ScatFuncFitParam[Z][5] + lgq*Sca 672 value = ScatFuncFitParam[Z][5] + lgq*ScatFuncFitParam[Z][6] + 673 lgq*lgq*ScatFuncFitParam[Z][7] + lgq*l 673 lgq*lgq*ScatFuncFitParam[Z][7] + lgq*lgq*lgq*ScatFuncFitParam[Z][8]; 674 } 674 } 675 value = G4Exp(value*ln10); 675 value = G4Exp(value*ln10); 676 } 676 } 677 return value; 677 return value; 678 } 678 } 679 679 680 680 681 //******************************************** 681 //**************************************************************************** 682 682 683 void 683 void 684 G4LowEPPolarizedComptonModel::InitialiseForEle 684 G4LowEPPolarizedComptonModel::InitialiseForElement(const G4ParticleDefinition*, 685 G4int Z) 685 G4int Z) 686 { 686 { 687 G4AutoLock l(&LowEPPolarizedComptonModelMute 687 G4AutoLock l(&LowEPPolarizedComptonModelMutex); 688 if(!data[Z]) { ReadData(Z); } 688 if(!data[Z]) { ReadData(Z); } 689 l.unlock(); 689 l.unlock(); 690 } 690 } 691 691 692 //******************************************** 692 //**************************************************************************** 693 693 694 //Fitting data to compute scattering function 694 //Fitting data to compute scattering function 695 const G4double G4LowEPPolarizedComptonModel::S 695 const G4double G4LowEPPolarizedComptonModel::ScatFuncFitParam[101][9] = { 696 { 0, 0., 0., 0., 0., 696 { 0, 0., 0., 0., 0., 0., 0., 0., 0.}, 697 { 1, 6.673, 1.49968E+08, -14.352, 1.999, -1 697 { 1, 6.673, 1.49968E+08, -14.352, 1.999, -143.374, 50.787, -5.951, 0.2304418}, 698 { 2, 6.500, 2.50035E+08, -14.215, 1.970, -5 698 { 2, 6.500, 2.50035E+08, -14.215, 1.970, -53.649, 13.892, -0.948, 0.006996759}, 699 { 3, 6.551, 3.99945E+08, -13.555, 1.993, -6 699 { 3, 6.551, 3.99945E+08, -13.555, 1.993, -62.090, 21.462, -2.453, 0.093416}, 700 { 4, 6.500, 5.00035E+08, -13.746, 1.998, -1 700 { 4, 6.500, 5.00035E+08, -13.746, 1.998, -127.906, 46.491, -5.614, 0.2262103}, 701 { 5, 6.500, 5.99791E+08, -13.800, 1.998, -1 701 { 5, 6.500, 5.99791E+08, -13.800, 1.998, -131.153, 47.132, -5.619, 0.2233819}, 702 { 6, 6.708, 6.99842E+08, -13.885, 1.999, -1 702 { 6, 6.708, 6.99842E+08, -13.885, 1.999, -128.143, 45.379, -5.325, 0.2083009}, 703 { 7, 6.685, 7.99834E+08, -13.885, 2.000, -1 703 { 7, 6.685, 7.99834E+08, -13.885, 2.000, -131.048, 46.314, -5.421, 0.2114925}, 704 { 8, 6.669, 7.99834E+08, -13.962, 2.001, -1 704 { 8, 6.669, 7.99834E+08, -13.962, 2.001, -128.225, 44.818, -5.183, 0.1997155}, 705 { 9, 6.711, 7.99834E+08, -13.999, 2.000, -1 705 { 9, 6.711, 7.99834E+08, -13.999, 2.000, -122.112, 42.103, -4.796, 0.1819099}, 706 { 10, 6.702, 7.99834E+08, -14.044, 1.999, -1 706 { 10, 6.702, 7.99834E+08, -14.044, 1.999, -110.143, 37.225, -4.143, 0.1532094}, 707 { 11, 6.425, 1.00000E+09, -13.423, 1.993, -4 707 { 11, 6.425, 1.00000E+09, -13.423, 1.993, -41.137, 12.313, -1.152, 0.03384553}, 708 { 12, 6.542, 1.00000E+09, -13.389, 1.997, -5 708 { 12, 6.542, 1.00000E+09, -13.389, 1.997, -53.549, 17.420, -1.840, 0.06431849}, 709 { 13, 6.570, 1.49968E+09, -13.401, 1.997, -6 709 { 13, 6.570, 1.49968E+09, -13.401, 1.997, -66.243, 22.297, -2.460, 0.09045854}, 710 { 14, 6.364, 1.49968E+09, -13.452, 1.999, -7 710 { 14, 6.364, 1.49968E+09, -13.452, 1.999, -78.271, 26.757, -3.008, 0.1128195}, 711 { 15, 6.500, 1.49968E+09, -13.488, 1.998, -8 711 { 15, 6.500, 1.49968E+09, -13.488, 1.998, -85.069, 29.164, -3.291, 0.1239113}, 712 { 16, 6.500, 1.49968E+09, -13.532, 1.998, -9 712 { 16, 6.500, 1.49968E+09, -13.532, 1.998, -93.640, 32.274, -3.665, 0.1388633}, 713 { 17, 6.500, 1.49968E+09, -13.584, 2.000, -9 713 { 17, 6.500, 1.49968E+09, -13.584, 2.000, -98.534, 33.958, -3.857, 0.1461557}, 714 { 18, 6.500, 1.49968E+09, -13.618, 1.999, -1 714 { 18, 6.500, 1.49968E+09, -13.618, 1.999, -100.077, 34.379, -3.891, 0.1468902}, 715 { 19, 6.500, 1.99986E+09, -13.185, 1.992, -5 715 { 19, 6.500, 1.99986E+09, -13.185, 1.992, -53.819, 17.528, -1.851, 0.0648722}, 716 { 20, 6.490, 1.99986E+09, -13.123, 1.993, -5 716 { 20, 6.490, 1.99986E+09, -13.123, 1.993, -52.221, 17.169, -1.832, 0.06502094}, 717 { 21, 6.498, 1.99986E+09, -13.157, 1.994, -5 717 { 21, 6.498, 1.99986E+09, -13.157, 1.994, -55.365, 18.276, -1.961, 0.07002778}, 718 { 22, 6.495, 1.99986E+09, -13.183, 1.994, -5 718 { 22, 6.495, 1.99986E+09, -13.183, 1.994, -57.412, 18.957, -2.036, 0.07278856}, 719 { 23, 6.487, 1.99986E+09, -13.216, 1.995, -5 719 { 23, 6.487, 1.99986E+09, -13.216, 1.995, -58.478, 19.270, -2.065, 0.07362722}, 720 { 24, 6.500, 1.99986E+09, -13.330, 1.997, -6 720 { 24, 6.500, 1.99986E+09, -13.330, 1.997, -62.192, 20.358, -2.167, 0.07666583}, 721 { 25, 6.488, 1.99986E+09, -13.277, 1.997, -5 721 { 25, 6.488, 1.99986E+09, -13.277, 1.997, -58.007, 18.924, -2.003, 0.0704305}, 722 { 26, 6.500, 5.00035E+09, -13.292, 1.997, -6 722 { 26, 6.500, 5.00035E+09, -13.292, 1.997, -61.176, 20.067, -2.141, 0.0760269}, 723 { 27, 6.500, 5.00035E+09, -13.321, 1.998, -6 723 { 27, 6.500, 5.00035E+09, -13.321, 1.998, -61.909, 20.271, -2.159, 0.07653559}, 724 { 28, 6.500, 5.00035E+09, -13.340, 1.998, -6 724 { 28, 6.500, 5.00035E+09, -13.340, 1.998, -62.402, 20.391, -2.167, 0.07664243}, 725 { 29, 6.500, 5.00035E+09, -13.439, 1.998, -6 725 { 29, 6.500, 5.00035E+09, -13.439, 1.998, -67.305, 21.954, -2.331, 0.0823267}, 726 { 30, 6.500, 5.00035E+09, -13.383, 1.999, -6 726 { 30, 6.500, 5.00035E+09, -13.383, 1.999, -62.064, 20.136, -2.122, 0.07437589}, 727 { 31, 6.500, 5.00035E+09, -13.349, 1.997, -6 727 { 31, 6.500, 5.00035E+09, -13.349, 1.997, -61.068, 19.808, -2.086, 0.07307488}, 728 { 32, 6.500, 5.00035E+09, -13.373, 1.999, -6 728 { 32, 6.500, 5.00035E+09, -13.373, 1.999, -63.126, 20.553, -2.175, 0.07660222}, 729 { 33, 6.500, 5.00035E+09, -13.395, 1.999, -6 729 { 33, 6.500, 5.00035E+09, -13.395, 1.999, -65.674, 21.445, -2.278, 0.08054694}, 730 { 34, 6.500, 5.00035E+09, -13.417, 1.999, -6 730 { 34, 6.500, 5.00035E+09, -13.417, 1.999, -69.457, 22.811, -2.442, 0.08709536}, 731 { 35, 6.500, 5.00035E+09, -13.442, 2.000, -7 731 { 35, 6.500, 5.00035E+09, -13.442, 2.000, -72.283, 23.808, -2.558, 0.09156808}, 732 { 36, 6.500, 5.00035E+09, -13.451, 1.998, -7 732 { 36, 6.500, 5.00035E+09, -13.451, 1.998, -74.696, 24.641, -2.653, 0.09516597}, 733 { 37, 6.500, 5.00035E+09, -13.082, 1.991, -4 733 { 37, 6.500, 5.00035E+09, -13.082, 1.991, -46.235, 14.519, -1.458, 0.04837659}, 734 { 38, 6.465, 5.00035E+09, -13.022, 1.993, -4 734 { 38, 6.465, 5.00035E+09, -13.022, 1.993, -41.784, 13.065, -1.300, 0.04267703}, 735 { 39, 6.492, 5.00035E+09, -13.043, 1.994, -4 735 { 39, 6.492, 5.00035E+09, -13.043, 1.994, -44.609, 14.114, -1.429, 0.0479348}, 736 { 40, 6.499, 5.00035E+09, -13.064, 1.994, -4 736 { 40, 6.499, 5.00035E+09, -13.064, 1.994, -47.142, 15.019, -1.536, 0.0521347}, 737 { 41, 6.384, 5.00035E+09, -13.156, 1.996, -5 737 { 41, 6.384, 5.00035E+09, -13.156, 1.996, -53.114, 17.052, -1.766, 0.06079426}, 738 { 42, 6.500, 5.00035E+09, -13.176, 1.996, -5 738 { 42, 6.500, 5.00035E+09, -13.176, 1.996, -54.590, 17.550, -1.822, 0.06290335}, 739 { 43, 6.500, 5.00035E+09, -13.133, 1.997, -5 739 { 43, 6.500, 5.00035E+09, -13.133, 1.997, -51.272, 16.423, -1.694, 0.05806108}, 740 { 44, 6.500, 5.00035E+09, -13.220, 1.996, -5 740 { 44, 6.500, 5.00035E+09, -13.220, 1.996, -58.314, 18.839, -1.969, 0.0684608}, 741 { 45, 6.500, 5.00035E+09, -13.246, 1.998, -5 741 { 45, 6.500, 5.00035E+09, -13.246, 1.998, -59.674, 19.295, -2.020, 0.07037294}, 742 { 46, 6.500, 5.00035E+09, -13.407, 1.999, -7 742 { 46, 6.500, 5.00035E+09, -13.407, 1.999, -72.228, 23.693, -2.532, 0.09017969}, 743 { 47, 6.500, 5.00035E+09, -13.277, 1.998, -6 743 { 47, 6.500, 5.00035E+09, -13.277, 1.998, -60.890, 19.647, -2.053, 0.07138694}, 744 { 48, 6.500, 5.00035E+09, -13.222, 1.998, -5 744 { 48, 6.500, 5.00035E+09, -13.222, 1.998, -56.152, 18.002, -1.863, 0.06410123}, 745 { 49, 6.500, 5.00035E+09, -13.199, 1.997, -5 745 { 49, 6.500, 5.00035E+09, -13.199, 1.997, -56.208, 18.052, -1.872, 0.06456884}, 746 { 50, 6.500, 5.00035E+09, -13.215, 1.998, -5 746 { 50, 6.500, 5.00035E+09, -13.215, 1.998, -58.478, 18.887, -1.973, 0.06860356}, 747 { 51, 6.500, 5.00035E+09, -13.230, 1.998, -6 747 { 51, 6.500, 5.00035E+09, -13.230, 1.998, -60.708, 19.676, -2.066, 0.07225841}, 748 { 52, 6.500, 7.99834E+09, -13.246, 1.998, -6 748 { 52, 6.500, 7.99834E+09, -13.246, 1.998, -63.341, 20.632, -2.180, 0.0767412}, 749 { 53, 6.500, 5.00035E+09, -13.262, 1.998, -6 749 { 53, 6.500, 5.00035E+09, -13.262, 1.998, -66.339, 21.716, -2.310, 0.08191981}, 750 { 54, 6.500, 7.99834E+09, -13.279, 1.998, -6 750 { 54, 6.500, 7.99834E+09, -13.279, 1.998, -67.649, 22.151, -2.357, 0.08357825}, 751 { 55, 6.500, 5.00035E+09, -12.951, 1.990, -4 751 { 55, 6.500, 5.00035E+09, -12.951, 1.990, -45.302, 14.219, -1.423, 0.04712317}, 752 { 56, 6.425, 5.00035E+09, -12.882, 1.992, -3 752 { 56, 6.425, 5.00035E+09, -12.882, 1.992, -39.825, 12.363, -1.214, 0.03931009}, 753 { 57, 6.466, 2.82488E+09, -12.903, 1.992, -3 753 { 57, 6.466, 2.82488E+09, -12.903, 1.992, -38.952, 11.982, -1.160, 0.03681554}, 754 { 58, 6.451, 5.00035E+09, -12.915, 1.993, -4 754 { 58, 6.451, 5.00035E+09, -12.915, 1.993, -41.959, 13.118, -1.302, 0.04271291}, 755 { 59, 6.434, 5.00035E+09, -12.914, 1.993, -4 755 { 59, 6.434, 5.00035E+09, -12.914, 1.993, -40.528, 12.555, -1.230, 0.03971407}, 756 { 60, 6.444, 5.00035E+09, -12.922, 1.992, -3 756 { 60, 6.444, 5.00035E+09, -12.922, 1.992, -39.986, 12.329, -1.200, 0.03843737}, 757 { 61, 6.414, 7.99834E+09, -12.930, 1.993, -4 757 { 61, 6.414, 7.99834E+09, -12.930, 1.993, -42.756, 13.362, -1.327, 0.0436124}, 758 { 62, 6.420, 7.99834E+09, -12.938, 1.992, -4 758 { 62, 6.420, 7.99834E+09, -12.938, 1.992, -42.682, 13.314, -1.319, 0.04322509}, 759 { 63, 6.416, 7.99834E+09, -12.946, 1.993, -4 759 { 63, 6.416, 7.99834E+09, -12.946, 1.993, -42.399, 13.185, -1.301, 0.04243861}, 760 { 64, 6.443, 7.99834E+09, -12.963, 1.993, -4 760 { 64, 6.443, 7.99834E+09, -12.963, 1.993, -43.226, 13.475, -1.335, 0.04377341}, 761 { 65, 6.449, 7.99834E+09, -12.973, 1.993, -4 761 { 65, 6.449, 7.99834E+09, -12.973, 1.993, -43.232, 13.456, -1.330, 0.04347536}, 762 { 66, 6.419, 7.99834E+09, -12.966, 1.993, -4 762 { 66, 6.419, 7.99834E+09, -12.966, 1.993, -42.047, 12.990, -1.270, 0.04095499}, 763 { 67, 6.406, 7.99834E+09, -12.976, 1.993, -4 763 { 67, 6.406, 7.99834E+09, -12.976, 1.993, -42.405, 13.106, -1.283, 0.04146024}, 764 { 68, 6.424, 7.99834E+09, -12.986, 1.993, -4 764 { 68, 6.424, 7.99834E+09, -12.986, 1.993, -41.974, 12.926, -1.259, 0.040435}, 765 { 69, 6.417, 7.99834E+09, -12.989, 1.993, -4 765 { 69, 6.417, 7.99834E+09, -12.989, 1.993, -42.132, 12.967, -1.262, 0.04048908}, 766 { 70, 6.405, 7.99834E+09, -13.000, 1.994, -4 766 { 70, 6.405, 7.99834E+09, -13.000, 1.994, -42.582, 13.122, -1.280, 0.04119599}, 767 { 71, 6.449, 7.99834E+09, -13.015, 1.994, -4 767 { 71, 6.449, 7.99834E+09, -13.015, 1.994, -42.586, 13.115, -1.278, 0.04107587}, 768 { 72, 6.465, 7.99834E+09, -13.030, 1.994, -4 768 { 72, 6.465, 7.99834E+09, -13.030, 1.994, -43.708, 13.509, -1.324, 0.04286491}, 769 { 73, 6.447, 7.99834E+09, -13.048, 1.996, -4 769 { 73, 6.447, 7.99834E+09, -13.048, 1.996, -44.838, 13.902, -1.369, 0.04457132}, 770 { 74, 6.452, 7.99834E+09, -13.073, 1.997, -4 770 { 74, 6.452, 7.99834E+09, -13.073, 1.997, -45.545, 14.137, -1.395, 0.04553459}, 771 { 75, 6.432, 7.99834E+09, -13.082, 1.997, -4 771 { 75, 6.432, 7.99834E+09, -13.082, 1.997, -46.426, 14.431, -1.428, 0.04678218}, 772 { 76, 6.439, 7.99834E+09, -13.100, 1.997, -4 772 { 76, 6.439, 7.99834E+09, -13.100, 1.997, -47.513, 14.806, -1.471, 0.04842566}, 773 { 77, 6.432, 7.99834E+09, -13.110, 1.997, -4 773 { 77, 6.432, 7.99834E+09, -13.110, 1.997, -48.225, 15.042, -1.497, 0.04938364}, 774 { 78, 6.500, 7.99834E+09, -13.185, 1.997, -5 774 { 78, 6.500, 7.99834E+09, -13.185, 1.997, -53.256, 16.739, -1.687, 0.05645173}, 775 { 79, 6.500, 7.99834E+09, -13.200, 1.997, -5 775 { 79, 6.500, 7.99834E+09, -13.200, 1.997, -53.900, 16.946, -1.709, 0.05723134}, 776 { 80, 6.500, 7.99834E+09, -13.156, 1.998, -4 776 { 80, 6.500, 7.99834E+09, -13.156, 1.998, -49.801, 15.536, -1.547, 0.05103522}, 777 { 81, 6.500, 7.99834E+09, -13.128, 1.997, -4 777 { 81, 6.500, 7.99834E+09, -13.128, 1.997, -49.651, 15.512, -1.548, 0.05123203}, 778 { 82, 6.500, 7.99834E+09, -13.134, 1.997, -5 778 { 82, 6.500, 7.99834E+09, -13.134, 1.997, -51.021, 16.018, -1.609, 0.05364831}, 779 { 83, 6.500, 7.99834E+09, -13.148, 1.998, -5 779 { 83, 6.500, 7.99834E+09, -13.148, 1.998, -52.693, 16.612, -1.679, 0.05638698}, 780 { 84, 6.500, 7.99834E+09, -13.161, 1.998, -5 780 { 84, 6.500, 7.99834E+09, -13.161, 1.998, -54.415, 17.238, -1.754, 0.05935566}, 781 { 85, 6.500, 7.99834E+09, -13.175, 1.998, -5 781 { 85, 6.500, 7.99834E+09, -13.175, 1.998, -56.083, 17.834, -1.824, 0.06206068}, 782 { 86, 6.500, 7.99834E+09, -13.189, 1.998, -5 782 { 86, 6.500, 7.99834E+09, -13.189, 1.998, -57.860, 18.463, -1.898, 0.0649633}, 783 { 87, 6.500, 7.99834E+09, -12.885, 1.990, -3 783 { 87, 6.500, 7.99834E+09, -12.885, 1.990, -39.973, 12.164, -1.162, 0.0364598}, 784 { 88, 6.417, 7.99834E+09, -12.816, 1.991, -3 784 { 88, 6.417, 7.99834E+09, -12.816, 1.991, -34.591, 10.338, -0.956, 0.0287409}, 785 { 89, 6.442, 7.99834E+09, -12.831, 1.992, -3 785 { 89, 6.442, 7.99834E+09, -12.831, 1.992, -36.002, 10.867, -1.021, 0.03136835}, 786 { 90, 6.463, 7.99834E+09, -12.850, 1.993, -3 786 { 90, 6.463, 7.99834E+09, -12.850, 1.993, -37.660, 11.475, -1.095, 0.03435334}, 787 { 91, 6.447, 7.99834E+09, -12.852, 1.993, -3 787 { 91, 6.447, 7.99834E+09, -12.852, 1.993, -37.268, 11.301, -1.071, 0.0330539}, 788 { 92, 6.439, 7.99834E+09, -12.858, 1.993, -3 788 { 92, 6.439, 7.99834E+09, -12.858, 1.993, -37.695, 11.438, -1.085, 0.03376669}, 789 { 93, 6.437, 1.00000E+10, -12.866, 1.993, -3 789 { 93, 6.437, 1.00000E+10, -12.866, 1.993, -39.010, 11.927, -1.146, 0.03630848}, 790 { 94, 6.432, 7.99834E+09, -12.862, 1.993, -3 790 { 94, 6.432, 7.99834E+09, -12.862, 1.993, -37.192, 11.229, -1.057, 0.0325621}, 791 { 95, 6.435, 7.99834E+09, -12.869, 1.993, -3 791 { 95, 6.435, 7.99834E+09, -12.869, 1.993, -37.589, 11.363, -1.072, 0.03312393}, 792 { 96, 6.449, 1.00000E+10, -12.886, 1.993, -3 792 { 96, 6.449, 1.00000E+10, -12.886, 1.993, -39.573, 12.095, -1.162, 0.03680527}, 793 { 97, 6.446, 1.00000E+10, -12.892, 1.993, -4 793 { 97, 6.446, 1.00000E+10, -12.892, 1.993, -40.007, 12.242, -1.178, 0.03737377}, 794 { 98, 6.421, 1.00000E+10, -12.887, 1.993, -3 794 { 98, 6.421, 1.00000E+10, -12.887, 1.993, -39.509, 12.041, -1.152, 0.03629023}, 795 { 99, 6.414, 1.00000E+10, -12.894, 1.993, -3 795 { 99, 6.414, 1.00000E+10, -12.894, 1.993, -39.939, 12.183, -1.168, 0.03690464}, 796 {100, 6.412, 1.00000E+10, -12.900, 1.993, -3 796 {100, 6.412, 1.00000E+10, -12.900, 1.993, -39.973, 12.180, -1.166, 0.036773} 797 }; 797 }; 798 798 799 //******************************************** 799 //**************************************************************************** 800 800 801 //Supporting functions for photon polarisation 801 //Supporting functions for photon polarisation effects 802 G4double G4LowEPPolarizedComptonModel::SetPhi( 802 G4double G4LowEPPolarizedComptonModel::SetPhi(G4double energyRate, 803 G4double sinT2) 803 G4double sinT2) 804 { 804 { 805 G4double rand1; 805 G4double rand1; 806 G4double rand2; 806 G4double rand2; 807 G4double phiProbability; 807 G4double phiProbability; 808 G4double phi; 808 G4double phi; 809 G4double a, b; 809 G4double a, b; 810 810 811 do 811 do 812 { 812 { 813 rand1 = G4UniformRand(); 813 rand1 = G4UniformRand(); 814 rand2 = G4UniformRand(); 814 rand2 = G4UniformRand(); 815 phiProbability=0.; 815 phiProbability=0.; 816 phi = twopi*rand1; 816 phi = twopi*rand1; 817 817 818 a = 2*sinT2; 818 a = 2*sinT2; 819 b = energyRate + 1/energyRate; 819 b = energyRate + 1/energyRate; 820 820 821 phiProbability = 1 - (a/b)*(std::cos(phi 821 phiProbability = 1 - (a/b)*(std::cos(phi)*std::cos(phi)); 822 } 822 } 823 while ( rand2 > phiProbability ); 823 while ( rand2 > phiProbability ); 824 return phi; 824 return phi; 825 } 825 } 826 826 827 //******************************************** 827 //**************************************************************************** 828 828 829 G4ThreeVector G4LowEPPolarizedComptonModel::Se 829 G4ThreeVector G4LowEPPolarizedComptonModel::SetPerpendicularVector(G4ThreeVector& a) 830 { 830 { 831 G4double dx = a.x(); 831 G4double dx = a.x(); 832 G4double dy = a.y(); 832 G4double dy = a.y(); 833 G4double dz = a.z(); 833 G4double dz = a.z(); 834 G4double x = dx < 0.0 ? -dx : dx; 834 G4double x = dx < 0.0 ? -dx : dx; 835 G4double y = dy < 0.0 ? -dy : dy; 835 G4double y = dy < 0.0 ? -dy : dy; 836 G4double z = dz < 0.0 ? -dz : dz; 836 G4double z = dz < 0.0 ? -dz : dz; 837 if (x < y) { 837 if (x < y) { 838 return x < z ? G4ThreeVector(-dy,dx,0) : G 838 return x < z ? G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy); 839 }else{ 839 }else{ 840 return y < z ? G4ThreeVector(dz,0,-dx) : G 840 return y < z ? G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0); 841 } 841 } 842 } 842 } 843 843 844 //******************************************** 844 //**************************************************************************** 845 845 846 G4ThreeVector G4LowEPPolarizedComptonModel::Ge 846 G4ThreeVector G4LowEPPolarizedComptonModel::GetRandomPolarization(G4ThreeVector& direction0) 847 { 847 { 848 G4ThreeVector d0 = direction0.unit(); 848 G4ThreeVector d0 = direction0.unit(); 849 G4ThreeVector a1 = SetPerpendicularVector(d0 849 G4ThreeVector a1 = SetPerpendicularVector(d0); //different orthogonal 850 G4ThreeVector a0 = a1.unit(); // unit vector 850 G4ThreeVector a0 = a1.unit(); // unit vector 851 851 852 G4double rand1 = G4UniformRand(); 852 G4double rand1 = G4UniformRand(); 853 853 854 G4double angle = twopi*rand1; // random pola 854 G4double angle = twopi*rand1; // random polar angle 855 G4ThreeVector b0 = d0.cross(a0); // cross pr 855 G4ThreeVector b0 = d0.cross(a0); // cross product 856 856 857 G4ThreeVector c; 857 G4ThreeVector c; 858 858 859 c.setX(std::cos(angle)*(a0.x())+std::sin(ang 859 c.setX(std::cos(angle)*(a0.x())+std::sin(angle)*b0.x()); 860 c.setY(std::cos(angle)*(a0.y())+std::sin(ang 860 c.setY(std::cos(angle)*(a0.y())+std::sin(angle)*b0.y()); 861 c.setZ(std::cos(angle)*(a0.z())+std::sin(ang 861 c.setZ(std::cos(angle)*(a0.z())+std::sin(angle)*b0.z()); 862 862 863 G4ThreeVector c0 = c.unit(); 863 G4ThreeVector c0 = c.unit(); 864 864 865 return c0; 865 return c0; 866 } 866 } 867 867 868 //******************************************** 868 //**************************************************************************** 869 869 870 G4ThreeVector G4LowEPPolarizedComptonModel::Ge 870 G4ThreeVector G4LowEPPolarizedComptonModel::GetPerpendicularPolarization 871 (const G4ThreeVector& photonDirection, const G 871 (const G4ThreeVector& photonDirection, const G4ThreeVector& photonPolarization) const 872 { 872 { 873 // 873 // 874 // The polarization of a photon is always pe 874 // The polarization of a photon is always perpendicular to its momentum direction. 875 // Therefore this function removes those vec 875 // Therefore this function removes those vector component of photonPolarization, which 876 // points in direction of photonDirection 876 // points in direction of photonDirection 877 // 877 // 878 // Mathematically we search the projection o 878 // Mathematically we search the projection of the vector a on the plane E, where n is the 879 // plains normal vector. 879 // plains normal vector. 880 // The basic equation can be found in each g 880 // The basic equation can be found in each geometry book (e.g. Bronstein): 881 // p = a - (a o n)/(n o n)*n 881 // p = a - (a o n)/(n o n)*n 882 882 883 return photonPolarization - photonPolarizati 883 return photonPolarization - photonPolarization.dot(photonDirection)/photonDirection.dot(photonDirection) * photonDirection; 884 } 884 } 885 885 886 //******************************************** 886 //**************************************************************************** 887 887 888 G4ThreeVector G4LowEPPolarizedComptonModel::Se 888 G4ThreeVector G4LowEPPolarizedComptonModel::SetNewPolarization(G4double LowEPPCepsilon, 889 G4double sinT2, 889 G4double sinT2, 890 G4double phi, 890 G4double phi, 891 G4double costheta) 891 G4double costheta) 892 { 892 { 893 G4double rand1; 893 G4double rand1; 894 G4double rand2; 894 G4double rand2; 895 G4double cosPhi = std::cos(phi); 895 G4double cosPhi = std::cos(phi); 896 G4double sinPhi = std::sin(phi); 896 G4double sinPhi = std::sin(phi); 897 G4double sinTheta = std::sqrt(sinT2); 897 G4double sinTheta = std::sqrt(sinT2); 898 G4double cosP2 = cosPhi*cosPhi; 898 G4double cosP2 = cosPhi*cosPhi; 899 G4double normalisation = std::sqrt(1. - cosP 899 G4double normalisation = std::sqrt(1. - cosP2*sinT2); 900 900 901 // Method based on: 901 // Method based on: 902 // D. Xu, Z. He and F. Zhang 902 // D. Xu, Z. He and F. Zhang 903 // "Detection of Gamma Ray Polarization Usin 903 // "Detection of Gamma Ray Polarization Using a 3-D Position Sensitive CdZnTe Detector" 904 // IEEE TNS, Vol. 52(4), 1160-1164, 2005. 904 // IEEE TNS, Vol. 52(4), 1160-1164, 2005. 905 905 906 // Determination of Theta 906 // Determination of Theta 907 907 908 G4double theta; 908 G4double theta; 909 rand1 = G4UniformRand(); 909 rand1 = G4UniformRand(); 910 rand2 = G4UniformRand(); 910 rand2 = G4UniformRand(); 911 911 912 if (rand1<(LowEPPCepsilon+1.0/LowEPPCepsilon 912 if (rand1<(LowEPPCepsilon+1.0/LowEPPCepsilon-2)/(2.0*(LowEPPCepsilon+1.0/LowEPPCepsilon)-4.0*sinT2*cosP2)) 913 { 913 { 914 if (rand2<0.5) 914 if (rand2<0.5) 915 theta = pi/2.0; 915 theta = pi/2.0; 916 else 916 else 917 theta = 3.0*pi/2.0; 917 theta = 3.0*pi/2.0; 918 } 918 } 919 else 919 else 920 { 920 { 921 if (rand2<0.5) 921 if (rand2<0.5) 922 theta = 0; 922 theta = 0; 923 else 923 else 924 theta = pi; 924 theta = pi; 925 } 925 } 926 G4double cosBeta = std::cos(theta); 926 G4double cosBeta = std::cos(theta); 927 G4double sinBeta = std::sqrt(1-cosBeta*cosBe 927 G4double sinBeta = std::sqrt(1-cosBeta*cosBeta); 928 928 929 G4ThreeVector photonPolarization1; 929 G4ThreeVector photonPolarization1; 930 930 931 G4double xParallel = normalisation*cosBeta; 931 G4double xParallel = normalisation*cosBeta; 932 G4double yParallel = -(sinT2*cosPhi*sinPhi)* 932 G4double yParallel = -(sinT2*cosPhi*sinPhi)*cosBeta/normalisation; 933 G4double zParallel = -(costheta*sinTheta*cos 933 G4double zParallel = -(costheta*sinTheta*cosPhi)*cosBeta/normalisation; 934 G4double xPerpendicular = 0.; 934 G4double xPerpendicular = 0.; 935 G4double yPerpendicular = (costheta)*sinBeta 935 G4double yPerpendicular = (costheta)*sinBeta/normalisation; 936 G4double zPerpendicular = -(sinTheta*sinPhi) 936 G4double zPerpendicular = -(sinTheta*sinPhi)*sinBeta/normalisation; 937 937 938 G4double xTotal = (xParallel + xPerpendicula 938 G4double xTotal = (xParallel + xPerpendicular); 939 G4double yTotal = (yParallel + yPerpendicula 939 G4double yTotal = (yParallel + yPerpendicular); 940 G4double zTotal = (zParallel + zPerpendicula 940 G4double zTotal = (zParallel + zPerpendicular); 941 941 942 photonPolarization1.setX(xTotal); 942 photonPolarization1.setX(xTotal); 943 photonPolarization1.setY(yTotal); 943 photonPolarization1.setY(yTotal); 944 photonPolarization1.setZ(zTotal); 944 photonPolarization1.setZ(zTotal); 945 945 946 return photonPolarization1; 946 return photonPolarization1; 947 947 948 } 948 } 949 949 950 //******************************************** 950 //**************************************************************************** 951 void G4LowEPPolarizedComptonModel::SystemOfRef 951 void G4LowEPPolarizedComptonModel::SystemOfRefChange(G4ThreeVector& direction0, 952 G4ThreeVector& direction1, 952 G4ThreeVector& direction1, 953 G4ThreeVector& polarization0, 953 G4ThreeVector& polarization0, 954 G4ThreeVector& polarization1) 954 G4ThreeVector& polarization1) 955 { 955 { 956 // direction0 is the original photon directi 956 // direction0 is the original photon direction ---> z 957 // polarization0 is the original photon pola 957 // polarization0 is the original photon polarization ---> x 958 // need to specify y axis in the real refere 958 // need to specify y axis in the real reference frame ---> y 959 G4ThreeVector Axis_Z0 = direction0.unit(); 959 G4ThreeVector Axis_Z0 = direction0.unit(); 960 G4ThreeVector Axis_X0 = polarization0.unit() 960 G4ThreeVector Axis_X0 = polarization0.unit(); 961 G4ThreeVector Axis_Y0 = (Axis_Z0.cross(Axis_ 961 G4ThreeVector Axis_Y0 = (Axis_Z0.cross(Axis_X0)).unit(); // to be confirmed; 962 962 963 G4double direction_x = direction1.getX(); 963 G4double direction_x = direction1.getX(); 964 G4double direction_y = direction1.getY(); 964 G4double direction_y = direction1.getY(); 965 G4double direction_z = direction1.getZ(); 965 G4double direction_z = direction1.getZ(); 966 966 967 direction1 = (direction_x*Axis_X0 + directio 967 direction1 = (direction_x*Axis_X0 + direction_y*Axis_Y0 + direction_z*Axis_Z0).unit(); 968 G4double polarization_x = polarization1.getX 968 G4double polarization_x = polarization1.getX(); 969 G4double polarization_y = polarization1.getY 969 G4double polarization_y = polarization1.getY(); 970 G4double polarization_z = polarization1.getZ 970 G4double polarization_z = polarization1.getZ(); 971 971 972 polarization1 = (polarization_x*Axis_X0 + po 972 polarization1 = (polarization_x*Axis_X0 + polarization_y*Axis_Y0 + polarization_z*Axis_Z0).unit(); 973 973 974 } 974 } 975 975 976 //******************************************** 976 //**************************************************************************** 977 void G4LowEPPolarizedComptonModel::SystemOfRef 977 void G4LowEPPolarizedComptonModel::SystemOfRefChangeElect(G4ThreeVector& pdirection, 978 G4ThreeVector& edirection, 978 G4ThreeVector& edirection, 979 G4ThreeVector& ppolarization) 979 G4ThreeVector& ppolarization) 980 { 980 { 981 // direction0 is the original photon directi 981 // direction0 is the original photon direction ---> z 982 // polarization0 is the original photon pola 982 // polarization0 is the original photon polarization ---> x 983 // need to specify y axis in the real refere 983 // need to specify y axis in the real reference frame ---> y 984 G4ThreeVector Axis_Z0 = pdirection.unit(); 984 G4ThreeVector Axis_Z0 = pdirection.unit(); 985 G4ThreeVector Axis_X0 = ppolarization.unit() 985 G4ThreeVector Axis_X0 = ppolarization.unit(); 986 G4ThreeVector Axis_Y0 = (Axis_Z0.cross(Axis_ 986 G4ThreeVector Axis_Y0 = (Axis_Z0.cross(Axis_X0)).unit(); // to be confirmed; 987 987 988 G4double direction_x = edirection.getX(); 988 G4double direction_x = edirection.getX(); 989 G4double direction_y = edirection.getY(); 989 G4double direction_y = edirection.getY(); 990 G4double direction_z = edirection.getZ(); 990 G4double direction_z = edirection.getZ(); 991 991 992 edirection = (direction_x*Axis_X0 + directio 992 edirection = (direction_x*Axis_X0 + direction_y*Axis_Y0 + direction_z*Axis_Z0).unit(); 993 } 993 } 994 994 995 995 996 996