File indexing completed on 2025-02-23 09:22:22
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
0034 #include "OpNovicePrimaryGeneratorAction.hh"
0035
0036 #include "OpNovicePrimaryGeneratorMessenger.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 OpNovicePrimaryGeneratorAction::OpNovicePrimaryGeneratorAction()
0047 : G4VUserPrimaryGeneratorAction(), fParticleGun(nullptr)
0048 {
0049 G4int n_particle = 1;
0050 fParticleGun = new G4ParticleGun(n_particle);
0051
0052 fGunMessenger = new OpNovicePrimaryGeneratorMessenger(this);
0053
0054
0055 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
0056 G4ParticleDefinition* particle = particleTable->FindParticle("e+");
0057
0058 fParticleGun->SetParticleDefinition(particle);
0059 fParticleGun->SetParticleTime(0.0 * ns);
0060 fParticleGun->SetParticlePosition(G4ThreeVector(0., 0., 0.));
0061 fParticleGun->SetParticleMomentumDirection(G4ThreeVector(1., 0., 0.));
0062 fParticleGun->SetParticleEnergy(500.0 * keV);
0063 }
0064
0065
0066 OpNovicePrimaryGeneratorAction::~OpNovicePrimaryGeneratorAction()
0067 {
0068 delete fParticleGun;
0069 delete fGunMessenger;
0070 }
0071
0072
0073 void OpNovicePrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
0074 {
0075 fParticleGun->GeneratePrimaryVertex(anEvent);
0076 }
0077
0078
0079 void OpNovicePrimaryGeneratorAction::SetOptPhotonPolar()
0080 {
0081 G4double angle = G4UniformRand() * 360.0 * deg;
0082 SetOptPhotonPolar(angle);
0083 }
0084
0085
0086 void OpNovicePrimaryGeneratorAction::SetOptPhotonPolar(G4double angle)
0087 {
0088 if (fParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") {
0089 G4ExceptionDescription ed;
0090 ed << "Warning: the particleGun is not an opticalphoton";
0091 G4Exception("OpNovicePrimaryGeneratorAction::SetOptPhotonPolar()", "OpNovice_010", JustWarning,
0092 ed);
0093 return;
0094 }
0095
0096 G4ThreeVector normal(1., 0., 0.);
0097 G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection();
0098 G4ThreeVector product = normal.cross(kphoton);
0099 G4double modul2 = product * product;
0100
0101 G4ThreeVector e_perpend(0., 0., 1.);
0102 if (modul2 > 0.) e_perpend = (1. / std::sqrt(modul2)) * product;
0103 G4ThreeVector e_paralle = e_perpend.cross(kphoton);
0104
0105 G4ThreeVector polar = std::cos(angle) * e_paralle + std::sin(angle) * e_perpend;
0106 fParticleGun->SetParticlePolarization(polar);
0107 }
0108