File indexing completed on 2025-02-23 09:22:40
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033 #include "PrimaryGeneratorAction.hh"
0034
0035 #include "DetectorConstruction.hh"
0036 #include "PrimaryGeneratorMessenger.hh"
0037
0038 #include "G4Event.hh"
0039 #include "G4ParticleDefinition.hh"
0040 #include "G4ParticleGun.hh"
0041 #include "G4ParticleTable.hh"
0042 #include "G4SystemOfUnits.hh"
0043 #include "Randomize.hh"
0044
0045
0046
0047 PrimaryGeneratorAction::PrimaryGeneratorAction(DetectorConstruction* det) : fDetector(det)
0048 {
0049 fParticleGun = new G4ParticleGun(1);
0050 SetDefaultKinematic();
0051
0052
0053 fGunMessenger = new PrimaryGeneratorMessenger(this);
0054 }
0055
0056
0057
0058 PrimaryGeneratorAction::~PrimaryGeneratorAction()
0059 {
0060 delete fParticleGun;
0061 delete fGunMessenger;
0062 }
0063
0064
0065
0066 void PrimaryGeneratorAction::SetDefaultKinematic()
0067 {
0068 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle("proton");
0069 fParticleGun->SetParticleDefinition(particle);
0070 fParticleGun->SetParticleEnergy(10 * MeV);
0071 fParticleGun->SetParticleMomentumDirection(G4ThreeVector(1., 0., 0.));
0072 G4double position = -0.45 * (fDetector->GetWorldSizeX());
0073 fParticleGun->SetParticlePosition(G4ThreeVector(position, 0. * cm, 0. * cm));
0074 }
0075
0076
0077 void PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
0078 {
0079
0080
0081
0082
0083 if (fRndmBeam > 0.) {
0084 G4ThreeVector oldPosition = fParticleGun->GetParticlePosition();
0085 if (fRndmBeam > fDetector->GetAbsorSizeYZ()) fRndmBeam = fDetector->GetAbsorSizeYZ();
0086 G4double rbeam = 0.5 * fRndmBeam;
0087 G4double x0 = oldPosition.x();
0088 G4double y0 = oldPosition.y() + (2 * G4UniformRand() - 1.) * rbeam;
0089 G4double z0 = oldPosition.z() + (2 * G4UniformRand() - 1.) * rbeam;
0090 fParticleGun->SetParticlePosition(G4ThreeVector(x0, y0, z0));
0091 fParticleGun->GeneratePrimaryVertex(anEvent);
0092 fParticleGun->SetParticlePosition(oldPosition);
0093 }
0094
0095 else
0096 fParticleGun->GeneratePrimaryVertex(anEvent);
0097
0098
0099
0100 G4double t0 = fTimeExposure * G4UniformRand();
0101 anEvent->GetPrimaryVertex()->SetT0(t0);
0102 }
0103
0104