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 #include "PrimaryGeneratorAction.hh" 27 28 #include "ChemistryWorld.hh" 29 #include "DetectorConstruction.hh" 30 31 #include "G4Electron.hh" 32 #include "G4Event.hh" 33 #include "G4ParticleDefinition.hh" 34 #include "G4ParticleGun.hh" 35 #include "G4ParticleTable.hh" 36 #include "G4SystemOfUnits.hh" 37 38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 39 40 PrimaryGeneratorAction::PrimaryGeneratorAction(DetectorConstruction* pDet) 41 : G4VUserPrimaryGeneratorAction(), fpDetector(pDet) 42 { 43 fpMessenger = std::make_unique<PrimaryGeneratorMessenger>(this); 44 fParticleGun = std::make_unique<G4SingleParticleSource>(); 45 G4ParticleDefinition* particle = G4Electron::Definition(); 46 fParticleGun->SetParticleDefinition(particle); 47 fParticleGun->SetNumberOfParticles(1000000); // by user 48 49 auto pPosDist = fParticleGun->GetPosDist(); 50 pPosDist->SetPosDisType("Plane"); 51 pPosDist->SetPosDisShape("Square"); 52 auto faceSiez = fpDetector->GetChemistryWorld()->GetChemistryBoundary()->halfSideLengthInY(); 53 pPosDist->SetCentreCoords(G4ThreeVector(0, 0, -faceSiez)); 54 pPosDist->SetHalfX(faceSiez); 55 pPosDist->SetHalfY(faceSiez); 56 57 auto pAngleDist = fParticleGun->GetAngDist(); 58 pAngleDist->SetParticleMomentumDirection(G4ThreeVector(0., 0., 1.)); 59 60 auto pEnergyDis = fParticleGun->GetEneDist(); 61 pEnergyDis->SetMonoEnergy(0.9999 * MeV); 62 } 63 64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 65 66 void PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent) 67 { 68 G4ParticleDefinition* particle = fParticleGun->GetParticleDefinition(); 69 auto NumberOfParticlesToBeGenerated = fParticleGun->GetNumberOfParticles(); 70 auto pPosDist = fParticleGun->GetPosDist(); 71 auto pAngleDist = fParticleGun->GetAngDist(); 72 auto pEnDist = fParticleGun->GetEneDist(); 73 auto rnd = fParticleGun->GetBiasRndm(); 74 auto charge = particle->GetPDGCharge(); 75 76 for (G4int i = 0; i < NumberOfParticlesToBeGenerated; i++) { 77 auto angle = pAngleDist->GenerateOne(); 78 auto energy = pEnDist->GenerateOne(particle); 79 auto pos = pPosDist->GenerateOne(); 80 auto mass = particle->GetPDGMass(); 81 auto p = new G4PrimaryParticle(particle); 82 auto vertex = new G4PrimaryVertex(pos, 0); 83 p->SetKineticEnergy(energy); 84 p->SetMass(mass); 85 p->SetMomentumDirection(angle); 86 p->SetCharge(charge); 87 p->SetPolarization(0., 0., 0.); 88 89 G4double weight = pEnDist->GetWeight() * rnd->GetBiasWeight(); 90 p->SetWeight(weight); 91 vertex->SetPrimary(p); 92 anEvent->AddPrimaryVertex(vertex); 93 } 94 } 95 96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 97