File indexing completed on 2025-02-23 09:22:04
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
0035
0036
0037
0038
0039
0040
0041
0042 #include "PrimaryGeneratorAction.hh"
0043
0044 #include "G4SystemOfUnits.hh"
0045
0046 #include "CommandLineParser.hh"
0047
0048 #include "G4Box.hh"
0049 #include "G4Event.hh"
0050 #include "G4LogicalVolume.hh"
0051 #include "G4LogicalVolumeStore.hh"
0052 #include "G4Orb.hh"
0053 #include "G4ParticleDefinition.hh"
0054 #include "G4ParticleGun.hh"
0055 #include "G4PhysicalConstants.hh"
0056 #include "G4PhysicalVolumeStore.hh"
0057 #include "G4Proton.hh"
0058 #include "G4VPhysicalVolume.hh"
0059 #include "Randomize.hh"
0060
0061 using namespace G4DNAPARSER;
0062
0063
0064
0065 PrimaryGeneratorAction::PrimaryGeneratorAction() : G4VUserPrimaryGeneratorAction()
0066 {
0067 G4int n_particle = 1;
0068 fpParticleGun = new G4ParticleGun(n_particle);
0069
0070 G4ParticleDefinition* particle = G4Proton::Proton();
0071 fpParticleGun->SetParticleDefinition(particle);
0072
0073 fpParticleGun->SetParticleEnergy(10. * MeV);
0074 fpParticleGun->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.));
0075 fpParticleGun->SetParticlePosition(G4ThreeVector(0., 0., 0.));
0076 }
0077
0078
0079
0080 PrimaryGeneratorAction::~PrimaryGeneratorAction()
0081 {
0082 delete fpParticleGun;
0083 }
0084
0085
0086
0087 void PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
0088 {
0089
0090
0091
0092
0093
0094
0095
0096
0097
0098
0099
0100
0101 G4double mediumRadius = 0.;
0102 G4double boundingXHalfLength = 0.;
0103 G4double boundingYHalfLength = 0.;
0104 G4double boundingZHalfLength = 0.;
0105 G4LogicalVolume* mediumLV = G4LogicalVolumeStore::GetInstance()->GetVolume("Medium");
0106 G4LogicalVolume* boundingLV = G4LogicalVolumeStore::GetInstance()->GetVolume("BoundingSlice");
0107 G4Orb* mediumSphere = 0;
0108 G4Box* boundingSlice = 0;
0109 if (mediumLV && boundingLV) {
0110 mediumSphere = dynamic_cast<G4Orb*>(mediumLV->GetSolid());
0111 boundingSlice = dynamic_cast<G4Box*>(boundingLV->GetSolid());
0112 }
0113 if (mediumSphere && boundingSlice) {
0114 mediumRadius = mediumSphere->GetRadius();
0115 boundingXHalfLength = boundingSlice->GetXHalfLength();
0116 boundingYHalfLength = boundingSlice->GetYHalfLength();
0117 boundingZHalfLength = boundingSlice->GetZHalfLength();
0118
0119
0120
0121
0122 if (CommandLineParser::GetParser()->GetCommandIfActive("-sXY")) {
0123
0124 fpParticleGun->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.));
0125
0126 fpParticleGun->SetParticlePosition(G4ThreeVector(
0127 CLHEP::RandFlat::shoot(-boundingXHalfLength, boundingXHalfLength),
0128 CLHEP::RandFlat::shoot(-boundingYHalfLength, boundingYHalfLength), -mediumRadius));
0129 }
0130
0131
0132
0133
0134 else if (CommandLineParser::GetParser()->GetCommandIfActive("-dXY")) {
0135 fpParticleGun->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.));
0136 G4double x0, y0, z0;
0137 x0 = 100. * mm;
0138 y0 = 100. * mm;
0139 z0 = -mediumRadius;
0140 while (!(std::sqrt(x0 * x0 + y0 * y0) <= mediumRadius)) {
0141 x0 = CLHEP::RandFlat::shoot(-mediumRadius, mediumRadius);
0142 y0 = CLHEP::RandFlat::shoot(-mediumRadius, mediumRadius);
0143 }
0144 fpParticleGun->SetParticlePosition(G4ThreeVector(x0, y0, z0));
0145 }
0146
0147
0148
0149
0150 else {
0151 G4double cosTheta = 2. * G4UniformRand() - 1;
0152 G4double sinTheta = std::sqrt(1. - cosTheta * cosTheta);
0153 G4double phi = twopi * G4UniformRand();
0154 G4ThreeVector positionStart(mediumRadius * sinTheta * std::cos(phi),
0155 mediumRadius * sinTheta * std::sin(phi), mediumRadius * cosTheta);
0156 fpParticleGun->SetParticlePosition(positionStart);
0157
0158 G4ThreeVector positionDir(boundingXHalfLength * (2. * G4UniformRand() - 1),
0159 boundingYHalfLength * (2. * G4UniformRand() - 1),
0160 boundingZHalfLength * (2. * G4UniformRand() - 1));
0161 fpParticleGun->SetParticleMomentumDirection((positionDir - positionStart).unit());
0162
0163 fGunArea = 4. * pi * mediumRadius * mediumRadius;
0164 }
0165 }
0166 else {
0167 G4cerr << "Bounding slice volume not found!" << G4endl;
0168 G4cerr << "Default particle kinematic used" << G4endl;
0169 }
0170
0171 fpParticleGun->GeneratePrimaryVertex(anEvent);
0172 }