Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
examples/extended/medical/GammaTherapy/src/PrimaryGeneratorAction.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: PrimaryGeneratorAction.cc 103662 2017-04-20 14:58:33Z gcosmo $
27 //
30 //
31 
32 //---------------------------------------------------------------------------
33 //
34 // ClassName: PrimaryGeneratorAction
35 //
36 // Description: Generate primary beam
37 //
38 // Authors: V.Grichine, V.Ivanchenko
39 //
40 // Modified:
41 //
42 //----------------------------------------------------------------------------
43 //
44 
45 #include "PrimaryGeneratorAction.hh"
46 #include "DetectorConstruction.hh"
47 #include "PrimaryGeneratorMessenger.hh"
48 #include "Randomize.hh"
49 #include "G4ParticleGun.hh"
50 #include "G4ParticleTable.hh"
51 #include "G4ParticleDefinition.hh"
52 #include "G4PhysicalConstants.hh"
53 #include "G4SystemOfUnits.hh"
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57 
59  fDetector(pDet)
60 {
61  InitializeMe();
62 }
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65 
67 {
71  fCounter = 0;
72  fX0 = 0.0;
73  fY0 = 0.0;
74  fZ0 = 0.0;
75  fSigmaX = 1.5*mm;
76  fSigmaY = 1.5*mm;
77  fSigmaZ = 0.0;
78  fSigmaE = 0.0;
79  fRMax2 = 2.5*2.5*mm*mm;
80  fSigmaTheta = 0.0;
81  // fSigmaTheta = 0.17*degree;
82  fMinCosTheta = 2.0;
83  SetBeamEnergy(50.0*MeV);
85  fDirection = G4ThreeVector(0.0,0.0,1.0);
86  fGauss = true;
87 }
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92 {
93  delete fParticleGun;
94  delete fMessenger;
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98 
100 {
101  fCounter++ ;
102 
103  // Simulation of beam position
104  G4double x = fX0;
105  G4double y = fY0;
107  do {
108  if(0.0 < fSigmaX) { x = G4RandGauss::shoot(fX0,fSigmaX); }
109  if(0.0 < fSigmaY) { y = G4RandGauss::shoot(fY0,fSigmaY); }
110  } while (x*x + y*y > fRMax2);
111 
112  fPosition = G4ThreeVector(x,y,z);
114 
115  // Simulation of beam direction
116  G4double ux = fDirection.x();
117  G4double uy = fDirection.y();
118  G4double uz = fDirection.z();
119 
120  // Beam particles are uniformly distributed over phi, cosTheta
121  if(1.0 > fMinCosTheta) {
122  uz = fMinCosTheta + (1.0 - fMinCosTheta)*G4UniformRand() ;
123  ux = std::sqrt((1.0 - uz)*(1.0 + uz)) ;
124  } else if (fSigmaTheta > 0.0) {
126  uz = std::sqrt((1.0 - ux)*(1.0 + ux));
127  }
128 
129  G4double phi = twopi*G4UniformRand() ;
130  uy = ux;
131  ux *= std::cos(phi) ;
132  uy *= std::sin(phi) ;
133  fDirection = G4ThreeVector(ux,uy,uz) ;
134 
136 
137  // Simulation of beam kinetic energy
138  G4double kinEnergy = fEnergy;
139 
140  if(fGauss == "flatE") {
141  kinEnergy = fEnergy - fSigmaE + 2.*fSigmaE*G4UniformRand();
142  } else if(0.0 < fSigmaE) {
143  kinEnergy = fEnergy + G4RandGauss::shoot(0.0,fSigmaE);
144  }
145  fParticleGun->SetParticleEnergy(kinEnergy);
146 
147  if(fVerbose > 0) {
149  G4String particleName = particle->GetParticleName() ;
150  G4cout << "Event# " << fCounter
151  << " Beam particle is generated by PrimaryGeneratorAction "
152  << G4endl;
153  G4cout << "ParticleName= " << particleName
154  << " PDGcode= " << particle->GetPDGEncoding()
155  << std::setprecision(5)
156  << " KinEnergy(GeV)= "
157  << fEnergy/GeV
158  << " x(mm)= "
159  << x/mm
160  << " y(mm)= "
161  << y/mm
162  << " z(mm)= "
163  << z/mm
164  << " ux= "
165  << ux
166  << " uy= "
167  << uy
168  << " uz= "
169  << uz
170  << G4endl;
171  }
172 
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177 
179 {
180  fEnergy = val;
181  if(fEnergy<fDetector->GetMaxEnergy()) fDetector->SetMaxEnergy(fEnergy);
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185 
186 
187 
Float_t x
Definition: compare.C:6
CLHEP::Hep3Vector G4ThreeVector
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double mm
Definition: G4SIunits.hh:115
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
Double_t z
const G4String & GetParticleName() const
double z() const
void SetParticlePosition(G4ThreeVector aPosition)
double G4double
Definition: G4Types.hh:76
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual void GeneratePrimaryVertex(G4Event *evt)
G4GLOB_DLL std::ostream G4cout
double x() const
void SetParticleEnergy(G4double aKineticEnergy)
double y() const
G4ParticleDefinition * GetParticleDefinition() const
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetParticleMomentumDirection(G4ParticleMomentum aMomentumDirection)
Simple detector construction with a box volume placed in a world.