Geant4 Cross Reference |
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 /// \file OpNovice/src/OpNovicePrimaryGeneratorAction.cc 27 /// \brief Implementation of the OpNovicePrimaryGeneratorAction class 28 // 29 // 30 // 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 33 34 #include "OpNovicePrimaryGeneratorAction.hh" 35 36 #include "OpNovicePrimaryGeneratorMessenger.hh" 37 38 #include "G4Event.hh" 39 #include "G4ParticleDefinition.hh" 40 #include "G4ParticleGun.hh" 41 #include "G4ParticleTable.hh" 42 #include "G4SystemOfUnits.hh" 43 #include "Randomize.hh" 44 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 46 OpNovicePrimaryGeneratorAction::OpNovicePrimaryGeneratorAction() 47 : G4VUserPrimaryGeneratorAction(), fParticleGun(nullptr) 48 { 49 G4int n_particle = 1; 50 fParticleGun = new G4ParticleGun(n_particle); 51 // create a messenger for this class 52 fGunMessenger = new OpNovicePrimaryGeneratorMessenger(this); 53 // default kinematic 54 // 55 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable(); 56 G4ParticleDefinition* particle = particleTable->FindParticle("e+"); 57 58 fParticleGun->SetParticleDefinition(particle); 59 fParticleGun->SetParticleTime(0.0 * ns); 60 fParticleGun->SetParticlePosition(G4ThreeVector(0., 0., 0.)); 61 fParticleGun->SetParticleMomentumDirection(G4ThreeVector(1., 0., 0.)); 62 fParticleGun->SetParticleEnergy(500.0 * keV); 63 } 64 65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 66 OpNovicePrimaryGeneratorAction::~OpNovicePrimaryGeneratorAction() 67 { 68 delete fParticleGun; 69 delete fGunMessenger; 70 } 71 72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 73 void OpNovicePrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent) 74 { 75 fParticleGun->GeneratePrimaryVertex(anEvent); 76 } 77 78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 79 void OpNovicePrimaryGeneratorAction::SetOptPhotonPolar() 80 { 81 G4double angle = G4UniformRand() * 360.0 * deg; 82 SetOptPhotonPolar(angle); 83 } 84 85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 86 void OpNovicePrimaryGeneratorAction::SetOptPhotonPolar(G4double angle) 87 { 88 if (fParticleGun->GetParticleDefinition()->GetParticleName() != "opticalphoton") { 89 G4ExceptionDescription ed; 90 ed << "Warning: the particleGun is not an opticalphoton"; 91 G4Exception("OpNovicePrimaryGeneratorAction::SetOptPhotonPolar()", "OpNovice_010", JustWarning, 92 ed); 93 return; 94 } 95 96 G4ThreeVector normal(1., 0., 0.); 97 G4ThreeVector kphoton = fParticleGun->GetParticleMomentumDirection(); 98 G4ThreeVector product = normal.cross(kphoton); 99 G4double modul2 = product * product; 100 101 G4ThreeVector e_perpend(0., 0., 1.); 102 if (modul2 > 0.) e_perpend = (1. / std::sqrt(modul2)) * product; 103 G4ThreeVector e_paralle = e_perpend.cross(kphoton); 104 105 G4ThreeVector polar = std::cos(angle) * e_paralle + std::sin(angle) * e_perpend; 106 fParticleGun->SetParticlePolarization(polar); 107 } 108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 109