Geant4_10
G4SPSAngDistribution.hh
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.hh
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 // 26/10/2004 F Lei
42 // Added a "focused" option to allow all primary particles pointing to
43 // a user specified focusing point.
44 //
45 // Version 1.0, 05/02/2004, Fan Lei, Created.
46 // Based on the G4GeneralParticleSource class in Geant4 v6.0
47 //
49 //
50 // Class Description:
51 //
52 // To generate the direction of a primary vertex according to the defined distribution
53 //
55 //
56 // MEMBER FUNCTIONS
57 // ----------------
58 //
59 // G4SPSAngDistribution ()
60 // Constructor: Initializes variables
61 //
62 // ~G4SPSAngDistribution ()
63 // Destructor:
64 //
65 // void SetAngDistType(G4String)
66 // Used to set the type of angular distribution wanted. Arguments
67 // are iso, cos, beam and user for isotropic, cosine-law, beam and user-defined
68 // respectively.
69 //
70 // void DefineAngRefAxes(G4String, G4ThreeVector)
71 // DefineAngRefAxes is used in a similar way as SetPosRot to
72 // define vectors, one x' and one in the plane x'y', to create
73 // a rotated set of axes for the angular distribution.
74 //
75 // void SetMinTheta(G4double)
76 // Sets the minimum value for the angle theta.
77 //
78 // void SetMinPhi(G4double)
79 // Sets the minimum value for phi.
80 //
81 // void SetMaxTheta(G4double)
82 // Sets the maximum value for theta.
83 //
84 // void SetMaxPhi(G4double)
85 // Sets the maximum value for phi.
86 //
87 // void UserDefAngTheta(G4ThreeVector)
88 // This method allows the user to define a histogram in Theta.
89 //
90 // void UserDefAngPhi(G4ThreeVector)
91 // This method allows the user to define a histogram in phi.
92 //
93 // void GenerateIsotropicFlux()
94 // This method generates momentum vectors for particles according
95 // to an isotropic distribution.
96 //
97 // void GenerateCosineLawFlux()
98 // This method generates momentum vectors for particles according
99 // to a cosine-law distribution.
100 //
101 // void GenerateFocusedFlux()
102 // This method generates momentum vectors for particles pointing to
103 // an user specified focusing point.
104 //
105 // void GenerateUserDefFlux()
106 // Controls generation of momentum vectors according to user-defined
107 // distributions.
108 //
109 // G4double GenerateUserDefTheta()
110 // Generates the theta angle according to a user-defined distribution.
111 //
112 // G4double GenerateUserDefPhi()
113 // Generates phi according to a user-defined distribution.
114 //
115 // void SetBeamSigmaInAngR(G4double);
116 // Sets the sigma for 1D beam
117 //
118 // void SetBeamSigmaInAngX(G4double);
119 // Sets the first sigma for 2D beam
120 //
121 // void SetBeamSigmaInAngY(G4double);
122 // Sets the second sigma for 2D beam
123 //
124 // void SetUserWRTSurface(G4bool)
125 // Allows user to have user-defined spectra either with respect to the
126 // co-ordinate system (default) or with respect to the surface normal.
127 //
128 // void SetPosDistribution(G4SPSPosDistribution* a) {posDist = a; };
129 // Sets the required position generator, required for determining the cosine-law distribution
130 //
131 // void SetBiasRndm (G4SPSRandomGenerator* a)
132 // Sets the biased random number generator
133 //
134 // G4ThreeVector GenerateOne();
135 // Generate one random direction
136 //
137 // void ReSetHist(G4String);
138 // Re-sets the histogram for user defined distribution
139 //
140 // void SetVerbosity(G4int)
141 // Sets the verbosity level.
142 //
144 //
145 #ifndef G4SPSAngDistribution_h
146 #define G4SPSAngDistribution_h 1
147 
149 #include "G4DataInterpolation.hh"
150 #include "G4ParticleMomentum.hh"
151 
152 #include "G4SPSPosDistribution.hh"
153 #include "G4SPSRandomGenerator.hh"
154 
156 {
157 public:
160 
161  // Angular Distribution Methods
162  void SetAngDistType(G4String);
164  void SetMinTheta(G4double);
165  void SetMinPhi(G4double);
166  void SetMaxTheta(G4double);
167  void SetMaxPhi(G4double);
174  inline void SetParticleMomentumDirection
175  (G4ParticleMomentum aMomentumDirection)
176  { particle_momentum_direction = aMomentumDirection.unit(); }
179  //
181  void SetBiasRndm(G4SPSRandomGenerator* a) {angRndm = a;}
182  // method to re-set the histograms
183  void ReSetHist(G4String);
184  //
185  // Set the verbosity level.
186  void SetVerbosity(G4int a) {verbosityLevel = a; }
187  // some get methods
188  G4String GetDistType() { return AngDistType;}
189  G4double GetMinTheta() { return MinTheta; }
190  G4double GetMaxTheta() { return MaxTheta; }
191  G4double GetMinPhi() { return MinPhi; }
192  G4double GetMaxPhi() { return MaxPhi; }
193  //
195 
196 private:
197  // These methods generate the momentum vectors for the particles.
198  void GenerateFocusedFlux();
199  void GenerateIsotropicFlux();
200  void GenerateCosineLawFlux();
201  void GenerateBeamFlux();
202  void GeneratePlanarFlux();
203  void GenerateUserDefFlux();
204  G4double GenerateUserDefTheta();
205  G4double GenerateUserDefPhi();
206 
207 private:
208 
209  // Angular distribution variables.
210  G4String AngDistType; // String to hold Ang dist type iso, cos, user
211  G4ThreeVector AngRef1, AngRef2, AngRef3; // Reference axes for ang dist
212  G4double MinTheta, MaxTheta, MinPhi, MaxPhi; // min/max theta/phi
213  G4double DR,DX,DY ; // Standard deviations for beam divergence
214  G4double Theta, Phi; // Store these for use with DEBUG
215  G4ThreeVector FocusPoint ; // the focusing point in mother coordinates
216  G4bool IPDFThetaExist, IPDFPhiExist; // tell whether IPDF histos exist
217  G4PhysicsOrderedFreeVector UDefThetaH; // Theta histo data
218  G4PhysicsOrderedFreeVector IPDFThetaH; //Cumulative Theta histogram.
219  G4PhysicsOrderedFreeVector UDefPhiH; // Phi histo bins
220  G4PhysicsOrderedFreeVector IPDFPhiH; // Cumulative phi histogram.
221  G4String UserDistType; //String to hold user distributions
222  G4bool UserWRTSurface; // G4bool to tell whether user wants distribution wrt
223  // surface normals or co-ordinate system
224  G4bool UserAngRef; // Set to true when user defines a new coordinates
225  //
226  G4ParticleMomentum particle_momentum_direction;
227  //
228  G4SPSPosDistribution* posDist; // need it here for the cosine-law distri
229  G4SPSRandomGenerator* angRndm; // biased random generator
230 
231  // Verbosity
232  G4int verbosityLevel;
233  //
234  G4PhysicsOrderedFreeVector ZeroPhysVector ; // for re-set only
235 };
236 
237 #endif
tuple a
Definition: test.py:11
void SetBeamSigmaInAngR(G4double)
void SetBeamSigmaInAngY(G4double)
int G4int
Definition: G4Types.hh:78
void UserDefAngPhi(G4ThreeVector)
void UserDefAngTheta(G4ThreeVector)
G4ParticleMomentum GenerateOne()
bool G4bool
Definition: G4Types.hh:79
void SetBiasRndm(G4SPSRandomGenerator *a)
void DefineAngRefAxes(G4String, G4ThreeVector)
void SetPosDistribution(G4SPSPosDistribution *a)
Hep3Vector unit() const
void SetFocusPoint(G4ThreeVector)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
double G4double
Definition: G4Types.hh:76
void SetBeamSigmaInAngX(G4double)