![]() |
|
|||
File indexing completed on 2025-02-23 09:21:11
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 /// \file eventgenerator/pythia/pythia8decayer/src/SingleParticleGun.cc 0028 /// \brief Implementation of the SingleParticleGun class 0029 /// 0030 /// \author J. Yarba; FNAL 0031 0032 #include "SingleParticleGun.hh" 0033 0034 #include "G4Event.hh" 0035 #include "G4ParticleDefinition.hh" 0036 #include "G4ParticleGun.hh" 0037 #include "G4ParticleTable.hh" 0038 #include "G4ThreeVector.hh" 0039 #include "Randomize.hh" 0040 0041 // SPECIAL CASE - needed for treatment of polarization for tau's 0042 #include "G4TauMinus.hh" 0043 #include "G4TauPlus.hh" 0044 0045 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0046 0047 SingleParticleGun::SingleParticleGun(const G4String& pname, const double pmom) 0048 : G4VUserPrimaryGeneratorAction(), fGun(nullptr), fMomentum(pmom) 0049 { 0050 int nparts = 1; 0051 fGun = new G4ParticleGun(nparts); 0052 0053 G4ParticleTable* pdt = G4ParticleTable::GetParticleTable(); 0054 G4ParticleDefinition* pd = pdt->FindParticle(pname); 0055 fGun->SetParticleDefinition(pd); 0056 fGun->SetParticlePosition(G4ThreeVector(0., 0., 0.)); 0057 double mass = pd->GetPDGMass(); 0058 double energy = std::sqrt(fMomentum * fMomentum + mass * mass); 0059 fGun->SetParticleEnergy(energy); 0060 0061 // NOTE: do not set momentum directions 0062 // as it'll be generated event-by-event 0063 } 0064 0065 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0066 0067 SingleParticleGun::~SingleParticleGun() 0068 { 0069 delete fGun; 0070 } 0071 0072 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0073 0074 void SingleParticleGun::GeneratePrimaries(G4Event* evt) 0075 { 0076 if (!fGun) { 0077 G4cout << " ParticleGun is NOT set up; BAIL OUT from this event " << G4endl; 0078 G4cout << " Check if you are using ctor " 0079 << "SinglePartGun(const G4string& partname, const double pmom ) " << G4endl; 0080 return; 0081 } 0082 0083 double phi = -1. * CLHEP::pi + CLHEP::twopi * G4UniformRand(); 0084 double theta = CLHEP::pi / 4. + CLHEP::halfpi * G4UniformRand(); 0085 0086 double x = std::sin(theta) * std::cos(phi); 0087 double y = std::sin(theta) * std::sin(phi); 0088 double z = std::cos(theta); 0089 0090 fGun->SetParticleMomentumDirection(G4ThreeVector(x, y, z)); 0091 0092 // SPECIAL CASE: for particle such as tau polarization should be -1 * p3 0093 // while for tau_bar it should be +1 0094 // 0095 // NOTE: one could probably compare by name but that'd mean comparing strings... 0096 // 0097 if (fGun->GetParticleDefinition() == G4TauMinus::TauMinus()) { 0098 fGun->SetParticlePolarization(G4ThreeVector(-1.0 * x, -1.0 * y, -1.0 * z)); 0099 } 0100 else if (fGun->GetParticleDefinition() == G4TauPlus::TauPlus()) { 0101 fGun->SetParticlePolarization(G4ThreeVector(x, y, z)); 0102 } 0103 0104 fGun->GeneratePrimaryVertex(evt); 0105 0106 return; 0107 } 0108 0109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |