Geant4  10.03
G4ParticleGun.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 //
26 //
27 // $Id: G4ParticleGun.hh 94950 2016-01-07 11:53:14Z gcosmo $
28 //
29 
30 #ifndef G4ParticleGun_h
31 #define G4ParticleGun_h 1
32 
33 
34 #include "globals.hh"
35 #include "G4VPrimaryGenerator.hh"
36 #include "G4ThreeVector.hh"
37 #include "G4ParticleDefinition.hh"
38 #include "G4PrimaryVertex.hh"
39 #include "G4ParticleMomentum.hh"
40 
41 class G4Event;
43 
44 // class description:
45 //
46 // This is a concrete class of G4VPrimaryGenerator. It shoots a particle of given type
47 // into a given direction with either a given kinetic energy or momentum.
48 // The position and time of the primary particle must be set by the corresponding
49 // set methods of G4VPrimaryGenerator base class, otherwise zero will be set.
50 //
51 // The FAQ to this class is for randomizing position/direction/kinetic energy of primary
52 // particle. But, G4ParticleGun does NOT have any way of randomization. Instead, the user's
53 // concrete implementation of G4VUserPrimaryGeneratorAction which transmits G4Event object
54 // to this particle gun can randomize these quantities and set to this particle gun before
55 // invoking GeneratePrimaryVertex() method.
56 // Note that, even if the particle gun shoots more than one particles at one invokation of
57 // GeneratePrimaryVertex() method, all particles have the same physical quantities. If the
58 // user wants to shoot two particles with different momentum, position, etc., invoke
59 // GeneratePrimaryVertex() method twice and set quantities on demand to the particle gun.
60 //
61 
63 {
64  public: // with description
65  G4ParticleGun();
66  G4ParticleGun(G4int numberofparticles);
67  G4ParticleGun(G4ParticleDefinition * particleDef,
68  G4int numberofparticles = 1);
69  // costructors. "numberofparticles" is number of particles to be shoot at one invokation
70  // of GeneratePrimaryVertex() method. All paricles are shot with the same physical
71  // quantities.
72 
73  public:
74  virtual ~G4ParticleGun();
75 
76  private:
77  G4ParticleGun(const G4ParticleGun&) = delete;
78  const G4ParticleGun & operator=(const G4ParticleGun&) = delete;
79  G4int operator==(const G4ParticleGun&) const = delete;
80  G4int operator!=(const G4ParticleGun&) const = delete;
81 
82  public: // with description
83  virtual void GeneratePrimaryVertex(G4Event* evt);
84  // Creates a primary vertex at the given point and put primary particles to it.
85  // Followings are set methods for the particle properties.
86  // SetParticleDefinition should be called first.
87  // By using SetParticleMomentum(), both particle_momentum_direction and
88  // particle_energy(Kinetic Energy) are set.
89  //
91  (G4ParticleDefinition * aParticleDefinition);
92  void SetParticleEnergy(G4double aKineticEnergy);
93  void SetParticleMomentum(G4double aMomentum);
96  (G4ParticleMomentum aMomentumDirection)
97  { particle_momentum_direction = aMomentumDirection.unit(); }
99  { particle_charge = aCharge; }
101  { particle_polarization = aVal; }
104 
105  public:
107  { return particle_definition; }
109  { return particle_momentum_direction; }
111  { return particle_energy; }
113  { return particle_momentum; }
115  { return particle_charge; }
117  { return particle_polarization; }
120 
121  protected:
122  virtual void SetInitialValues();
123 
131 
132  private:
134 };
135 
136 #endif
137 
138 
139 
140 
141 
142 
143 
void SetParticleMomentum(G4double aMomentum)
CLHEP::Hep3Vector G4ThreeVector
G4double GetParticleMomentum() const
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
G4double particle_momentum
int G4int
Definition: G4Types.hh:78
virtual void GeneratePrimaryVertex(G4Event *evt)
G4ParticleDefinition * particle_definition
G4ParticleMomentum GetParticleMomentumDirection() const
G4int operator==(const G4ParticleGun &) const =delete
void SetParticlePolarization(G4ThreeVector aVal)
G4ThreeVector GetParticlePolarization() const
G4ParticleMomentum particle_momentum_direction
void SetParticleCharge(G4double aCharge)
virtual void SetInitialValues()
G4double GetParticleCharge() const
const G4ParticleGun & operator=(const G4ParticleGun &)=delete
void SetNumberOfParticles(G4int i)
G4double particle_charge
G4ThreeVector particle_polarization
G4int GetNumberOfParticles() const
void SetParticleEnergy(G4double aKineticEnergy)
G4int NumberOfParticlesToBeGenerated
G4double particle_energy
G4ParticleGunMessenger * theMessenger
G4int operator!=(const G4ParticleGun &) const =delete
G4ParticleDefinition * GetParticleDefinition() const
virtual ~G4ParticleGun()
double G4double
Definition: G4Types.hh:76
G4ThreeVector G4ParticleMomentum
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GetParticleEnergy() const