File indexing completed on 2026-06-01 07:55:10
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 #include "PrimaryGeneratorAction.hh"
0030
0031 #include "PrimaryGeneratorMessenger.hh"
0032
0033 #include "G4Event.hh"
0034 #include "G4OpticalPhoton.hh"
0035 #include "G4ParticleDefinition.hh"
0036 #include "G4ParticleGun.hh"
0037 #include "G4ParticleTable.hh"
0038 #include "G4SystemOfUnits.hh"
0039 #include "Randomize.hh"
0040
0041
0042
0043 PrimaryGeneratorAction::PrimaryGeneratorAction()
0044 : G4VUserPrimaryGeneratorAction(), fParticleGun(nullptr)
0045 {
0046 G4int n_particle = 1;
0047 fParticleGun = new G4ParticleGun(n_particle);
0048
0049
0050 fGunMessenger = new PrimaryGeneratorMessenger(this);
0051
0052 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
0053 G4ParticleDefinition* particle = particleTable->FindParticle("e+");
0054
0055 fParticleGun->SetParticleDefinition(particle);
0056 fParticleGun->SetParticleTime(0.0 * ns);
0057 fParticleGun->SetParticlePosition(G4ThreeVector(0.0 * cm, 0.0 * cm, 0.0 * cm));
0058 fParticleGun->SetParticleMomentumDirection(G4ThreeVector(1., 0., 0.));
0059 fParticleGun->SetParticleEnergy(500.0 * keV);
0060 }
0061
0062
0063
0064 PrimaryGeneratorAction::~PrimaryGeneratorAction()
0065 {
0066 delete fParticleGun;
0067 delete fGunMessenger;
0068 }
0069
0070
0071
0072 void PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
0073 {
0074 if (fRandomDirection) {
0075 G4double theta = CLHEP::halfpi * G4UniformRand();
0076 G4double phi = CLHEP::twopi * G4UniformRand();
0077 G4double x = std::cos(theta);
0078 G4double y = std::sin(theta) * std::sin(phi);
0079 G4double z = std::sin(theta) * std::cos(phi);
0080 G4ThreeVector dir(x, y, z);
0081 fParticleGun->SetParticleMomentumDirection(dir);
0082 }
0083 if (fParticleGun->GetParticleDefinition() == G4OpticalPhoton::OpticalPhotonDefinition()) {
0084 if (fPolarized)
0085 SetOptPhotonPolar(fPolarization);
0086 else
0087 SetOptPhotonPolar();
0088 }
0089 fParticleGun->GeneratePrimaryVertex(anEvent);
0090 }
0091
0092
0093
0094 void PrimaryGeneratorAction::SetOptPhotonPolar()
0095 {
0096 G4double angle = G4UniformRand() * 360.0 * deg;
0097 SetOptPhotonPolar(angle);
0098 }
0099
0100
0101
0102 void PrimaryGeneratorAction::SetOptPhotonPolar(G4double angle)
0103 {
0104 if (fParticleGun->GetParticleDefinition() != G4OpticalPhoton::OpticalPhotonDefinition()) {
0105 G4ExceptionDescription ed;
0106 ed << "The particleGun is not an opticalphoton.";
0107 G4Exception("PrimaryGeneratorAction::SetOptPhotonPolar", "OpNovice2_004", JustWarning, ed);
0108 return;
0109 }
0110
0111 fPolarized = true;
0112 fPolarization = angle;
0113
0114 G4ThreeVector normal(1., 0., 0.);
0115 G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection();
0116 G4ThreeVector product = normal.cross(kphoton);
0117 G4double modul2 = product * product;
0118
0119 G4ThreeVector e_perpend(0., 0., 1.);
0120 if (modul2 > 0.) e_perpend = (1. / std::sqrt(modul2)) * product;
0121 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
0122
0123 G4ThreeVector polar = std::cos(angle) * e_paralle + std::sin(angle) * e_perpend;
0124 fParticleGun->SetParticlePolarization(polar);
0125 }
0126
0127
0128 void PrimaryGeneratorAction::SetRandomDirection(G4bool val)
0129 {
0130 fRandomDirection = val;
0131 }