Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GammaRayTelPrimaryGeneratorAction Class Reference

#include <GammaRayTelPrimaryGeneratorAction.hh>

Inheritance diagram for GammaRayTelPrimaryGeneratorAction:
Collaboration diagram for GammaRayTelPrimaryGeneratorAction:

Public Member Functions

 GammaRayTelPrimaryGeneratorAction ()
 
 ~GammaRayTelPrimaryGeneratorAction ()
 
void GeneratePrimaries (G4Event *)
 
void SetRndmFlag (G4String val)
 
void SetSourceType (G4int val)
 
void SetSpectrumType (G4int val)
 
void SetVertexRadius (G4double val)
 
void SetSourceGen (G4bool val)
 
- Public Member Functions inherited from G4VUserPrimaryGeneratorAction
 G4VUserPrimaryGeneratorAction ()
 
virtual ~G4VUserPrimaryGeneratorAction ()
 

Detailed Description

Definition at line 58 of file GammaRayTelPrimaryGeneratorAction.hh.

Constructor & Destructor Documentation

GammaRayTelPrimaryGeneratorAction::GammaRayTelPrimaryGeneratorAction ( )

Definition at line 58 of file GammaRayTelPrimaryGeneratorAction.cc.

59  :rndmFlag("off"),nSourceType(0),nSpectrumType(0),sourceGun(false)
60 {
61  GammaRayTelDetector = static_cast<const GammaRayTelDetectorConstruction*>
63 
64  //create a messenger for this class
65 
66  gunMessenger = new GammaRayTelPrimaryGeneratorMessenger(this);
67 
68  G4int n_particle = 1;
69 
70  particleGun = new G4ParticleGun(n_particle);
71  // default particle kinematic
72 
74  G4String particleName;
75  G4ParticleDefinition* particle
76  = particleTable->FindParticle(particleName="e-");
77  particleGun->SetParticleDefinition(particle);
78  particleGun->SetParticleMomentumDirection(G4ThreeVector(0.,0.,-1.));
79  particleGun->SetParticleEnergy(30.*MeV);
80  G4double position = 0.5*(GammaRayTelDetector->GetWorldSizeZ());
81  particleGun->SetParticlePosition(G4ThreeVector(0.*cm,0.*cm,position));
82  particleSource = new G4GeneralParticleSource();
83 
84 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
CLHEP::Hep3Vector G4ThreeVector
const G4VUserDetectorConstruction * GetUserDetectorConstruction() const
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
int G4int
Definition: G4Types.hh:78
void SetParticlePosition(G4ThreeVector aPosition)
static constexpr double cm
Definition: G4SIunits.hh:119
void SetParticleEnergy(G4double aKineticEnergy)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
static G4ParticleTable * GetParticleTable()
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)

Here is the call graph for this function:

GammaRayTelPrimaryGeneratorAction::~GammaRayTelPrimaryGeneratorAction ( )

Definition at line 88 of file GammaRayTelPrimaryGeneratorAction.cc.

89 {
90 
91  delete particleGun;
92  delete particleSource;
93 
94  delete gunMessenger;
95 }

Member Function Documentation

void GammaRayTelPrimaryGeneratorAction::GeneratePrimaries ( G4Event anEvent)
virtual

Implements G4VUserPrimaryGeneratorAction.

Definition at line 99 of file GammaRayTelPrimaryGeneratorAction.cc.

100 {
101  if (sourceGun)
102  {
103 
104  //this function is called at the begining of event
105  //
106  G4double z0 = 0.5*(GammaRayTelDetector->GetWorldSizeZ());
107  G4double x0 = 0.*cm, y0 = 0.*cm;
108 
109  G4ThreeVector pos0;
110  G4ThreeVector dir0;
111  G4ThreeVector vertex0 = G4ThreeVector(x0,y0,z0);
112 
113  dir0 = G4ThreeVector(0.,0.,-1.);
114 
115  G4double theta, phi, y, f;
116  G4double theta0=0.;
117  G4double phi0=0.;
118 
119  switch(nSourceType) {
120  case 0:
121  particleGun->SetParticlePosition(vertex0);
122  particleGun->SetParticleMomentumDirection(dir0);
123  break;
124  case 1:
125  // GS: Generate random position on the 4PIsphere to create a unif. distrib.
126  // GS: on the sphere
127  phi = G4UniformRand() * twopi;
128  do {
129  y = G4UniformRand()*1.0;
130  theta = G4UniformRand() * pi;
131  f = std::sin(theta);
132  } while (y > f);
133  vertex0 = G4ThreeVector(1.,0.,0.);
134  vertex0.setMag(dVertexRadius);
135  vertex0.setTheta(theta);
136  vertex0.setPhi(phi);
137  particleGun->SetParticlePosition(vertex0);
138 
139  dir0 = G4ThreeVector(1.,0.,0.);
140  do {
141  phi = G4UniformRand() * twopi;
142  do {
143  y = G4UniformRand()*1.0;
144  theta = G4UniformRand() * pi;
145  f = std::sin(theta);
146  } while (y > f);
147  dir0.setPhi(phi);
148  dir0.setTheta(theta);
149  } while (vertex0.dot(dir0) >= -0.7 * vertex0.mag());
151 
152  break;
153  case 2:
154  // GS: Generate random position on the upper semi-sphere z>0 to create a unif. distrib.
155  // GS: on a plane
156  phi = G4UniformRand() * twopi;
157  do {
158  y = G4UniformRand()*1.0;
159  theta = G4UniformRand() * halfpi;
160  f = std::sin(theta) * std::cos(theta);
161  } while (y > f);
162  vertex0 = G4ThreeVector(1.,0.,0.);
163 
164  G4double xy = GammaRayTelDetector->GetWorldSizeXY();
165  G4double z = GammaRayTelDetector->GetWorldSizeZ();
166 
167  if (dVertexRadius > xy*0.5)
168  {
169  G4cout << "vertexRadius too big " << G4endl;
170  G4cout << "vertexRadius setted to " << xy*0.45 << G4endl;
171  dVertexRadius = xy*0.45;
172  }
173 
174  if (dVertexRadius > z*0.5)
175  {
176  G4cout << "vertexRadius too high " << G4endl;
177  G4cout << "vertexRadius setted to " << z*0.45 << G4endl;
178  dVertexRadius = z*0.45;
179  }
180 
181 
182  vertex0.setMag(dVertexRadius);
183  vertex0.setTheta(theta);
184  vertex0.setPhi(phi);
185 
186  // GS: Get the user defined direction for the primaries and
187  // GS: Rotate the random position according to the user defined direction for the particle
188 
189  dir0 = particleGun->GetParticleMomentumDirection();
190  if (dir0.mag() > 0.001)
191  {
192  theta0 = dir0.theta();
193  phi0 = dir0.phi();
194  }
195 
196  if (theta0!=0.)
197  {
198  G4ThreeVector rotationAxis(1.,0.,0.);
199  rotationAxis.setPhi(phi0+halfpi);
200  vertex0.rotate(theta0+pi,rotationAxis);
201  }
202  particleGun->SetParticlePosition(vertex0);
203  break;
204  }
205 
206 
207  G4double pEnergy;
208 
209  switch(nSpectrumType) {
210  case 0:
211  break;
212  case 1:
213  break;
214  case 2:
215  do {
216  y = G4UniformRand()*100000.0;
217  pEnergy = G4UniformRand() * 10. * GeV;
218  f = std::pow(pEnergy * (1/GeV), -4.);
219  } while (y > f);
220 
221  particleGun->SetParticleEnergy(pEnergy);
222 
223  break;
224  case 3:
225  break;
226  }
227  G4cout << particleGun->GetParticleDefinition()->GetParticleName() << G4endl;
228  particleGun->GeneratePrimaryVertex(anEvent);
229  }
230  else
231  {
232  particleSource->GeneratePrimaryVertex(anEvent);
233  }
234 
235 }
void setPhi(double)
CLHEP::Hep3Vector G4ThreeVector
double dot(const Hep3Vector &) const
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
virtual void GeneratePrimaryVertex(G4Event *evt)
const G4String & GetParticleName() const
static constexpr double twopi
Definition: G4SIunits.hh:76
G4ParticleMomentum GetParticleMomentumDirection() const
void SetParticlePosition(G4ThreeVector aPosition)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
double phi() const
void setTheta(double)
double theta() const
void SetParticleEnergy(G4double aKineticEnergy)
G4ParticleDefinition * GetParticleDefinition() const
static constexpr double GeV
Definition: G4SIunits.hh:217
tuple z
Definition: test.py:28
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
static constexpr double halfpi
Definition: G4SIunits.hh:77
void setMag(double)
Definition: ThreeVector.cc:25
double G4double
Definition: G4Types.hh:76
double mag() const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28

Here is the call graph for this function:

void GammaRayTelPrimaryGeneratorAction::SetRndmFlag ( G4String  val)
inline

Definition at line 67 of file GammaRayTelPrimaryGeneratorAction.hh.

67 { rndmFlag = val;}

Here is the caller graph for this function:

void GammaRayTelPrimaryGeneratorAction::SetSourceGen ( G4bool  val)
inline

Definition at line 71 of file GammaRayTelPrimaryGeneratorAction.hh.

71 { sourceGun = val;}

Here is the caller graph for this function:

void GammaRayTelPrimaryGeneratorAction::SetSourceType ( G4int  val)
inline

Definition at line 68 of file GammaRayTelPrimaryGeneratorAction.hh.

68 { nSourceType = val;}

Here is the caller graph for this function:

void GammaRayTelPrimaryGeneratorAction::SetSpectrumType ( G4int  val)
inline

Definition at line 69 of file GammaRayTelPrimaryGeneratorAction.hh.

69 { nSpectrumType = val;}

Here is the caller graph for this function:

void GammaRayTelPrimaryGeneratorAction::SetVertexRadius ( G4double  val)
inline

Definition at line 70 of file GammaRayTelPrimaryGeneratorAction.hh.

70 { dVertexRadius = val;}

Here is the caller graph for this function:


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