Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
B5PrimaryGeneratorAction.cc
이 파일의 문서화 페이지로 가기
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 101036 2016-11-04 09:00:23Z 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(nullptr), fMessenger(nullptr),
46  fPositron(nullptr), fMuon(nullptr), fPion(nullptr),
47  fKaon(nullptr), fProton(nullptr),
48  fMomentum(1000.*MeV),
49  fSigmaMomentum(50.*MeV),
50  fSigmaAngle(2.*deg),
51  fRandomizePrimary(true)
52 {
53  G4int nofParticles = 1;
54  fParticleGun = new G4ParticleGun(nofParticles);
55 
56  auto particleTable = G4ParticleTable::GetParticleTable();
57  fPositron = particleTable->FindParticle("e+");
58  fMuon = particleTable->FindParticle("mu+");
59  fPion = particleTable->FindParticle("pi+");
60  fKaon = particleTable->FindParticle("kaon+");
61  fProton = particleTable->FindParticle("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  if (fRandomizePrimary) {
85  G4int i = (int)(5.*G4UniformRand());
86  switch(i) {
87  case 0:
88  particle = fPositron;
89  break;
90  case 1:
91  particle = fMuon;
92  break;
93  case 2:
94  particle = fPion;
95  break;
96  case 3:
97  particle = fKaon;
98  break;
99  default:
100  particle = fProton;
101  break;
102  }
104  }
105  else {
106  particle = fParticleGun->GetParticleDefinition();
107  }
108 
109  auto pp = fMomentum + (G4UniformRand()-0.5)*fSigmaMomentum;
110  auto mass = particle->GetPDGMass();
111  auto ekin = std::sqrt(pp*pp+mass*mass)-mass;
113 
114  auto angle = (G4UniformRand()-0.5)*fSigmaAngle;
116  G4ThreeVector(std::sin(angle),0.,std::cos(angle)));
117 
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
124 {
125  // Define /B5/generator command directory using generic messenger class
126  fMessenger
127  = new G4GenericMessenger(this,
128  "/B5/generator/",
129  "Primary generator control");
130 
131  // momentum command
132  auto& momentumCmd
133  = fMessenger->DeclarePropertyWithUnit("momentum", "GeV", fMomentum,
134  "Mean momentum of primaries.");
135  momentumCmd.SetParameterName("p", true);
136  momentumCmd.SetRange("p>=0.");
137  momentumCmd.SetDefaultValue("1.");
138  // ok
139  //momentumCmd.SetParameterName("p", true);
140  //momentumCmd.SetRange("p>=0.");
141 
142  // sigmaMomentum command
143  auto& sigmaMomentumCmd
144  = fMessenger->DeclarePropertyWithUnit("sigmaMomentum",
145  "MeV", fSigmaMomentum, "Sigma momentum of primaries.");
146  sigmaMomentumCmd.SetParameterName("sp", true);
147  sigmaMomentumCmd.SetRange("sp>=0.");
148  sigmaMomentumCmd.SetDefaultValue("50.");
149 
150  // sigmaAngle command
151  auto& sigmaAngleCmd
152  = fMessenger->DeclarePropertyWithUnit("sigmaAngle", "deg", fSigmaAngle,
153  "Sigma angle divergence of primaries.");
154  sigmaAngleCmd.SetParameterName("t", true);
155  sigmaAngleCmd.SetRange("t>=0.");
156  sigmaAngleCmd.SetDefaultValue("2.");
157 
158  // randomizePrimary command
159  auto& randomCmd
160  = fMessenger->DeclareProperty("randomizePrimary", fRandomizePrimary);
161  G4String guidance
162  = "Boolean flag for randomizing primary particle types.\n";
163  guidance
164  += "In case this flag is false, you can select the primary particle\n";
165  guidance += " with /gun/particle command.";
166  randomCmd.SetGuidance(guidance);
167  randomCmd.SetParameterName("flg", true);
168  randomCmd.SetDefaultValue("true");
169 }
170 
171 //..oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void GeneratePrimaries(G4Event *)
This class is generic messenger.
CLHEP::Hep3Vector G4ThreeVector
Command & DeclarePropertyWithUnit(const G4String &name, const G4String &defaultUnit, const G4AnyType &variable, const G4String &doc="")
static G4ParticleTable * GetParticleTable()
static constexpr double MeV
Definition: G4SIunits.hh:214
Definition of the B5PrimaryGeneratorAction class.
Command & SetDefaultValue(const G4String &)
void SetParticlePosition(G4ThreeVector aPosition)
G4double GetPDGMass() const
static constexpr double m
Definition: G4SIunits.hh:129
void SetParticleDefinition(G4ParticleDefinition *aParticleDefinition)
static constexpr double deg
Definition: G4SIunits.hh:152
#define G4UniformRand()
Definition: Randomize.hh:53
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
Command & SetRange(const G4String &range)
G4ParticleDefinition * fProton
Command & SetParameterName(const G4String &, G4bool, G4bool=false)
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * fPositron
virtual void GeneratePrimaryVertex(G4Event *evt)
void SetParticleEnergy(G4double aKineticEnergy)
G4ParticleDefinition * GetParticleDefinition() const
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
Command & DeclareProperty(const G4String &name, const G4AnyType &variable, const G4String &doc="")
Declare Methods.