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 // 27 /// \file eventgenerator/pythia/pythia8decayer/src/SingleParticleGun.cc 28 /// \brief Implementation of the SingleParticleGun class 29 /// 30 /// \author J. Yarba; FNAL 31 32 #include "SingleParticleGun.hh" 33 34 #include "G4Event.hh" 35 #include "G4ParticleDefinition.hh" 36 #include "G4ParticleGun.hh" 37 #include "G4ParticleTable.hh" 38 #include "G4ThreeVector.hh" 39 #include "Randomize.hh" 40 41 // SPECIAL CASE - needed for treatment of polarization for tau's 42 #include "G4TauMinus.hh" 43 #include "G4TauPlus.hh" 44 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 46 47 SingleParticleGun::SingleParticleGun(const G4String& pname, const double pmom) 48 : G4VUserPrimaryGeneratorAction(), fGun(nullptr), fMomentum(pmom) 49 { 50 int nparts = 1; 51 fGun = new G4ParticleGun(nparts); 52 53 G4ParticleTable* pdt = G4ParticleTable::GetParticleTable(); 54 G4ParticleDefinition* pd = pdt->FindParticle(pname); 55 fGun->SetParticleDefinition(pd); 56 fGun->SetParticlePosition(G4ThreeVector(0., 0., 0.)); 57 double mass = pd->GetPDGMass(); 58 double energy = std::sqrt(fMomentum * fMomentum + mass * mass); 59 fGun->SetParticleEnergy(energy); 60 61 // NOTE: do not set momentum directions 62 // as it'll be generated event-by-event 63 } 64 65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 66 67 SingleParticleGun::~SingleParticleGun() 68 { 69 delete fGun; 70 } 71 72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 73 74 void SingleParticleGun::GeneratePrimaries(G4Event* evt) 75 { 76 if (!fGun) { 77 G4cout << " ParticleGun is NOT set up; BAIL OUT from this event " << G4endl; 78 G4cout << " Check if you are using ctor " 79 << "SinglePartGun(const G4string& partname, const double pmom ) " << G4endl; 80 return; 81 } 82 83 double phi = -1. * CLHEP::pi + CLHEP::twopi * G4UniformRand(); 84 double theta = CLHEP::pi / 4. + CLHEP::halfpi * G4UniformRand(); 85 86 double x = std::sin(theta) * std::cos(phi); 87 double y = std::sin(theta) * std::sin(phi); 88 double z = std::cos(theta); 89 90 fGun->SetParticleMomentumDirection(G4ThreeVector(x, y, z)); 91 92 // SPECIAL CASE: for particle such as tau polarization should be -1 * p3 93 // while for tau_bar it should be +1 94 // 95 // NOTE: one could probably compare by name but that'd mean comparing strings... 96 // 97 if (fGun->GetParticleDefinition() == G4TauMinus::TauMinus()) { 98 fGun->SetParticlePolarization(G4ThreeVector(-1.0 * x, -1.0 * y, -1.0 * z)); 99 } 100 else if (fGun->GetParticleDefinition() == G4TauPlus::TauPlus()) { 101 fGun->SetParticlePolarization(G4ThreeVector(x, y, z)); 102 } 103 104 fGun->GeneratePrimaryVertex(evt); 105 106 return; 107 } 108 109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 110