GGS(GenericGEANT4Simulation)Software  2.7.1
 All Data Structures Namespaces Files Functions Variables Typedefs Macros
Public Member Functions
GGSGunGeneratorAction Class Reference

A single-particle generator. More...

#include <GGSGunGeneratorAction.h>

Inheritance diagram for GGSGunGeneratorAction:
Inheritance graph
[legend]
Collaboration diagram for GGSGunGeneratorAction:
Collaboration graph
[legend]

Public Member Functions

 GGSGunGeneratorAction ()
 Constructor.
 
 ~GGSGunGeneratorAction ()
 Destructor.
 
void GeneratePrimaries (G4Event *)
 GeneratePrimaries. More...
 
void SetAcceptanceFlag (G4bool val)
 Set the acceptance control. More...
 
void SetRndmMinEnergyFlag (G4bool val)
 Sets the min energy flag. More...
 
void SetRndmMaxEnergyFlag (G4bool val)
 Sets the max energy flag. More...
 
void SetSpectralIndexFlag (G4bool val)
 Sets the spectral index geeration. More...
 
void SetEnergy (G4double val)
 Sets the particle's kinetic energy. More...
 
void SetMinEnergy (G4double val)
 Sets the particle's min kinetic energy. More...
 
void SetMaxEnergy (G4double val)
 Sets the particle's max kinetic energy. More...
 
void SetSpectralIndex (G4double val)
 Sets the spectral index. More...
 
void SetPosition (const G4ThreeVector &val)
 Sets the particle's generation position. More...
 
void SetMinPosition (const G4ThreeVector &val)
 Sets the particle's minimum generation position. More...
 
void SetMaxPosition (const G4ThreeVector &val)
 Sets the particle's maximum generation position. More...
 
void SetSphereCenter (const G4ThreeVector &val)
 Sets the center of the spherical surface. More...
 
void SetSphereCapPosition (const G4ThreeVector &val)
 Sets the position of the spherical cap on the sphere. More...
 
void SetSphereCapExtension (G4double val)
 Sets the extension of the spherical cap. More...
 
void SetTheta (G4double val)
 Sets the value for theta. More...
 
void SetMinTheta (G4double val)
 Sets the minimum random value for theta. More...
 
void SetMaxTheta (G4double val)
 Sets the maximum random value for theta. More...
 
void SetPhi (G4double val)
 Sets the value for phi. More...
 
void SetMinPhi (G4double val)
 Sets the minimum value for phi. More...
 
void SetMaxPhi (G4double val)
 Sets the maximum value for phi. More...
 
void SetGunParticle (const G4String &particle)
 Sets the particle to be used by G4ParticleGun. More...
 
GGSParameters GetParameters () const
 Gets the generation parameters. More...
 
- Public Member Functions inherited from GGSGeneratorAction
 GGSGeneratorAction ()
 Constructor.
 
 ~GGSGeneratorAction ()
 Destructor.
 
unsigned int GetNDiscarded () const
 Returns the number of discarded events for the current event generation. More...
 

Additional Inherited Members

- Protected Attributes inherited from GGSGeneratorAction
unsigned long int _nDiscarded
 

Detailed Description

A single-particle generator.

This is a wrapper for the G4ParticleGun generator, with some custom features.

Definition at line 26 of file GGSGunGeneratorAction.h.

Member Function Documentation

void GGSGunGeneratorAction::GeneratePrimaries ( G4Event *  anEvent)

GeneratePrimaries.

This class defines the parameters of the projectile primary particles It is called at the beginning of each event.

Definition at line 208 of file GGSGunGeneratorAction.cpp.

208  {
209 
210  static const std::string routineName("GGSGunGeneratorAction::GeneratePrimaries");
211  static const G4ThreeVector xDirection(1, 0, 0);
212  static const G4ThreeVector yDirection(0, 1, 0);
213 
214  _nDiscarded = 0;
215 
216  // Initialize shooting parameters with default values
217  G4ThreeVector position(_position);
218  G4ThreeVector direction(_direction);
219  G4double energy(_energy);
220 
221  bool acceptanceCheck = true;
222  do { // Acceptance check loop
223 
224  // ++++++++++++++ Randomize direction ++++++++++++++
225  // The random direction computed here is always referred to the shooting frame.
226  // Since the shooting frame may depend on shooting point (eg., shooting
227  // from a sphere) the transformation from shooting frame to absolute
228  // frame is done after random position have been chosen.
229  if (_rndmDirFlag) {
230  G4double cosPhi(_cosPhi);
231  G4double sinPhi(_sinPhi);
232  if (_rndmMinPhiFlag && _rndmMaxPhiFlag) {
233  double phi = CLHEP::RandFlat::shoot(_minPhi, _maxPhi); // (0,2pi) or (-2pi,2pi) if shooting along -Z or -X
234  cosPhi = cos(phi);
235  sinPhi = sin(phi);
236  }
237 
238  G4double sinTheta(_sinTheta);
239  G4double cosTheta(_cosTheta);
240  if (_rndmMinThetaFlag && _rndmMaxThetaFlag) {
241  if (_rndmFromFlat || _rndmFromSphere) {
242  // To generate an isotropic flux from a flat surface, the random angle must be uniform in cos^2 theta
243  float cos2Theta = CLHEP::RandFlat::shoot(_minCos2ThetaRand, _maxCos2ThetaRand);
244  cosTheta = sqrt(cos2Theta);
245  sinTheta = sqrt(1 - cos2Theta);
246  } else {
247  // This will sample a uniform theta distribution on the sphere. This formula comes from the inverse CDF
248  // method. See eg. http://en.wikibooks.org/wiki/Mathematica/Uniform_Spherical_Distribution
249  // The form reported here is slightly different but equivalent (it can be easily seen using the
250  // identity (sin x/2)^2 = (1 - cos x)/2 ).
251  cosTheta = 1 - 2 * CLHEP::RandFlat::shoot(_minThetaRand, _maxThetaRand);
252  sinTheta = sqrt(1 - cosTheta * cosTheta);
253  }
254  }
255 
256  // Compute position std::vector and direction versor (Z is shoot axis in the shooting frame)
257  direction[0] = -sinTheta * cosPhi; // momentum direction (unitary vectors)
258  direction[1] = -sinTheta * sinPhi;
259  direction[2] = -cosTheta;
260  }
261 
262  // ++++++++++++++ Randomize position ++++++++++++++
263  if (_rndmPosFlag) {
264  if (_rndmMinPositionFlag && _rndmMaxPositionFlag) {
265  // This works for all "flat" domains (point, line, plane, parallelepiped)
266  position[0] = CLHEP::RandFlat::shoot(_minPosition(0), _maxPosition(0));
267  position[1] = CLHEP::RandFlat::shoot(_minPosition(1), _maxPosition(1));
268  position[2] = CLHEP::RandFlat::shoot(_minPosition(2), _maxPosition(2));
269  if (_rndmFromFlat) { // Transform shoot direction to absolute reference frame
270  // Transform direction if shoot axis is not the Z axis.
271  // 3rd coordinate is relative to shooting axis so it is x when shooting along X,
272  // so the coordinate transformation to go to absolute reference frame is:
273  // (y, z, x) -> (x, y, z) when shooting along X
274  // (z, x, y) -> (x, y, z) when shooting along X
275  if (_shootAxis == xDirection) {
276  double buf = direction[0];
277  direction[0] = direction[2];
278  direction[2] = direction[1];
279  direction[1] = buf;
280  } else if (_shootAxis == yDirection) {
281  double buf = direction[0];
282  direction[0] = direction[1];
283  direction[1] = direction[2];
284  direction[2] = buf;
285  }
286  }
287  } else if (_rndmFromSphere) {
288  // Randomize position on the spherical cap (in the reference frame of the generation sphere).
289  // Gamma and delta are the polar and azimuth angle of the position vector w.r.t. the
290  // sphere center.
291  float cosGamma = 1 - 2 * CLHEP::RandFlat::shoot(0., _rndmSphereCapExtRand);
292  float delta = CLHEP::RandFlat::shoot(0., 2 * CLHEP::pi);
293  float sinGamma = sqrt(1 - cosGamma * cosGamma);
294  float sinDelta = sin(delta);
295  float cosDelta = cos(delta);
296 
297  // Position vector in cap frame
298  position[0] = _rndmSphereRadius * sinGamma * cosDelta;
299  position[1] = _rndmSphereRadius * sinGamma * sinDelta;
300  position[2] = _rndmSphereRadius * cosGamma;
301 
302  // Position vector in absolute (world) frame
303  position *= _rndmSphereCapRotMat;
304  position += _rndmSphereCenter;
305 
306  // Direction vector in absolute (world) frame
307  // At this point, the direction vector is built so that the distribution of its polar angle with
308  // respect to the normal incidence is flat in cos^2 (isotropic flux). The vector coordinates have
309  // been constructed in the reference frame of the generation point on the sphere, so they must be
310  // rotated to the reference frame of the cap center and then to the world reference frame, in this
311  // order.
312  direction *= _ChangeOfBasisMatrix(cosGamma, sinGamma, cosDelta, sinDelta); // Direction
313  direction *= _rndmSphereCapRotMat;
314  }
315  }
316 
317  // ++++++++++++++ Check acceptance ++++++++++++++
318  if (_acceptanceFlag && (_rndmPosFlag || _rndmDirFlag)) {
319  GGSDetectorConstruction *detConst =
320  (GGSDetectorConstruction *)(GGSRunManager::GetRunManager()->GetUserDetectorConstruction());
321  acceptanceCheck = detConst->GetGeometry()->IsInsideAcceptance(position, direction);
322  if (!acceptanceCheck)
323  _nDiscarded++;
324  }
325 
326  } while (!acceptanceCheck);
327 
328  // ++++++++++++++ Randomize energy ++++++++++++++
329  if (_rndmMinEnergyFlag && _rndmMaxEnergyFlag) {
330  if (_spectralIndexFlag) {
331  energy = _GenSpectrumPowerLaw(_minEnergy, _maxEnergy, _spectralIndex);
332  } else {
333  energy = CLHEP::RandFlat::shoot(_minEnergy, _maxEnergy);
334  }
335  }
336 
337  COUT(DEEPDEB) << "Generated primary:" << ENDL;
338  CCOUT(DEEPDEB) << " PDG code: " << _gunGenerator->GetParticleDefinition()->GetPDGEncoding() << ENDL;
339  CCOUT(DEEPDEB) << " position: " << position << ENDL;
340  CCOUT(DEEPDEB) << " direction: " << direction << ENDL;
341  CCOUT(DEEPDEB) << " energy: " << energy << ENDL;
342 
343  //
344  // Set values for primary particle and shoot
345  //
346  _gunGenerator->SetParticlePosition(position);
347  _gunGenerator->SetParticleMomentumDirection(direction);
348  _gunGenerator->SetParticleEnergy(energy);
349 
350  // Shoot
351  _gunGenerator->GeneratePrimaryVertex(anEvent);
352 }
#define ENDL
Definition: GGSSmartLog.h:105
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:76
Class for GGS detector construction.
const GGSVGeometryConstruction * GetGeometry() const
Returns a pointer to the geometry.
#define CCOUT(level)
Smart log utility which prints no header at the beginning of the line.
Definition: GGSSmartLog.h:100
static GGSRunManager * GetRunManager()
Static getter function the run manager.
virtual bool IsInsideAcceptance(const G4ThreeVector &, const G4ThreeVector &) const
Checks if a particle is inside the instrument&#39;s acceptance.
GGSParameters GGSGunGeneratorAction::GetParameters ( ) const
virtual

Gets the generation parameters.

The parameters exported by the GGSGunGeneratorAction are:

General

Name Type Description Exported
generator string The name of the generator "gun" Always
acceptanceCheck bool Flag for acceptance check Always
particlePdgCode int the PDG code of the generated particle Always

Generation position

Name Type Description Exported
genPosition string Generation position: "point", "plane", "box", "sphere" Always
minPosition vector<double> Coordinates of the minimum generation point [cm] When genPosition is "plane" or "box"
maxPosition vector<double> Coordinates of the maximum generation point [cm] When genPosition is "plane" or "box"
sphereCenter vector<double> Center of the generation sphere cap [cm] When genPosition is "sphere"
sphereCapPosition vector<double> Position of the sphere cap w.r.t. the sphere center [cm] When genPosition is "sphere"
sphereCapExtensiondouble Angular extension of the sphere [rad] When genPosition is "sphere"
position vector<double> Coordinates of the generation point [cm] When genPosition is "point"

Generation direction

Name Type Description Exported
genDirection string Generation direction: "fixed", "uniform", "isotropic" Always
minPhi double Minimum azimuth angle [rad] When genDirection is "uniform" or "isotropic"
maxPhi double Maximum azimuth angle [rad] When genDirection is "uniform" or "isotropic"
minTheta double Minimum polar angle [rad] When genDirection is "uniform" or "isotropic"
maxTheta double Maximum polar angle [rad] When genDirection is "uniform" or "isotropic"
phi double Azimuth angle [rad] When minPhi and maxPhi are not exported
theta double Polar angle [rad] When minTheta and maxTheta are not exported

Notes: for generations from surfaces (i.e. genPosition == "point" or genPosition == "sphere") the polar and azimuth angles are always referred to the normal direction to the surface in the generation point.

Generation energy

Name Type Description Exported
genSpectrum string Generation spectrum: "monochromatic", "flat", "powerLaw" Always
minEnergy double Minimum kinetic energy [GeV] When genSpectrum is "flat" or "powerLaw"
maxEnergy double Maximum kinetic energy [GeV] When genSpectrum is "flat" or "powerLaw"
spectralIndex double Spectral index When genSpectrum is "powerLaw"
energy double Kinetic energy [GeV] When genSpectrum is "monochromatic"
Returns
The generation parameters object;

Reimplemented from GGSGeneratorAction.

Definition at line 354 of file GGSGunGeneratorAction.cpp.

354  {
355  GGSParameters params;
356  params.SetParam("generator", std::string("gun"));
357 
358  params.SetParam("acceptanceCheck", _acceptanceFlag);
359 
360  // Particle
361  params.SetParam("particlePdgCode", _gunGenerator->GetParticleDefinition()->GetPDGEncoding());
362 
363  // Position
364  if (_rndmPosFlag) {
365  if (_rndmMinPositionFlag && _rndmMaxPositionFlag) {
366  params.SetParam("minPosition", std::vector<double>{_minPosition[0] / CLHEP::cm, _minPosition[1] / CLHEP::cm,
367  _minPosition[2] / CLHEP::cm});
368  params.SetParam("maxPosition", std::vector<double>{_maxPosition[0] / CLHEP::cm, _maxPosition[1] / CLHEP::cm,
369  _maxPosition[2] / CLHEP::cm});
370  if (_rndmFromFlat) {
371  params.SetParam("genPosition", std::string("plane"));
372  } else {
373  params.SetParam("genPosition", std::string("box"));
374  }
375  } else if (_rndmFromSphere) {
376  params.SetParam("genPosition", std::string("sphere"));
377  params.SetParam("sphereCenter",
378  std::vector<double>{_rndmSphereCenter[0] / CLHEP::cm, _rndmSphereCenter[1] / CLHEP::cm,
379  _rndmSphereCenter[2] / CLHEP::cm});
380  params.SetParam("sphereCapPosition",
381  std::vector<double>{_rndmSphereCapPosition[0] / CLHEP::cm, _rndmSphereCapPosition[1] / CLHEP::cm,
382  _rndmSphereCapPosition[2] / CLHEP::cm});
383  params.SetParam("sphereCapExtension", acos(-2. * _rndmSphereCapExtRand + 1) / CLHEP::rad);
384  }
385  } else {
386  params.SetParam("genPosition", std::string("point"));
387  params.SetParam("position",
388  std::vector<double>{_position[0] / CLHEP::cm, _position[1] / CLHEP::cm, _position[2] / CLHEP::cm});
389  }
390 
391  // Direction
392  if (_rndmDirFlag) {
393  bool isotropic = false;
394  if (_rndmMinPhiFlag && _rndmMaxPhiFlag) {
395  params.SetParam("minPhi", _minPhi / CLHEP::rad);
396  params.SetParam("maxPhi", _maxPhi / CLHEP::rad);
397  }
398  if (_rndmMinThetaFlag && _rndmMaxThetaFlag) {
399  params.SetParam("minTheta", acos(-2. * _minThetaRand + 1) / CLHEP::rad);
400  params.SetParam("maxTheta", acos(-2. * _maxThetaRand + 1) / CLHEP::rad);
401  if (_rndmFromFlat || _rndmFromSphere) {
402  isotropic = true;
403  }
404  }
405  if (isotropic) {
406  params.SetParam("genDirection", std::string("isotropic"));
407  } else {
408  params.SetParam("genDirection", std::string("uniform"));
409  }
410  } else {
411  params.SetParam("genDirection", std::string("fixed"));
412  params.SetParam("theta", acos(_cosTheta) / CLHEP::rad);
413  params.SetParam("phi", acos(_cosPhi) / CLHEP::rad);
414  }
415 
416  // Energy
417  if (_rndmMinEnergyFlag && _rndmMaxEnergyFlag) {
418  params.SetParam("minEnergy", _minEnergy / CLHEP::GeV);
419  params.SetParam("maxEnergy", _maxEnergy / CLHEP::GeV);
420  if (_spectralIndexFlag) {
421  params.SetParam("genSpectrum", std::string("powerLaw"));
422  params.SetParam("spectralIndex", _spectralIndex);
423  } else {
424  params.SetParam("genSpectrum", std::string("flat"));
425  }
426  } else {
427  params.SetParam("genSpectrum", std::string("monochromatic"));
428  params.SetParam("energy", _energy / CLHEP::GeV);
429  }
430 
431  return params;
432 }
void SetParam(const std::string &name, T value)
Sets a parameter.
Class for storing parameters.
Definition: GGSParameters.h:29
void GGSGunGeneratorAction::SetAcceptanceFlag ( G4bool  val)
inline

Set the acceptance control.

Acceptance must defined by the implementation of GGSVGeometryConstruction::IsInsideAcceptance in the current geometry. Events discarded by the acceptance check can be retrieved by means of the GetNDiscarded method defined in the mother class.

Parameters
valIf true, acceptance check will be performed for random events.
See Also
GGSGeneratorAction::GetNDiscarded

Definition at line 51 of file GGSGunGeneratorAction.h.

51 { _acceptanceFlag = val; }
void GGSGunGeneratorAction::SetEnergy ( G4double  val)
inline

Sets the particle's kinetic energy.

This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe particle's energy.

Definition at line 93 of file GGSGunGeneratorAction.h.

93 { _energy = val; }
void GGSGunGeneratorAction::SetGunParticle ( const G4String &  particle)

Sets the particle to be used by G4ParticleGun.

Parameters
particleThe particle to shoot. A std::string value like "e+" or "proton".

Definition at line 62 of file GGSGunGeneratorAction.cpp.

62  {
63  if (particle == "C12") {
64  // G4ParticleTable* pTable = G4ParticleTable::GetParticleTable();
65  G4ParticleDefinition *c12 = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(6, 12, 0.);
66  ((G4ParticleGun *)_gunGenerator)->SetParticleDefinition(c12);
67  ((G4ParticleGun *)_gunGenerator)->SetParticleCharge(6 * eplus);
68  } else if (particle == "Fe56") {
69  // G4ParticleTable* pTable = G4ParticleTable::GetParticleTable();
70  G4ParticleDefinition *fe56 = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(26, 56, 0.);
71  ((G4ParticleGun *)_gunGenerator)->SetParticleDefinition(fe56);
72  ((G4ParticleGun *)_gunGenerator)->SetParticleCharge(26 * eplus);
73  } else
74  G4UImanager::GetUIpointer()->ApplyCommand(G4String("/gun/particle ").append(particle));
75 }
void GGSGunGeneratorAction::SetMaxEnergy ( G4double  val)
inline

Sets the particle's max kinetic energy.

Sets the maximum value of randomized energy spectrum. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe particle's max energy.

Definition at line 113 of file GGSGunGeneratorAction.h.

113 { _maxEnergy = val; }
void GGSGunGeneratorAction::SetMaxPhi ( G4double  val)
inline

Sets the maximum value for phi.

For isotropic generation, this sets the maximum random value for phi. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe maximum random value for phi.

Definition at line 310 of file GGSGunGeneratorAction.h.

310  {
311  _maxPhi = val;
312  _rndmMaxPhiFlag = true;
313  _SetRndmDirFlag();
314  }
void GGSGunGeneratorAction::SetMaxPosition ( const G4ThreeVector &  val)
inline

Sets the particle's maximum generation position.

Random position generation is done inside a volume defined by specifying the "minimum" and "maximum" coordinates vectors. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe particle's max position.

Definition at line 166 of file GGSGunGeneratorAction.h.

166  {
167  _maxPosition = val;
168  _rndmMaxPositionFlag = true;
169  _SetRndmPosFlag();
170  _SetRndmFromSurface(1);
171  }
void GGSGunGeneratorAction::SetMaxTheta ( G4double  val)
inline

Sets the maximum random value for theta.

For isotropic generation, this sets the maximum random value for theta. From the theta value passed as parameter it computes and store the lower limit of the uniform distribution from which a random number will be extracted and subsequently transformed into theta by means of an inverse CDF technique.

Parameters
valThe maximum random value for theta.

Definition at line 263 of file GGSGunGeneratorAction.h.

263  {
264  _maxThetaRand = (1. - cos(val)) / 2.;
265  _minCos2ThetaRand = cos(val);
266  _minCos2ThetaRand *= _minCos2ThetaRand;
267  _rndmMaxThetaFlag = true;
268  _SetRndmDirFlag();
269  }
void GGSGunGeneratorAction::SetMinEnergy ( G4double  val)
inline

Sets the particle's min kinetic energy.

Sets the minimum value of randomized energy spectrum. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe particle's min energy.

Definition at line 103 of file GGSGunGeneratorAction.h.

103 { _minEnergy = val; }
void GGSGunGeneratorAction::SetMinPhi ( G4double  val)
inline

Sets the minimum value for phi.

For isotropic generation, this sets the minimum random value for phi. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe minimum random value for phi.

Definition at line 296 of file GGSGunGeneratorAction.h.

296  {
297  _minPhi = val;
298  _rndmMinPhiFlag = true;
299  _SetRndmDirFlag();
300  }
void GGSGunGeneratorAction::SetMinPosition ( const G4ThreeVector &  val)
inline

Sets the particle's minimum generation position.

Random position generation is done inside a volume defined by specifying the "minimum" and "maximum" coordinates vectors. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe particle's min position.

Definition at line 150 of file GGSGunGeneratorAction.h.

150  {
151  _minPosition = val;
152  _rndmMinPositionFlag = true;
153  _SetRndmPosFlag();
154  _SetRndmFromSurface(1);
155  }
void GGSGunGeneratorAction::SetMinTheta ( G4double  val)
inline

Sets the minimum random value for theta.

For isotropic generation, this sets the minimum random value for theta. From the theta value passed as parameter it computes and store the lower limit of the uniform distribution from which a random number will be extracted and subsequently transformed into theta by means of an inverse CDF technique. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe minimum random value for theta.

Definition at line 247 of file GGSGunGeneratorAction.h.

247  {
248  _minThetaRand = (1. - cos(val)) / 2.;
249  _maxCos2ThetaRand = cos(val);
250  _maxCos2ThetaRand *= _maxCos2ThetaRand;
251  _rndmMinThetaFlag = true;
252  _SetRndmDirFlag();
253  }
void GGSGunGeneratorAction::SetPhi ( G4double  val)
inline

Sets the value for phi.

Phi is the azimuth angle, defined as zero for particles traveling along the X axis towards negative direction. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe value for phi.

Definition at line 280 of file GGSGunGeneratorAction.h.

280  {
281  _sinPhi = sin(val);
282  _cosPhi = cos(val);
283  _rndmMinPhiFlag = _rndmMaxPhiFlag = false;
284  _SetRndmDirFlag();
285  _SetShootingDir();
286  }
void GGSGunGeneratorAction::SetPosition ( const G4ThreeVector &  val)
inline

Sets the particle's generation position.

This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe particle's position.

Definition at line 133 of file GGSGunGeneratorAction.h.

133  {
134  _position = val;
135  // Reset random position generation
136  _rndmMinPositionFlag = _rndmMaxPositionFlag = false;
137  _SetRndmPosFlag();
138  _SetRndmFromSurface(0);
139  }
void GGSGunGeneratorAction::SetRndmMaxEnergyFlag ( G4bool  val)
inline

Sets the max energy flag.

If both min energy and max energy flags are set to true, energy will be randomized. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe flag value.

Definition at line 73 of file GGSGunGeneratorAction.h.

73 { _rndmMaxEnergyFlag = val; }
void GGSGunGeneratorAction::SetRndmMinEnergyFlag ( G4bool  val)
inline

Sets the min energy flag.

If both min energy and max energy flags are set to true, energy will be randomized. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe flag value.

Definition at line 62 of file GGSGunGeneratorAction.h.

62 { _rndmMinEnergyFlag = val; }
void GGSGunGeneratorAction::SetSpectralIndex ( G4double  val)
inline

Sets the spectral index.

If spectral index flag is set, this method specifies the value of the index to generate a power law spectrum. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe spectral index value.

Definition at line 124 of file GGSGunGeneratorAction.h.

124 { _spectralIndex = val; }
void GGSGunGeneratorAction::SetSpectralIndexFlag ( G4bool  val)
inline

Sets the spectral index geeration.

If set to true, the eventually randomized energies will be extracted to generate a power law spectrum. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe flag value.

Definition at line 84 of file GGSGunGeneratorAction.h.

84 { _spectralIndexFlag = val; }
void GGSGunGeneratorAction::SetSphereCapExtension ( G4double  val)
inline

Sets the extension of the spherical cap.

The angular width of the spherical cap is set by this method. By "angular width" we mean the angle between the cap position std::vector (set by SetSphereCapPosition) and the position std::vector of a point on the cap's border (both the vectors are relative to the sphere center).

Parameters
valThe angular extension of the spherical cap.

Definition at line 215 of file GGSGunGeneratorAction.h.

215  {
216  _rndmSphereCapExtRand = (1. - cos(val)) / 2.;
217  _SetRndmFromSurface(2);
218  }
void GGSGunGeneratorAction::SetSphereCapPosition ( const G4ThreeVector &  val)
inline

Sets the position of the spherical cap on the sphere.

The position of the spherical cap from which the particles will be randomly generated can be set by means of this method. The position std::vector is relative to the sphere center and not to the absolute frame. The width of the cap can be set using SetSphereCapExtension. The modulus of the specified position std::vector of the spherical cap will be used as the sphere radius.

Parameters
val

Definition at line 198 of file GGSGunGeneratorAction.h.

198  {
199  _rndmSphereCapPosition = val;
200  _rndmSphereRadius = val.mag();
201  _rndmSphereCapAlpha = acos(val[2] / _rndmSphereRadius);
202  _rndmSphereCapBeta = atan2(val[1], val[0]);
203  _SetRndmFromSurface(2);
204  }
void GGSGunGeneratorAction::SetSphereCenter ( const G4ThreeVector &  val)
inline

Sets the center of the spherical surface.

This method sets the center of the sphere from which particles will be randomly generated. Actually, the generation is not done on the whole sphere but on a cap.

See Also
SetSphereCapPosition
Parameters
valThe center of the sphere.

Definition at line 182 of file GGSGunGeneratorAction.h.

182  {
183  _rndmSphereCenter = val;
184  _SetRndmFromSurface(2);
185  }
void GGSGunGeneratorAction::SetTheta ( G4double  val)
inline

Sets the value for theta.

Theta is the altitude angle, defined as zero for particles traveling along the Z axis towards negative direction. This function is called by the messenger. Has effect only on gun generator.

Parameters
valThe value for theta.

Definition at line 229 of file GGSGunGeneratorAction.h.

229  {
230  _sinTheta = sin(val);
231  _cosTheta = cos(val);
232  _rndmMinThetaFlag = _rndmMaxThetaFlag = false;
233  _SetRndmDirFlag();
234  _SetShootingDir();
235  }

The documentation for this class was generated from the following files: