[816] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | /////////////////////////////////////////////////////////////////////////////// |
---|
| 27 | // |
---|
| 28 | // MODULE: G4SPSAngDistribution.cc |
---|
| 29 | // |
---|
| 30 | // Version: 1.0 |
---|
| 31 | // Date: 5/02/04 |
---|
| 32 | // Author: Fan Lei |
---|
| 33 | // Organisation: QinetiQ ltd. |
---|
| 34 | // Customer: ESA/ESTEC |
---|
| 35 | // |
---|
| 36 | /////////////////////////////////////////////////////////////////////////////// |
---|
| 37 | // |
---|
| 38 | // |
---|
| 39 | // CHANGE HISTORY |
---|
| 40 | // -------------- |
---|
| 41 | // |
---|
| 42 | // |
---|
| 43 | // Version 1.0, 05/02/2004, Fan Lei, Created. |
---|
| 44 | // Based on the G4GeneralParticleSource class in Geant4 v6.0 |
---|
| 45 | // |
---|
| 46 | /////////////////////////////////////////////////////////////////////////////// |
---|
| 47 | // |
---|
| 48 | #include "Randomize.hh" |
---|
| 49 | #include "G4SPSAngDistribution.hh" |
---|
| 50 | |
---|
| 51 | G4SPSAngDistribution::G4SPSAngDistribution() |
---|
| 52 | { |
---|
| 53 | // Angular distribution Variables |
---|
| 54 | G4ThreeVector zero; |
---|
| 55 | particle_momentum_direction = G4ParticleMomentum(0,0,-1); |
---|
| 56 | |
---|
| 57 | AngDistType = "planar"; |
---|
| 58 | AngRef1 = CLHEP::HepXHat; |
---|
| 59 | AngRef2 = CLHEP::HepYHat; |
---|
| 60 | AngRef3 = CLHEP::HepZHat; |
---|
| 61 | MinTheta = 0.; |
---|
| 62 | MaxTheta = pi; |
---|
| 63 | MinPhi = 0.; |
---|
| 64 | MaxPhi = twopi; |
---|
| 65 | DR = 0.; |
---|
| 66 | DX = 0.; |
---|
| 67 | DY = 0.; |
---|
| 68 | FocusPoint = G4ThreeVector(0., 0., 0.); |
---|
| 69 | UserDistType = "NULL"; |
---|
| 70 | UserWRTSurface = true; |
---|
| 71 | UserAngRef = false; |
---|
| 72 | IPDFThetaExist = false; |
---|
| 73 | IPDFPhiExist = false; |
---|
| 74 | verbosityLevel = 0 ; |
---|
| 75 | } |
---|
| 76 | |
---|
| 77 | G4SPSAngDistribution::~G4SPSAngDistribution() |
---|
| 78 | {} |
---|
| 79 | |
---|
| 80 | // |
---|
| 81 | void G4SPSAngDistribution::SetAngDistType(G4String atype) |
---|
| 82 | { |
---|
| 83 | if(atype != "iso" && atype != "cos" && atype != "user" && atype != "planar" |
---|
| 84 | && atype != "beam1d" && atype != "beam2d" && atype != "focused") |
---|
| 85 | G4cout << "Error, distribution must be iso, cos, planar, beam1d, beam2d, focused or user" << G4endl; |
---|
| 86 | else |
---|
| 87 | AngDistType = atype; |
---|
| 88 | if (AngDistType == "cos") MaxTheta = pi/2. ; |
---|
| 89 | if (AngDistType == "user") { |
---|
| 90 | UDefThetaH = IPDFThetaH = ZeroPhysVector ; |
---|
| 91 | IPDFThetaExist = false ; |
---|
| 92 | UDefPhiH = IPDFPhiH = ZeroPhysVector ; |
---|
| 93 | IPDFPhiExist = false ; |
---|
| 94 | } |
---|
| 95 | } |
---|
| 96 | |
---|
| 97 | void G4SPSAngDistribution::DefineAngRefAxes(G4String refname, G4ThreeVector ref) |
---|
| 98 | { |
---|
| 99 | if(refname == "angref1") |
---|
| 100 | AngRef1 = ref.unit(); // x' |
---|
| 101 | else if(refname == "angref2") |
---|
| 102 | AngRef2 = ref.unit(); // vector in x'y' plane |
---|
| 103 | |
---|
| 104 | // User defines x' (AngRef1) and a vector in the x'y' |
---|
| 105 | // plane (AngRef2). Then, AngRef1 x AngRef2 = AngRef3 |
---|
| 106 | // the z' vector. Then, AngRef3 x AngRef1 = AngRef2 |
---|
| 107 | // which will now be y'. |
---|
| 108 | |
---|
| 109 | AngRef3 = AngRef1.cross(AngRef2); // z' |
---|
| 110 | AngRef2 = AngRef3.cross(AngRef1); // y' |
---|
| 111 | UserAngRef = true ; |
---|
| 112 | if(verbosityLevel == 2) |
---|
| 113 | { |
---|
| 114 | G4cout << "Angular distribution rotation axes " << AngRef1 << " " << AngRef2 << " " << AngRef3 << G4endl; |
---|
| 115 | } |
---|
| 116 | } |
---|
| 117 | |
---|
| 118 | void G4SPSAngDistribution::SetMinTheta(G4double mint) |
---|
| 119 | { |
---|
| 120 | MinTheta = mint; |
---|
| 121 | } |
---|
| 122 | |
---|
| 123 | void G4SPSAngDistribution::SetMinPhi(G4double minp) |
---|
| 124 | { |
---|
| 125 | MinPhi = minp; |
---|
| 126 | } |
---|
| 127 | |
---|
| 128 | void G4SPSAngDistribution::SetMaxTheta(G4double maxt) |
---|
| 129 | { |
---|
| 130 | MaxTheta = maxt; |
---|
| 131 | } |
---|
| 132 | |
---|
| 133 | void G4SPSAngDistribution::SetMaxPhi(G4double maxp) |
---|
| 134 | { |
---|
| 135 | MaxPhi = maxp; |
---|
| 136 | } |
---|
| 137 | |
---|
| 138 | void G4SPSAngDistribution::SetBeamSigmaInAngR(G4double r) |
---|
| 139 | { |
---|
| 140 | DR = r; |
---|
| 141 | } |
---|
| 142 | |
---|
| 143 | void G4SPSAngDistribution::SetBeamSigmaInAngX(G4double r) |
---|
| 144 | { |
---|
| 145 | DX = r; |
---|
| 146 | } |
---|
| 147 | |
---|
| 148 | void G4SPSAngDistribution::SetBeamSigmaInAngY(G4double r) |
---|
| 149 | { |
---|
| 150 | DY = r; |
---|
| 151 | } |
---|
| 152 | |
---|
| 153 | void G4SPSAngDistribution::UserDefAngTheta(G4ThreeVector input) |
---|
| 154 | { |
---|
| 155 | if(UserDistType == "NULL") UserDistType = "theta"; |
---|
| 156 | if(UserDistType == "phi") UserDistType = "both"; |
---|
| 157 | G4double thi, val; |
---|
| 158 | thi = input.x(); |
---|
| 159 | val = input.y(); |
---|
| 160 | if(verbosityLevel >= 1) |
---|
| 161 | G4cout << "In UserDefAngTheta" << G4endl; |
---|
| 162 | UDefThetaH.InsertValues(thi, val); |
---|
| 163 | } |
---|
| 164 | |
---|
| 165 | void G4SPSAngDistribution::UserDefAngPhi(G4ThreeVector input) |
---|
| 166 | { |
---|
| 167 | if(UserDistType == "NULL") UserDistType = "phi"; |
---|
| 168 | if(UserDistType == "theta") UserDistType = "both"; |
---|
| 169 | G4double phhi, val; |
---|
| 170 | phhi = input.x(); |
---|
| 171 | val = input.y(); |
---|
| 172 | if(verbosityLevel >= 1) |
---|
| 173 | G4cout << "In UserDefAngPhi" << G4endl; |
---|
| 174 | UDefPhiH.InsertValues(phhi, val); |
---|
| 175 | } |
---|
| 176 | |
---|
| 177 | void G4SPSAngDistribution::SetFocusPoint(G4ThreeVector input) |
---|
| 178 | { |
---|
| 179 | FocusPoint = input; |
---|
| 180 | } |
---|
| 181 | |
---|
| 182 | void G4SPSAngDistribution::SetUserWRTSurface(G4bool wrtSurf) |
---|
| 183 | { |
---|
| 184 | // This is only applied in user mode? |
---|
| 185 | // if UserWRTSurface = true then the user wants momenta with respect |
---|
| 186 | // to the surface normals. |
---|
| 187 | // When doing this theta has to be 0-90 only otherwise there will be |
---|
| 188 | // errors, which currently are flagged anywhere. |
---|
| 189 | UserWRTSurface = wrtSurf; |
---|
| 190 | } |
---|
| 191 | |
---|
| 192 | void G4SPSAngDistribution::SetUseUserAngAxis(G4bool userang) |
---|
| 193 | { |
---|
| 194 | // if UserAngRef = true the angular distribution is defined wrt |
---|
| 195 | // the user defined co-ordinates |
---|
| 196 | UserAngRef = userang; |
---|
| 197 | } |
---|
| 198 | |
---|
| 199 | void G4SPSAngDistribution::GenerateBeamFlux() |
---|
| 200 | { |
---|
| 201 | G4double theta, phi; |
---|
| 202 | G4double px, py, pz; |
---|
| 203 | if (AngDistType == "beam1d") |
---|
| 204 | { |
---|
| 205 | theta = G4RandGauss::shoot(0.0,DR); |
---|
| 206 | phi = twopi * G4UniformRand(); |
---|
| 207 | } |
---|
| 208 | else |
---|
| 209 | { |
---|
| 210 | px = G4RandGauss::shoot(0.0,DX); |
---|
| 211 | py = G4RandGauss::shoot(0.0,DY); |
---|
| 212 | theta = std::sqrt (px*px + py*py); |
---|
| 213 | if (theta != 0.) { |
---|
| 214 | phi = std::acos(px/theta); |
---|
| 215 | if ( py < 0.) phi = -phi; |
---|
| 216 | } |
---|
| 217 | else |
---|
| 218 | { |
---|
| 219 | phi = 0.0; |
---|
| 220 | } |
---|
| 221 | } |
---|
| 222 | px = -std::sin(theta) * std::cos(phi); |
---|
| 223 | py = -std::sin(theta) * std::sin(phi); |
---|
| 224 | pz = -std::cos(theta); |
---|
| 225 | G4double finx, finy, finz ; |
---|
| 226 | finx = px, finy =py, finz =pz; |
---|
| 227 | if (UserAngRef){ |
---|
| 228 | // Apply Angular Rotation Matrix |
---|
| 229 | // x * AngRef1, y * AngRef2 and z * AngRef3 |
---|
| 230 | finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x()); |
---|
| 231 | finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y()); |
---|
| 232 | finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z()); |
---|
| 233 | G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz)); |
---|
| 234 | finx = finx/ResMag; |
---|
| 235 | finy = finy/ResMag; |
---|
| 236 | finz = finz/ResMag; |
---|
| 237 | } |
---|
| 238 | particle_momentum_direction.setX(finx); |
---|
| 239 | particle_momentum_direction.setY(finy); |
---|
| 240 | particle_momentum_direction.setZ(finz); |
---|
| 241 | |
---|
| 242 | // particle_momentum_direction now holds unit momentum vector. |
---|
| 243 | if(verbosityLevel >= 1) |
---|
| 244 | G4cout << "Generating beam vector: " << particle_momentum_direction << G4endl; |
---|
| 245 | } |
---|
| 246 | |
---|
| 247 | void G4SPSAngDistribution::GenerateFocusedFlux() |
---|
| 248 | { |
---|
| 249 | particle_momentum_direction = (FocusPoint - posDist->particle_position).unit(); |
---|
| 250 | // |
---|
| 251 | // particle_momentum_direction now holds unit momentum vector. |
---|
| 252 | if(verbosityLevel >= 1) |
---|
| 253 | G4cout << "Generating focused vector: " << particle_momentum_direction << G4endl; |
---|
| 254 | } |
---|
| 255 | |
---|
| 256 | void G4SPSAngDistribution::GenerateIsotropicFlux() |
---|
| 257 | { |
---|
| 258 | // generates isotropic flux. |
---|
| 259 | // No vectors are needed. |
---|
| 260 | G4double rndm, rndm2; |
---|
| 261 | G4double px, py, pz; |
---|
| 262 | |
---|
| 263 | // |
---|
| 264 | G4double sintheta, sinphi,costheta,cosphi; |
---|
| 265 | rndm = angRndm->GenRandTheta(); |
---|
| 266 | costheta = std::cos(MinTheta) - rndm * (std::cos(MinTheta) - std::cos(MaxTheta)); |
---|
| 267 | sintheta = std::sqrt(1. - costheta*costheta); |
---|
| 268 | |
---|
| 269 | rndm2 = angRndm->GenRandPhi(); |
---|
| 270 | Phi = MinPhi + (MaxPhi - MinPhi) * rndm2; |
---|
| 271 | sinphi = std::sin(Phi); |
---|
| 272 | cosphi = std::cos(Phi); |
---|
| 273 | |
---|
| 274 | px = -sintheta * cosphi; |
---|
| 275 | py = -sintheta * sinphi; |
---|
| 276 | pz = -costheta; |
---|
| 277 | |
---|
| 278 | // for volume and ponit source use mother or user defined co-ordinates |
---|
| 279 | // for plane and surface source user surface-normal or userdefined co-ordinates |
---|
| 280 | // |
---|
| 281 | G4double finx, finy, finz; |
---|
| 282 | if (posDist->SourcePosType == "Point" || posDist->SourcePosType == "Volume") { |
---|
| 283 | if (UserAngRef){ |
---|
| 284 | // Apply Rotation Matrix |
---|
| 285 | // x * AngRef1, y * AngRef2 and z * AngRef3 |
---|
| 286 | finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x()); |
---|
| 287 | finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y()); |
---|
| 288 | finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z()); |
---|
| 289 | } else { |
---|
| 290 | finx = px; |
---|
| 291 | finy = py; |
---|
| 292 | finz = pz; |
---|
| 293 | } |
---|
| 294 | } else { // for plane and surface source |
---|
| 295 | if (UserAngRef){ |
---|
| 296 | // Apply Rotation Matrix |
---|
| 297 | // x * AngRef1, y * AngRef2 and z * AngRef3 |
---|
| 298 | finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x()); |
---|
| 299 | finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y()); |
---|
| 300 | finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z()); |
---|
| 301 | } else { |
---|
| 302 | finx = (px*posDist->SideRefVec1.x()) + (py*posDist->SideRefVec2.x()) + (pz*posDist->SideRefVec3.x()); |
---|
| 303 | finy = (px*posDist->SideRefVec1.y()) + (py*posDist->SideRefVec2.y()) + (pz*posDist->SideRefVec3.y()); |
---|
| 304 | finz = (px*posDist->SideRefVec1.z()) + (py*posDist->SideRefVec2.z()) + (pz*posDist->SideRefVec3.z()); |
---|
| 305 | } |
---|
| 306 | } |
---|
| 307 | G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz)); |
---|
| 308 | finx = finx/ResMag; |
---|
| 309 | finy = finy/ResMag; |
---|
| 310 | finz = finz/ResMag; |
---|
| 311 | |
---|
| 312 | particle_momentum_direction.setX(finx); |
---|
| 313 | particle_momentum_direction.setY(finy); |
---|
| 314 | particle_momentum_direction.setZ(finz); |
---|
| 315 | |
---|
| 316 | // particle_momentum_direction now holds unit momentum vector. |
---|
| 317 | if(verbosityLevel >= 1) |
---|
| 318 | G4cout << "Generating isotropic vector: " << particle_momentum_direction << G4endl; |
---|
| 319 | } |
---|
| 320 | |
---|
| 321 | void G4SPSAngDistribution::GenerateCosineLawFlux() |
---|
| 322 | { |
---|
| 323 | // Method to generate flux distributed with a cosine law |
---|
| 324 | G4double px, py, pz; |
---|
| 325 | G4double rndm, rndm2; |
---|
| 326 | // |
---|
| 327 | G4double sintheta, sinphi,costheta,cosphi; |
---|
| 328 | rndm = angRndm->GenRandTheta(); |
---|
| 329 | sintheta = std::sqrt( rndm * (std::sin(MaxTheta)*std::sin(MaxTheta) - std::sin(MinTheta)*std::sin(MinTheta) ) |
---|
| 330 | +std::sin(MinTheta)*std::sin(MinTheta) ); |
---|
| 331 | costheta = std::sqrt(1. -sintheta*sintheta); |
---|
| 332 | |
---|
| 333 | rndm2 = angRndm->GenRandPhi(); |
---|
| 334 | Phi = MinPhi + (MaxPhi - MinPhi) * rndm2; |
---|
| 335 | sinphi = std::sin(Phi); |
---|
| 336 | cosphi = std::cos(Phi); |
---|
| 337 | |
---|
| 338 | px = -sintheta * cosphi; |
---|
| 339 | py = -sintheta * sinphi; |
---|
| 340 | pz = -costheta; |
---|
| 341 | |
---|
| 342 | // for volume and ponit source use mother or user defined co-ordinates |
---|
| 343 | // for plane and surface source user surface-normal or userdefined co-ordinates |
---|
| 344 | // |
---|
| 345 | G4double finx, finy, finz; |
---|
| 346 | if (posDist->SourcePosType == "Point" || posDist->SourcePosType == "Volume") { |
---|
| 347 | if (UserAngRef){ |
---|
| 348 | // Apply Rotation Matrix |
---|
| 349 | finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x()); |
---|
| 350 | finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y()); |
---|
| 351 | finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z()); |
---|
| 352 | } else { |
---|
| 353 | finx = px; |
---|
| 354 | finy = py; |
---|
| 355 | finz = pz; |
---|
| 356 | } |
---|
| 357 | } else { // for plane and surface source |
---|
| 358 | if (UserAngRef){ |
---|
| 359 | // Apply Rotation Matrix |
---|
| 360 | finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x()); |
---|
| 361 | finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y()); |
---|
| 362 | finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z()); |
---|
| 363 | } else { |
---|
| 364 | finx = (px*posDist->SideRefVec1.x()) + (py*posDist->SideRefVec2.x()) + (pz*posDist->SideRefVec3.x()); |
---|
| 365 | finy = (px*posDist->SideRefVec1.y()) + (py*posDist->SideRefVec2.y()) + (pz*posDist->SideRefVec3.y()); |
---|
| 366 | finz = (px*posDist->SideRefVec1.z()) + (py*posDist->SideRefVec2.z()) + (pz*posDist->SideRefVec3.z()); |
---|
| 367 | } |
---|
| 368 | } |
---|
| 369 | G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz)); |
---|
| 370 | finx = finx/ResMag; |
---|
| 371 | finy = finy/ResMag; |
---|
| 372 | finz = finz/ResMag; |
---|
| 373 | |
---|
| 374 | particle_momentum_direction.setX(finx); |
---|
| 375 | particle_momentum_direction.setY(finy); |
---|
| 376 | particle_momentum_direction.setZ(finz); |
---|
| 377 | |
---|
| 378 | // particle_momentum_direction now contains unit momentum vector. |
---|
| 379 | if(verbosityLevel >= 1) |
---|
| 380 | { |
---|
| 381 | G4cout << "Resultant cosine-law unit momentum vector " << particle_momentum_direction << G4endl; |
---|
| 382 | } |
---|
| 383 | } |
---|
| 384 | |
---|
| 385 | void G4SPSAngDistribution::GeneratePlanarFlux() |
---|
| 386 | { |
---|
| 387 | // particle_momentum_direction now contains unit momentum vector. |
---|
| 388 | // nothing need be done here as the m-directions have been set directly |
---|
| 389 | // under this option |
---|
| 390 | if(verbosityLevel >= 1) |
---|
| 391 | { |
---|
| 392 | G4cout << "Resultant Planar wave momentum vector " << particle_momentum_direction << G4endl; |
---|
| 393 | } |
---|
| 394 | } |
---|
| 395 | |
---|
| 396 | void G4SPSAngDistribution::GenerateUserDefFlux() |
---|
| 397 | { |
---|
| 398 | G4double rndm, px, py, pz, pmag; |
---|
| 399 | |
---|
| 400 | if(UserDistType == "NULL") |
---|
| 401 | G4cout << "Error: UserDistType undefined" << G4endl; |
---|
| 402 | else if(UserDistType == "theta") { |
---|
| 403 | Theta = 10.; |
---|
| 404 | while(Theta > MaxTheta || Theta < MinTheta) |
---|
| 405 | Theta = GenerateUserDefTheta(); |
---|
| 406 | Phi = 10.; |
---|
| 407 | while(Phi > MaxPhi || Phi < MinPhi) { |
---|
| 408 | rndm = angRndm->GenRandPhi(); |
---|
| 409 | Phi = twopi * rndm; |
---|
| 410 | } |
---|
| 411 | } |
---|
| 412 | else if(UserDistType == "phi") { |
---|
| 413 | Theta = 10.; |
---|
| 414 | while(Theta > MaxTheta || Theta < MinTheta) |
---|
| 415 | { |
---|
| 416 | rndm = angRndm->GenRandTheta(); |
---|
| 417 | Theta = std::acos(1. - (2. * rndm)); |
---|
| 418 | } |
---|
| 419 | Phi = 10.; |
---|
| 420 | while(Phi > MaxPhi || Phi < MinPhi) |
---|
| 421 | Phi = GenerateUserDefPhi(); |
---|
| 422 | } |
---|
| 423 | else if(UserDistType == "both") |
---|
| 424 | { |
---|
| 425 | Theta = 10.; |
---|
| 426 | while(Theta > MaxTheta || Theta < MinTheta) |
---|
| 427 | Theta = GenerateUserDefTheta(); |
---|
| 428 | Phi = 10.; |
---|
| 429 | while(Phi > MaxPhi || Phi < MinPhi) |
---|
| 430 | Phi = GenerateUserDefPhi(); |
---|
| 431 | } |
---|
| 432 | px = -std::sin(Theta) * std::cos(Phi); |
---|
| 433 | py = -std::sin(Theta) * std::sin(Phi); |
---|
| 434 | pz = -std::cos(Theta); |
---|
| 435 | |
---|
| 436 | pmag = std::sqrt((px*px) + (py*py) + (pz*pz)); |
---|
| 437 | |
---|
| 438 | if(!UserWRTSurface) { |
---|
| 439 | G4double finx, finy, finz; |
---|
| 440 | if (UserAngRef) { |
---|
| 441 | // Apply Rotation Matrix |
---|
| 442 | // x * AngRef1, y * AngRef2 and z * AngRef3 |
---|
| 443 | finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x()); |
---|
| 444 | finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y()); |
---|
| 445 | finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z()); |
---|
| 446 | } else { // use mother co-ordinates |
---|
| 447 | finx = px; |
---|
| 448 | finy = py; |
---|
| 449 | finz = pz; |
---|
| 450 | } |
---|
| 451 | G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz)); |
---|
| 452 | finx = finx/ResMag; |
---|
| 453 | finy = finy/ResMag; |
---|
| 454 | finz = finz/ResMag; |
---|
| 455 | |
---|
| 456 | particle_momentum_direction.setX(finx); |
---|
| 457 | particle_momentum_direction.setY(finy); |
---|
| 458 | particle_momentum_direction.setZ(finz); |
---|
| 459 | } |
---|
| 460 | else { // UserWRTSurface = true |
---|
| 461 | G4double pxh = px/pmag; |
---|
| 462 | G4double pyh = py/pmag; |
---|
| 463 | G4double pzh = pz/pmag; |
---|
| 464 | if(verbosityLevel > 1) { |
---|
| 465 | G4cout <<"SideRefVecs " <<posDist->SideRefVec1<<posDist->SideRefVec2<<posDist->SideRefVec3<<G4endl; |
---|
| 466 | G4cout <<"Raw Unit vector "<<pxh<<","<<pyh<<","<<pzh<<G4endl; |
---|
| 467 | } |
---|
| 468 | G4double resultx = (pxh*posDist->SideRefVec1.x()) + (pyh*posDist->SideRefVec2.x()) + |
---|
| 469 | (pzh*posDist->SideRefVec3.x()); |
---|
| 470 | |
---|
| 471 | G4double resulty = (pxh*posDist->SideRefVec1.y()) + (pyh*posDist->SideRefVec2.y()) + |
---|
| 472 | (pzh*posDist->SideRefVec3.y()); |
---|
| 473 | |
---|
| 474 | G4double resultz = (pxh*posDist->SideRefVec1.z()) + (pyh*posDist->SideRefVec2.z()) + |
---|
| 475 | (pzh*posDist->SideRefVec3.z()); |
---|
| 476 | |
---|
| 477 | G4double ResMag = std::sqrt((resultx*resultx) + (resulty*resulty) + (resultz*resultz)); |
---|
| 478 | resultx = resultx/ResMag; |
---|
| 479 | resulty = resulty/ResMag; |
---|
| 480 | resultz = resultz/ResMag; |
---|
| 481 | |
---|
| 482 | particle_momentum_direction.setX(resultx); |
---|
| 483 | particle_momentum_direction.setY(resulty); |
---|
| 484 | particle_momentum_direction.setZ(resultz); |
---|
| 485 | } |
---|
| 486 | |
---|
| 487 | // particle_momentum_direction now contains unit momentum vector. |
---|
| 488 | if(verbosityLevel > 0 ) |
---|
| 489 | { |
---|
| 490 | G4cout << "Final User Defined momentum vector " << particle_momentum_direction << G4endl; |
---|
| 491 | } |
---|
| 492 | } |
---|
| 493 | |
---|
| 494 | G4double G4SPSAngDistribution::GenerateUserDefTheta() |
---|
| 495 | { |
---|
| 496 | // Create cumulative histogram if not already done so. Then use RandFlat |
---|
| 497 | //::shoot to generate the output Theta value. |
---|
| 498 | if(UserDistType == "NULL" || UserDistType == "phi") |
---|
| 499 | { |
---|
| 500 | // No user defined theta distribution |
---|
| 501 | G4cout << "Error ***********************" << G4endl; |
---|
| 502 | G4cout << "UserDistType = " << UserDistType << G4endl; |
---|
| 503 | return (0.); |
---|
| 504 | } |
---|
| 505 | else |
---|
| 506 | { |
---|
| 507 | // UserDistType = theta or both and so a theta distribution |
---|
| 508 | // is defined. This should be integrated if not already done. |
---|
| 509 | if(IPDFThetaExist == false) |
---|
| 510 | { |
---|
| 511 | // IPDF has not been created, so create it |
---|
| 512 | G4double bins[1024],vals[1024], sum; |
---|
| 513 | G4int ii; |
---|
| 514 | G4int maxbin = G4int(UDefThetaH.GetVectorLength()); |
---|
| 515 | bins[0] = UDefThetaH.GetLowEdgeEnergy(size_t(0)); |
---|
| 516 | vals[0] = UDefThetaH(size_t(0)); |
---|
| 517 | sum = vals[0]; |
---|
| 518 | for(ii=1;ii<maxbin;ii++) |
---|
| 519 | { |
---|
| 520 | bins[ii] = UDefThetaH.GetLowEdgeEnergy(size_t(ii)); |
---|
| 521 | vals[ii] = UDefThetaH(size_t(ii)) + vals[ii-1]; |
---|
| 522 | sum = sum + UDefThetaH(size_t(ii)); |
---|
| 523 | } |
---|
| 524 | for(ii=0;ii<maxbin;ii++) |
---|
| 525 | { |
---|
| 526 | vals[ii] = vals[ii]/sum; |
---|
| 527 | IPDFThetaH.InsertValues(bins[ii], vals[ii]); |
---|
| 528 | } |
---|
| 529 | // Make IPDFThetaExist = true |
---|
| 530 | IPDFThetaExist = true; |
---|
| 531 | } |
---|
| 532 | // IPDF has been create so carry on |
---|
| 533 | G4double rndm = G4UniformRand(); |
---|
| 534 | return(IPDFThetaH.GetEnergy(rndm)); |
---|
| 535 | } |
---|
| 536 | } |
---|
| 537 | |
---|
| 538 | G4double G4SPSAngDistribution::GenerateUserDefPhi() |
---|
| 539 | { |
---|
| 540 | // Create cumulative histogram if not already done so. Then use RandFlat |
---|
| 541 | //::shoot to generate the output Theta value. |
---|
| 542 | |
---|
| 543 | if(UserDistType == "NULL" || UserDistType == "theta") |
---|
| 544 | { |
---|
| 545 | // No user defined phi distribution |
---|
| 546 | G4cout << "Error ***********************" << G4endl; |
---|
| 547 | G4cout << "UserDistType = " << UserDistType << G4endl; |
---|
| 548 | return(0.); |
---|
| 549 | } |
---|
| 550 | else |
---|
| 551 | { |
---|
| 552 | // UserDistType = phi or both and so a phi distribution |
---|
| 553 | // is defined. This should be integrated if not already done. |
---|
| 554 | if(IPDFPhiExist == false) |
---|
| 555 | { |
---|
| 556 | // IPDF has not been created, so create it |
---|
| 557 | G4double bins[1024],vals[1024], sum; |
---|
| 558 | G4int ii; |
---|
| 559 | G4int maxbin = G4int(UDefPhiH.GetVectorLength()); |
---|
| 560 | bins[0] = UDefPhiH.GetLowEdgeEnergy(size_t(0)); |
---|
| 561 | vals[0] = UDefPhiH(size_t(0)); |
---|
| 562 | sum = vals[0]; |
---|
| 563 | for(ii=1;ii<maxbin;ii++) |
---|
| 564 | { |
---|
| 565 | bins[ii] = UDefPhiH.GetLowEdgeEnergy(size_t(ii)); |
---|
| 566 | vals[ii] = UDefPhiH(size_t(ii)) + vals[ii-1]; |
---|
| 567 | sum = sum + UDefPhiH(size_t(ii)); |
---|
| 568 | } |
---|
| 569 | |
---|
| 570 | for(ii=0;ii<maxbin;ii++) |
---|
| 571 | { |
---|
| 572 | vals[ii] = vals[ii]/sum; |
---|
| 573 | IPDFPhiH.InsertValues(bins[ii], vals[ii]); |
---|
| 574 | } |
---|
| 575 | // Make IPDFPhiExist = true |
---|
| 576 | IPDFPhiExist = true; |
---|
| 577 | } |
---|
| 578 | // IPDF has been create so carry on |
---|
| 579 | G4double rndm = G4UniformRand(); |
---|
| 580 | return(IPDFPhiH.GetEnergy(rndm)); |
---|
| 581 | } |
---|
| 582 | } |
---|
| 583 | // |
---|
| 584 | void G4SPSAngDistribution::ReSetHist(G4String atype) |
---|
| 585 | { |
---|
| 586 | if (atype == "theta") { |
---|
| 587 | UDefThetaH = IPDFThetaH = ZeroPhysVector ; |
---|
| 588 | IPDFThetaExist = false ;} |
---|
| 589 | else if (atype == "phi"){ |
---|
| 590 | UDefPhiH = IPDFPhiH = ZeroPhysVector ; |
---|
| 591 | IPDFPhiExist = false ;} |
---|
| 592 | else { |
---|
| 593 | G4cout << "Error, histtype not accepted " << G4endl; |
---|
| 594 | } |
---|
| 595 | } |
---|
| 596 | |
---|
| 597 | |
---|
| 598 | G4ParticleMomentum G4SPSAngDistribution::GenerateOne() |
---|
| 599 | { |
---|
| 600 | // Angular stuff |
---|
| 601 | if(AngDistType == "iso") |
---|
| 602 | GenerateIsotropicFlux(); |
---|
| 603 | else if(AngDistType == "cos") |
---|
| 604 | GenerateCosineLawFlux(); |
---|
| 605 | else if(AngDistType == "planar") |
---|
| 606 | GeneratePlanarFlux(); |
---|
| 607 | else if(AngDistType == "beam1d" || AngDistType == "beam2d" ) |
---|
| 608 | GenerateBeamFlux(); |
---|
| 609 | else if(AngDistType == "user") |
---|
| 610 | GenerateUserDefFlux(); |
---|
| 611 | else if(AngDistType == "focused") |
---|
| 612 | GenerateFocusedFlux(); |
---|
| 613 | else |
---|
| 614 | G4cout << "Error: AngDistType has unusual value" << G4endl; |
---|
| 615 | return particle_momentum_direction; |
---|
| 616 | } |
---|
| 617 | |
---|
| 618 | |
---|
| 619 | |
---|
| 620 | |
---|
| 621 | |
---|
| 622 | |
---|
| 623 | |
---|
| 624 | |
---|
| 625 | |
---|