Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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// G4ParticleGum
27//
28// Class description:
29//
30// This is a concrete class of G4VPrimaryGenerator. It shoots a particle of
31// given type into a given direction with either a given kinetic energy or
32// momentum.
33// The position and time of the primary particle must be set by the
34// corresponding set methods of the G4VPrimaryGenerator base class, otherwise
35// zero will be set.
36//
37// The FAQ to this class is for randomizing position/direction/kinetic energy
38// of the primary particle. But, G4ParticleGun does NOT have any way of
39// randomization. Instead, the user's concrete implementation of
40// G4VUserPrimaryGeneratorAction which transfer the G4Event object
41// to this particle gun can randomize these quantities and set to this
42// particle gun before invoking GeneratePrimaryVertex() method.
43// Note that, even if the particle gun shoots more than one particles at one
44// invokation of the GeneratePrimaryVertex() method, all particles have the
45// same physical quantities. If the user wants to shoot two particles with
46// different momentum, position, etc., should invoke GeneratePrimaryVertex()
47// method twice and set quantities on demand to the particle gun.
48
49// Author: Makoto Asai, 1997
50// --------------------------------------------------------------------
51#ifndef G4ParticleGun_hh
52#define G4ParticleGun_hh 1
53
54#include "globals.hh"
56#include "G4ThreeVector.hh"
58#include "G4PrimaryVertex.hh"
59#include "G4ParticleMomentum.hh"
60
61class G4Event;
63
65{
66 public:
67
69 explicit G4ParticleGun(G4int numberofparticles);
70 explicit G4ParticleGun(G4ParticleDefinition* particleDef,
71 G4int numberofparticles = 1);
72 // Costructors. "numberofparticles" is the number of particles to be
73 // shot at one invokation of GeneratePrimaryVertex() method.
74 // All particles are shot with the same physical quantities.
75
76 ~G4ParticleGun() override;
77
78 G4ParticleGun(const G4ParticleGun&) = delete;
79 const G4ParticleGun& operator=(const G4ParticleGun&) = delete;
80 G4bool operator==(const G4ParticleGun&) const = delete;
81 G4bool operator!=(const G4ParticleGun&) const = delete;
82
83 void GeneratePrimaryVertex(G4Event* evt) override;
84 // Creates a primary vertex at the given point
85 // and put primary particles to it.
86
87 // Followings are the Set methods for the particle properties.
88 // SetParticleDefinition() should be called first.
89 // By using SetParticleMomentum(), both particle_momentum_direction and
90 // particle_energy(Kinetic Energy) are set.
91 //
92 void SetParticleDefinition(G4ParticleDefinition* aParticleDefinition);
93 void SetParticleEnergy(G4double aKineticEnergy);
94 void SetParticleMomentum(G4double aMomentum);
97 { particle_momentum_direction = aMomDirection.unit(); }
98 inline void SetParticleCharge(G4double aCharge)
99 { particle_charge = aCharge; }
101 { particle_polarization = aVal; }
104
106 { return particle_definition; }
110 { return particle_energy; }
112 { return particle_momentum; }
114 { return particle_charge; }
116 { return particle_polarization; }
119
120 protected:
121
122 virtual void SetInitialValues();
123
131
132 private:
133
134 G4ParticleGunMessenger* theMessenger = nullptr;
135};
136
137#endif
138
139
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Hep3Vector unit() const
G4double particle_momentum
virtual void SetInitialValues()
const G4ParticleGun & operator=(const G4ParticleGun &)=delete
G4ThreeVector GetParticlePolarization() const
G4ParticleMomentum GetParticleMomentumDirection() const
G4double particle_energy
G4bool operator!=(const G4ParticleGun &) const =delete
G4double GetParticleMomentum() const
G4ParticleGun(const G4ParticleGun &)=delete
G4bool operator==(const G4ParticleGun &) const =delete
void SetNumberOfParticles(G4int i)
void SetParticlePolarization(G4ThreeVector aVal)
G4ParticleDefinition * GetParticleDefinition() const
void GeneratePrimaryVertex(G4Event *evt) override
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
G4ThreeVector particle_polarization
G4ParticleMomentum particle_momentum_direction
G4int NumberOfParticlesToBeGenerated
G4double particle_charge
~G4ParticleGun() override
G4int GetNumberOfParticles() const
void SetParticleEnergy(G4double aKineticEnergy)
void SetParticleMomentumDirection(G4ParticleMomentum aMomDirection)
G4double GetParticleCharge() const
G4ParticleDefinition * particle_definition
void SetParticleMomentum(G4double aMomentum)
G4double GetParticleEnergy() const
void SetParticleCharge(G4double aCharge)