Geant4  10.01
B5PrimaryGeneratorAction.cc
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 // $Id: B5PrimaryGeneratorAction.cc 77781 2013-11-28 07:54:07Z gcosmo $
27 //
30 
32 
33 #include "G4Event.hh"
34 #include "G4ParticleGun.hh"
35 #include "G4ParticleTable.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4GenericMessenger.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "Randomize.hh"
40 
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42 
45  fParticleGun(0), fMessenger(0),
46  fPositron(0), fMuon(0), fPion(0), fKaon(0), fProton(0),
47  fMomentum(1000.*MeV),
48  fSigmaMomentum(50.*MeV),
49  fSigmaAngle(2.*deg),
50  fRandomizePrimary(true)
51 {
52  G4int n_particle = 1;
53  fParticleGun = new G4ParticleGun(n_particle);
54 
56  G4String particleName;
57  fPositron = particleTable->FindParticle(particleName="e+");
58  fMuon = particleTable->FindParticle(particleName="mu+");
59  fPion = particleTable->FindParticle(particleName="pi+");
60  fKaon = particleTable->FindParticle(particleName="kaon+");
61  fProton = particleTable->FindParticle(particleName="proton");
62 
63  // default particle kinematics
66 
67  // define commands for this class
69 }
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 {
75  delete fParticleGun;
76  delete fMessenger;
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80 
82 {
83  G4ParticleDefinition* particle;
84 
86  {
87  G4int i = (int)(5.*G4UniformRand());
88  switch(i)
89  {
90  case 0:
91  particle = fPositron;
92  break;
93  case 1:
94  particle = fMuon;
95  break;
96  case 2:
97  particle = fPion;
98  break;
99  case 3:
100  particle = fKaon;
101  break;
102  default:
103  particle = fProton;
104  break;
105  }
107  }
108  else
109  {
110  particle = fParticleGun->GetParticleDefinition();
111  }
112 
114  G4double mass = particle->GetPDGMass();
115  G4double Ekin = std::sqrt(pp*pp+mass*mass)-mass;
117 
118  G4double angle = (G4UniformRand()-0.5)*fSigmaAngle;
120  std::cos(angle)));
121 
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126 
128 {
129  // Define /B5/generator command directory using generic messenger class
130  fMessenger
131  = new G4GenericMessenger(this,
132  "/B5/generator/",
133  "Primary generator control");
134 
135  // momentum command
136  G4GenericMessenger::Command& momentumCmd
137  = fMessenger->DeclarePropertyWithUnit("momentum", "GeV", fMomentum,
138  "Mean momentum of primaries.");
139  momentumCmd.SetParameterName("p", true);
140  momentumCmd.SetRange("p>=0.");
141  momentumCmd.SetDefaultValue("1.");
142  // ok
143  //momentumCmd.SetParameterName("p", true);
144  //momentumCmd.SetRange("p>=0.");
145 
146  // sigmaMomentum command
147  G4GenericMessenger::Command& sigmaMomentumCmd
148  = fMessenger->DeclarePropertyWithUnit("sigmaMomentum",
149  "MeV", fSigmaMomentum, "Sigma momentum of primaries.");
150  sigmaMomentumCmd.SetParameterName("sp", true);
151  sigmaMomentumCmd.SetRange("sp>=0.");
152  sigmaMomentumCmd.SetDefaultValue("50.");
153 
154  // sigmaAngle command
155  G4GenericMessenger::Command& sigmaAngleCmd
156  = fMessenger->DeclarePropertyWithUnit("sigmaAngle", "deg", fSigmaAngle,
157  "Sigma angle divergence of primaries.");
158  sigmaAngleCmd.SetParameterName("t", true);
159  sigmaAngleCmd.SetRange("t>=0.");
160  sigmaAngleCmd.SetDefaultValue("2.");
161 
162  // randomizePrimary command
163  G4GenericMessenger::Command& randomCmd
164  = fMessenger->DeclareProperty("randomizePrimary", fRandomizePrimary);
165  G4String guidance
166  = "Boolean flag for randomizing primary particle types.\n";
167  guidance
168  += "In case this flag is false, you can select the primary particle\n";
169  guidance += " with /gun/particle command.";
170  randomCmd.SetGuidance(guidance);
171  randomCmd.SetParameterName("flg", true);
172  randomCmd.SetDefaultValue("true");
173 }
174 
175 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static const double MeV
Definition: G4SIunits.hh:193
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
CLHEP::Hep3Vector G4ThreeVector
This class is generic messenger.
Command & DeclareProperty(const G4String &name, const G4AnyType &variable, const G4String &doc="")
Declare Methods.
virtual void GeneratePrimaries(G4Event *)
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
int G4int
Definition: G4Types.hh:78
virtual void GeneratePrimaryVertex(G4Event *evt)
Command & SetDefaultValue(const G4String &)
void SetParticlePosition(G4ThreeVector aPosition)
#define G4UniformRand()
Definition: Randomize.hh:95
Command & SetGuidance(const G4String &s0)
G4ParticleDefinition * fPositron
static const double deg
Definition: G4SIunits.hh:133
Command & SetRange(const G4String &range)
void SetParticleEnergy(G4double aKineticEnergy)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
Definition of the B5PrimaryGeneratorAction class.
G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * fProton
Command & SetParameterName(const G4String &, G4bool, G4bool=false)
static const double m
Definition: G4SIunits.hh:110
double G4double
Definition: G4Types.hh:76
Command & DeclarePropertyWithUnit(const G4String &name, const G4String &defaultUnit, const G4AnyType &variable, const G4String &doc="")
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)