GGS(GenericGEANT4Simulation)Software  2.6.0
 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...
 
- 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 227 of file GGSGunGeneratorAction.cpp.

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

52  {
53  _acceptanceFlag = val;
54  }
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 102 of file GGSGunGeneratorAction.h.

102  {
103  _energy = val;
104  }
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 63 of file GGSGunGeneratorAction.cpp.

63  {
64  if (particle == "C12") {
65  //G4ParticleTable* pTable = G4ParticleTable::GetParticleTable();
66  G4ParticleDefinition* c12 = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(6, 12, 0.);
67  ((G4ParticleGun*) _gunGenerator)->SetParticleDefinition(c12);
68  ((G4ParticleGun*) _gunGenerator)->SetParticleCharge(6 * eplus);
69  }
70  else if (particle == "Fe56") {
71  //G4ParticleTable* pTable = G4ParticleTable::GetParticleTable();
72  G4ParticleDefinition* fe56 = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(26, 56, 0.);
73  ((G4ParticleGun*) _gunGenerator)->SetParticleDefinition(fe56);
74  ((G4ParticleGun*) _gunGenerator)->SetParticleCharge(26 * eplus);
75  }
76  else
77  G4UImanager::GetUIpointer()->ApplyCommand(G4String("/gun/particle ").append(particle));
78 }
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 126 of file GGSGunGeneratorAction.h.

126  {
127  _maxEnergy = val;
128  }
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 327 of file GGSGunGeneratorAction.h.

327  {
328  _maxPhi = val;
329  _rndmMaxPhiFlag = true;
330  _SetRndmDirFlag();
331  }
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 183 of file GGSGunGeneratorAction.h.

183  {
184  _maxPosition = val;
185  _rndmMaxPositionFlag = true;
186  _SetRndmPosFlag();
187  _SetRndmFromSurface(1);
188  }
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 280 of file GGSGunGeneratorAction.h.

280  {
281  _maxThetaRand = (1. - cos(val)) / 2.;
282  _minCos2ThetaRand = cos(val);
283  _minCos2ThetaRand *= _minCos2ThetaRand;
284  _rndmMaxThetaFlag = true;
285  _SetRndmDirFlag();
286  }
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 114 of file GGSGunGeneratorAction.h.

114  {
115  _minEnergy = val;
116  }
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 313 of file GGSGunGeneratorAction.h.

313  {
314  _minPhi = val;
315  _rndmMinPhiFlag = true;
316  _SetRndmDirFlag();
317  }
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 167 of file GGSGunGeneratorAction.h.

167  {
168  _minPosition = val;
169  _rndmMinPositionFlag = true;
170  _SetRndmPosFlag();
171  _SetRndmFromSurface(1);
172  }
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 264 of file GGSGunGeneratorAction.h.

264  {
265  _minThetaRand = (1. - cos(val)) / 2.;
266  _maxCos2ThetaRand = cos(val);
267  _maxCos2ThetaRand *= _maxCos2ThetaRand;
268  _rndmMinThetaFlag = true;
269  _SetRndmDirFlag();
270  }
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 297 of file GGSGunGeneratorAction.h.

297  {
298  _sinPhi = sin(val);
299  _cosPhi = cos(val);
300  _rndmMinPhiFlag = _rndmMaxPhiFlag = false;
301  _SetRndmDirFlag();
302  _SetShootingDir();
303  }
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 150 of file GGSGunGeneratorAction.h.

150  {
151  _position = val;
152  // Reset random position generation
153  _rndmMinPositionFlag = _rndmMaxPositionFlag = false;
154  _SetRndmPosFlag();
155  _SetRndmFromSurface(0);
156  }
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 78 of file GGSGunGeneratorAction.h.

78  {
79  _rndmMaxEnergyFlag = val;
80  }
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 65 of file GGSGunGeneratorAction.h.

65  {
66  _rndmMinEnergyFlag = val;
67  }
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 139 of file GGSGunGeneratorAction.h.

139  {
140  _spectralIndex = val;
141  }
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 91 of file GGSGunGeneratorAction.h.

91  {
92  _spectralIndexFlag = val;
93  }
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 231 of file GGSGunGeneratorAction.h.

231  {
232  _rndmSphereCapExtRand = (1. - cos(val)) / 2.;
233  _SetRndmFromSurface(2);
234  }
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 215 of file GGSGunGeneratorAction.h.

215  {
216  _rndmSphereRadius = val.mag();
217  _rndmSphereCapAlpha = acos(val[2] / _rndmSphereRadius);
218  _rndmSphereCapBeta = atan2(val[1], val[0]);
219  _SetRndmFromSurface(2);
220  }
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 199 of file GGSGunGeneratorAction.h.

199  {
200  _rndmSphereCenter = val;
201  _SetRndmFromSurface(2);
202  }
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 245 of file GGSGunGeneratorAction.h.

245  {
246  _sinTheta = sin(val);
247  _cosTheta = cos(val);
248  _rndmMinThetaFlag = _rndmMaxThetaFlag = false;
249  _SetRndmDirFlag();
250  _SetShootingDir();
251 
252  }

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