Geant4  10.01
G4SPSAngDistribution.cc
Go to the documentation of this file.
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 //
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 //
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 //
47 //
48 
49 #include "G4SPSAngDistribution.hh"
50 
51 #include "Randomize.hh"
52 #include "G4PhysicalConstants.hh"
53 
54 
56  : Theta(0.), Phi(0.)
57 {
58  // Angular distribution Variables
59  G4ThreeVector zero;
61 
62  AngDistType = "planar";
66  MinTheta = 0.;
67  MaxTheta = pi;
68  MinPhi = 0.;
69  MaxPhi = twopi;
70  DR = 0.;
71  DX = 0.;
72  DY = 0.;
73  FocusPoint = G4ThreeVector(0., 0., 0.);
74  UserDistType = "NULL";
75  UserWRTSurface = true;
76  UserAngRef = false;
77  IPDFThetaExist = false;
78  IPDFPhiExist = false;
79  verbosityLevel = 0 ;
80 
82 }
83 
85 {
87 }
88 
89 //
91 {
92  G4AutoLock l(&mutex);
93  if(atype != "iso" && atype != "cos" && atype != "user" && atype != "planar"
94  && atype != "beam1d" && atype != "beam2d" && atype != "focused")
95  G4cout << "Error, distribution must be iso, cos, planar, beam1d, beam2d, focused or user" << G4endl;
96  else
97  AngDistType = atype;
98  if (AngDistType == "cos") MaxTheta = pi/2. ;
99  if (AngDistType == "user") {
101  IPDFThetaExist = false ;
103  IPDFPhiExist = false ;
104  }
105 }
106 
108 {
109  G4AutoLock l(&mutex);
110  if(refname == "angref1")
111  AngRef1 = ref.unit(); // x'
112  else if(refname == "angref2")
113  AngRef2 = ref.unit(); // vector in x'y' plane
114 
115  // User defines x' (AngRef1) and a vector in the x'y'
116  // plane (AngRef2). Then, AngRef1 x AngRef2 = AngRef3
117  // the z' vector. Then, AngRef3 x AngRef1 = AngRef2
118  // which will now be y'.
119 
120  AngRef3 = AngRef1.cross(AngRef2); // z'
121  AngRef2 = AngRef3.cross(AngRef1); // y'
122  UserAngRef = true ;
123  if(verbosityLevel == 2)
124  {
125  G4cout << "Angular distribution rotation axes " << AngRef1 << " " << AngRef2 << " " << AngRef3 << G4endl;
126  }
127 }
128 
130 {
131  G4AutoLock l(&mutex);
132  MinTheta = mint;
133 }
134 
136 {
137  G4AutoLock l(&mutex);
138  MinPhi = minp;
139 }
140 
142 {
143  G4AutoLock l(&mutex);
144  MaxTheta = maxt;
145 }
146 
148 {
149  G4AutoLock l(&mutex);
150  MaxPhi = maxp;
151 }
152 
154 {
155  G4AutoLock l(&mutex);
156  DR = r;
157 }
158 
160 {
161  G4AutoLock l(&mutex);
162  DX = r;
163 }
164 
166 {
167  G4AutoLock l(&mutex);
168  DY = r;
169 }
170 
172 {
173  G4AutoLock l(&mutex);
174  particle_momentum_direction = aMomentumDirection.unit();
175 
176 }
177 
179 {
180  G4AutoLock l(&mutex);
181  posDist = a;
182 }
183 
185 {
186  G4AutoLock l(&mutex);
187  angRndm = a;
188 }
189 
191 {
192  G4AutoLock l(&mutex);
193  verbosityLevel = a;
194 }
195 
196 
198 {
199  G4AutoLock l(&mutex);
200  if(UserDistType == "NULL") UserDistType = "theta";
201  if(UserDistType == "phi") UserDistType = "both";
202  G4double thi, val;
203  thi = input.x();
204  val = input.y();
205  if(verbosityLevel >= 1)
206  G4cout << "In UserDefAngTheta" << G4endl;
207  UDefThetaH.InsertValues(thi, val);
208 }
209 
215 
217 {
218  G4AutoLock l(&mutex);
219  if(UserDistType == "NULL") UserDistType = "phi";
220  if(UserDistType == "theta") UserDistType = "both";
221  G4double phhi, val;
222  phhi = input.x();
223  val = input.y();
224  if(verbosityLevel >= 1)
225  G4cout << "In UserDefAngPhi" << G4endl;
226  UDefPhiH.InsertValues(phhi, val);
227 }
228 
230 {
231  G4AutoLock l(&mutex);
232  FocusPoint = input;
233 }
234 
236 {
237  G4AutoLock l(&mutex);
238  // This is only applied in user mode?
239  // if UserWRTSurface = true then the user wants momenta with respect
240  // to the surface normals.
241  // When doing this theta has to be 0-90 only otherwise there will be
242  // errors, which currently are flagged anywhere.
243  UserWRTSurface = wrtSurf;
244 }
245 
247 {
248  G4AutoLock l(&mutex);
249  // if UserAngRef = true the angular distribution is defined wrt
250  // the user defined co-ordinates
251  UserAngRef = userang;
252 }
253 
255 {
256  G4double theta, phi;
257  G4double px, py, pz;
258  if (AngDistType == "beam1d")
259  {
260  theta = G4RandGauss::shoot(0.0,DR);
261  phi = twopi * G4UniformRand();
262  }
263  else
264  {
265  px = G4RandGauss::shoot(0.0,DX);
266  py = G4RandGauss::shoot(0.0,DY);
267  theta = std::sqrt (px*px + py*py);
268  if (theta != 0.) {
269  phi = std::acos(px/theta);
270  if ( py < 0.) phi = -phi;
271  } else {
272  phi = 0.0;
273  }
274  }
275  px = -std::sin(theta) * std::cos(phi);
276  py = -std::sin(theta) * std::sin(phi);
277  pz = -std::cos(theta);
278  G4double finx, finy, finz ;
279  finx = px, finy =py, finz =pz;
280  if (UserAngRef){
281  // Apply Angular Rotation Matrix
282  // x * AngRef1, y * AngRef2 and z * AngRef3
283  finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x());
284  finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y());
285  finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z());
286  G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz));
287  finx = finx/ResMag;
288  finy = finy/ResMag;
289  finz = finz/ResMag;
290  }
291  mom.setX(finx);
292  mom.setY(finy);
293  mom.setZ(finz);
294 
295  // particle_momentum_direction now holds unit momentum vector.
296  if(verbosityLevel >= 1)
297  G4cout << "Generating beam vector: " << mom << G4endl;
298 }
299 
301 {
302  mom = (FocusPoint - posDist->GetParticlePos()).unit();
303  //
304  // particle_momentum_direction now holds unit momentum vector.
305  if(verbosityLevel >= 1)
306  G4cout << "Generating focused vector: " << mom << G4endl;
307 }
308 
310 {
311  // generates isotropic flux.
312  // No vectors are needed.
313  G4double rndm, rndm2;
314  G4double px, py, pz;
315 
316  //
317  G4double sintheta, sinphi,costheta,cosphi;
318  rndm = angRndm->GenRandTheta();
319  costheta = std::cos(MinTheta) - rndm * (std::cos(MinTheta) - std::cos(MaxTheta));
320  sintheta = std::sqrt(1. - costheta*costheta);
321 
322  rndm2 = angRndm->GenRandPhi();
323  Phi = MinPhi + (MaxPhi - MinPhi) * rndm2;
324  sinphi = std::sin(Phi);
325  cosphi = std::cos(Phi);
326 
327  px = -sintheta * cosphi;
328  py = -sintheta * sinphi;
329  pz = -costheta;
330 
331  // for volume and ponit source use mother or user defined co-ordinates
332  // for plane and surface source user surface-normal or userdefined co-ordinates
333  //
334  G4double finx, finy, finz;
335  if (posDist->GetSourcePosType() == "Point" || posDist->GetSourcePosType() == "Volume") {
336  if (UserAngRef){
337  // Apply Rotation Matrix
338  // x * AngRef1, y * AngRef2 and z * AngRef3
339  finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x());
340  finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y());
341  finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z());
342  } else {
343  finx = px;
344  finy = py;
345  finz = pz;
346  }
347  } else { // for plane and surface source
348  if (UserAngRef){
349  // Apply Rotation Matrix
350  // x * AngRef1, y * AngRef2 and z * AngRef3
351  finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x());
352  finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y());
353  finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z());
354  } else {
355  finx = (px*posDist->GetSideRefVec1().x()) + (py*posDist->GetSideRefVec2().x()) + (pz*posDist->GetSideRefVec3().x());
356  finy = (px*posDist->GetSideRefVec1().y()) + (py*posDist->GetSideRefVec2().y()) + (pz*posDist->GetSideRefVec3().y());
357  finz = (px*posDist->GetSideRefVec1().z()) + (py*posDist->GetSideRefVec2().z()) + (pz*posDist->GetSideRefVec3().z());
358  }
359  }
360  G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz));
361  finx = finx/ResMag;
362  finy = finy/ResMag;
363  finz = finz/ResMag;
364 
365  mom.setX(finx);
366  mom.setY(finy);
367  mom.setZ(finz);
368 
369  // particle_momentum_direction now holds unit momentum vector.
370  if(verbosityLevel >= 1)
371  G4cout << "Generating isotropic vector: " << mom << G4endl;
372 }
373 
375 {
376  // Method to generate flux distributed with a cosine law
377  G4double px, py, pz;
378  G4double rndm, rndm2;
379  //
380  G4double sintheta, sinphi,costheta,cosphi;
381  rndm = angRndm->GenRandTheta();
382  sintheta = std::sqrt( rndm * (std::sin(MaxTheta)*std::sin(MaxTheta) - std::sin(MinTheta)*std::sin(MinTheta) )
383  +std::sin(MinTheta)*std::sin(MinTheta) );
384  costheta = std::sqrt(1. -sintheta*sintheta);
385 
386  rndm2 = angRndm->GenRandPhi();
387  Phi = MinPhi + (MaxPhi - MinPhi) * rndm2;
388  sinphi = std::sin(Phi);
389  cosphi = std::cos(Phi);
390 
391  px = -sintheta * cosphi;
392  py = -sintheta * sinphi;
393  pz = -costheta;
394 
395  // for volume and ponit source use mother or user defined co-ordinates
396  // for plane and surface source user surface-normal or userdefined co-ordinates
397  //
398  G4double finx, finy, finz;
399  if (posDist->GetSourcePosType() == "Point" || posDist->GetSourcePosType() == "Volume") {
400  if (UserAngRef){
401  // Apply Rotation Matrix
402  finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x());
403  finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y());
404  finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z());
405  } else {
406  finx = px;
407  finy = py;
408  finz = pz;
409  }
410  } else { // for plane and surface source
411  if (UserAngRef){
412  // Apply Rotation Matrix
413  finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x());
414  finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y());
415  finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z());
416  } else {
417  finx = (px*posDist->GetSideRefVec1().x()) + (py*posDist->GetSideRefVec2().x()) + (pz*posDist->GetSideRefVec3().x());
418  finy = (px*posDist->GetSideRefVec1().y()) + (py*posDist->GetSideRefVec2().y()) + (pz*posDist->GetSideRefVec3().y());
419  finz = (px*posDist->GetSideRefVec1().z()) + (py*posDist->GetSideRefVec2().z()) + (pz*posDist->GetSideRefVec3().z());
420  }
421  }
422  G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz));
423  finx = finx/ResMag;
424  finy = finy/ResMag;
425  finz = finz/ResMag;
426 
427  mom.setX(finx);
428  mom.setY(finy);
429  mom.setZ(finz);
430 
431  // particle_momentum_direction now contains unit momentum vector.
432  if(verbosityLevel >= 1)
433  {
434  G4cout << "Resultant cosine-law unit momentum vector " << mom << G4endl;
435  }
436 }
437 
439 {
440  // particle_momentum_direction now contains unit momentum vector.
441  // nothing need be done here as the m-directions have been set directly
442  // under this option
443  if(verbosityLevel >= 1)
444  {
445  G4cout << "Resultant Planar wave momentum vector " << mom << G4endl;
446  }
447 }
448 
450 {
451  G4double rndm, px, py, pz, pmag;
452 
453  if(UserDistType == "NULL")
454  G4cout << "Error: UserDistType undefined" << G4endl;
455  else if(UserDistType == "theta") {
456  Theta = 10.;
457  while(Theta > MaxTheta || Theta < MinTheta)
459  Phi = 10.;
460  while(Phi > MaxPhi || Phi < MinPhi) {
461  rndm = angRndm->GenRandPhi();
462  Phi = twopi * rndm;
463  }
464  }
465  else if(UserDistType == "phi") {
466  Theta = 10.;
467  while(Theta > MaxTheta || Theta < MinTheta)
468  {
469  rndm = angRndm->GenRandTheta();
470  Theta = std::acos(1. - (2. * rndm));
471  }
472  Phi = 10.;
473  while(Phi > MaxPhi || Phi < MinPhi)
475  }
476  else if(UserDistType == "both")
477  {
478  Theta = 10.;
479  while(Theta > MaxTheta || Theta < MinTheta)
481  Phi = 10.;
482  while(Phi > MaxPhi || Phi < MinPhi)
484  }
485  px = -std::sin(Theta) * std::cos(Phi);
486  py = -std::sin(Theta) * std::sin(Phi);
487  pz = -std::cos(Theta);
488 
489  pmag = std::sqrt((px*px) + (py*py) + (pz*pz));
490 
491  if(!UserWRTSurface) {
492  G4double finx, finy, finz;
493  if (UserAngRef) {
494  // Apply Rotation Matrix
495  // x * AngRef1, y * AngRef2 and z * AngRef3
496  finx = (px * AngRef1.x()) + (py * AngRef2.x()) + (pz * AngRef3.x());
497  finy = (px * AngRef1.y()) + (py * AngRef2.y()) + (pz * AngRef3.y());
498  finz = (px * AngRef1.z()) + (py * AngRef2.z()) + (pz * AngRef3.z());
499  } else { // use mother co-ordinates
500  finx = px;
501  finy = py;
502  finz = pz;
503  }
504  G4double ResMag = std::sqrt((finx*finx) + (finy*finy) + (finz*finz));
505  finx = finx/ResMag;
506  finy = finy/ResMag;
507  finz = finz/ResMag;
508 
509  mom.setX(finx);
510  mom.setY(finy);
511  mom.setZ(finz);
512  }
513  else { // UserWRTSurface = true
514  G4double pxh = px/pmag;
515  G4double pyh = py/pmag;
516  G4double pzh = pz/pmag;
517  if(verbosityLevel > 1) {
519  G4cout <<"Raw Unit vector "<<pxh<<","<<pyh<<","<<pzh<<G4endl;
520  }
521  G4double resultx = (pxh*posDist->GetSideRefVec1().x()) + (pyh*posDist->GetSideRefVec2().x()) +
522  (pzh*posDist->GetSideRefVec3().x());
523 
524  G4double resulty = (pxh*posDist->GetSideRefVec1().y()) + (pyh*posDist->GetSideRefVec2().y()) +
525  (pzh*posDist->GetSideRefVec3().y());
526 
527  G4double resultz = (pxh*posDist->GetSideRefVec1().z()) + (pyh*posDist->GetSideRefVec2().z()) +
528  (pzh*posDist->GetSideRefVec3().z());
529 
530  G4double ResMag = std::sqrt((resultx*resultx) + (resulty*resulty) + (resultz*resultz));
531  resultx = resultx/ResMag;
532  resulty = resulty/ResMag;
533  resultz = resultz/ResMag;
534 
535  mom.setX(resultx);
536  mom.setY(resulty);
537  mom.setZ(resultz);
538  }
539 
540  // particle_momentum_direction now contains unit momentum vector.
541  if(verbosityLevel > 0 )
542  {
543  G4cout << "Final User Defined momentum vector " << particle_momentum_direction << G4endl;
544  }
545 }
546 
548 {
549  // Create cumulative histogram if not already done so. Then use RandFlat
550  //::shoot to generate the output Theta value.
551  if(UserDistType == "NULL" || UserDistType == "phi")
552  {
553  // No user defined theta distribution
554  G4cout << "Error ***********************" << G4endl;
555  G4cout << "UserDistType = " << UserDistType << G4endl;
556  return (0.);
557  }
558  else
559  {
560  // UserDistType = theta or both and so a theta distribution
561  // is defined. This should be integrated if not already done.
562  G4AutoLock l(&mutex);
563  if(IPDFThetaExist == false)
564  {
565  // IPDF has not been created, so create it
566  G4double bins[1024],vals[1024], sum;
567  G4int ii;
569  bins[0] = UDefThetaH.GetLowEdgeEnergy(size_t(0));
570  vals[0] = UDefThetaH(size_t(0));
571  sum = vals[0];
572  for(ii=1;ii<maxbin;ii++)
573  {
574  bins[ii] = UDefThetaH.GetLowEdgeEnergy(size_t(ii));
575  vals[ii] = UDefThetaH(size_t(ii)) + vals[ii-1];
576  sum = sum + UDefThetaH(size_t(ii));
577  }
578  for(ii=0;ii<maxbin;ii++)
579  {
580  vals[ii] = vals[ii]/sum;
581  IPDFThetaH.InsertValues(bins[ii], vals[ii]);
582  }
583  // Make IPDFThetaExist = true
584  IPDFThetaExist = true;
585  }
586  l.unlock();
587  // IPDF has been create so carry on
588  G4double rndm = G4UniformRand();
589  return(IPDFThetaH.GetEnergy(rndm));
590  }
591 }
592 
594 {
595  // Create cumulative histogram if not already done so. Then use RandFlat
596  //::shoot to generate the output Theta value.
597 
598  if(UserDistType == "NULL" || UserDistType == "theta")
599  {
600  // No user defined phi distribution
601  G4cout << "Error ***********************" << G4endl;
602  G4cout << "UserDistType = " << UserDistType << G4endl;
603  return(0.);
604  }
605  else
606  {
607  // UserDistType = phi or both and so a phi distribution
608  // is defined. This should be integrated if not already done.
609  G4AutoLock l(&mutex);
610  if(IPDFPhiExist == false)
611  {
612  // IPDF has not been created, so create it
613  G4double bins[1024],vals[1024], sum;
614  G4int ii;
615  G4int maxbin = G4int(UDefPhiH.GetVectorLength());
616  bins[0] = UDefPhiH.GetLowEdgeEnergy(size_t(0));
617  vals[0] = UDefPhiH(size_t(0));
618  sum = vals[0];
619  for(ii=1;ii<maxbin;ii++)
620  {
621  bins[ii] = UDefPhiH.GetLowEdgeEnergy(size_t(ii));
622  vals[ii] = UDefPhiH(size_t(ii)) + vals[ii-1];
623  sum = sum + UDefPhiH(size_t(ii));
624  }
625  for(ii=0;ii<maxbin;ii++)
626  {
627  vals[ii] = vals[ii]/sum;
628  IPDFPhiH.InsertValues(bins[ii], vals[ii]);
629  }
630  // Make IPDFPhiExist = true
631  IPDFPhiExist = true;
632  }
633  l.unlock();
634  // IPDF has been create so carry on
635  G4double rndm = G4UniformRand();
636  return(IPDFPhiH.GetEnergy(rndm));
637  }
638 }
639 //
641 {
642  G4AutoLock l(&mutex);
643  if (atype == "theta") {
645  IPDFThetaExist = false ;}
646  else if (atype == "phi"){
648  IPDFPhiExist = false ;}
649  else {
650  G4cout << "Error, histtype not accepted " << G4endl;
651  }
652 }
653 
654 
656 {
657  //Local copy for thread safety
659  // Angular stuff
660  if(AngDistType == "iso")
661  GenerateIsotropicFlux(localM);
662  else if(AngDistType == "cos")
663  GenerateCosineLawFlux(localM);
664  else if(AngDistType == "planar")
665  GeneratePlanarFlux(localM);
666  else if(AngDistType == "beam1d" || AngDistType == "beam2d" )
667  GenerateBeamFlux(localM);
668  else if(AngDistType == "user")
669  GenerateUserDefFlux(localM);
670  else if(AngDistType == "focused")
671  GenerateFocusedFlux(localM);
672  else
673  G4cout << "Error: AngDistType has unusual value" << G4endl;
674  return localM;
675 }
676 
677 
678 
679 
680 
681 
682 
683 
684 
ThreeVector shoot(const G4int Ap, const G4int Af)
G4ThreeVector GetSideRefVec2() const
void GenerateBeamFlux(G4ParticleMomentum &outputMom)
CLHEP::Hep3Vector G4ThreeVector
G4PhysicsOrderedFreeVector UDefPhiH
void GenerateCosineLawFlux(G4ParticleMomentum &outputMom)
void SetBeamSigmaInAngR(G4double)
G4PhysicsOrderedFreeVector IPDFPhiH
const Hep3Vector HepYHat(0.0, 1.0, 0.0)
const G4double pi
#define G4MUTEXINIT(mutex)
Definition: G4Threading.hh:165
G4PhysicsOrderedFreeVector IPDFThetaH
void SetBeamSigmaInAngY(G4double)
void InsertValues(G4double energy, G4double value)
G4double a
Definition: TRTMaterials.hh:39
G4SPSPosDistribution * posDist
size_t GetVectorLength() const
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4SPSRandomGenerator * angRndm
const Hep3Vector HepXHat(1.0, 0.0, 0.0)
void UserDefAngPhi(G4ThreeVector)
void UserDefAngTheta(G4ThreeVector)
#define G4UniformRand()
Definition: Randomize.hh:95
G4GLOB_DLL std::ostream G4cout
G4PhysicsOrderedFreeVector UDefThetaH
G4ParticleMomentum particle_momentum_direction
G4ParticleMomentum GenerateOne()
bool G4bool
Definition: G4Types.hh:79
void SetBiasRndm(G4SPSRandomGenerator *a)
G4PhysicsOrderedFreeVector ZeroPhysVector
void DefineAngRefAxes(G4String, G4ThreeVector)
G4double GetEnergy(G4double aValue)
G4ThreeVector GetSideRefVec3() const
void GenerateFocusedFlux(G4ParticleMomentum &outputMom)
void SetPosDistribution(G4SPSPosDistribution *a)
void SetFocusPoint(G4ThreeVector)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
const Hep3Vector HepZHat(0.0, 0.0, 1.0)
G4ThreeVector GetSideRefVec1() const
void GenerateUserDefFlux(G4ParticleMomentum &outputMom)
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector GetParticlePos() const
void GenerateIsotropicFlux(G4ParticleMomentum &outputMom)
double G4double
Definition: G4Types.hh:76
G4ThreeVector G4ParticleMomentum
#define G4MUTEXDESTROY(mutex)
Definition: G4Threading.hh:166
void GeneratePlanarFlux(G4ParticleMomentum &outputMom)
void SetBeamSigmaInAngX(G4double)