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 // ------------------------------------------------------------------- 28 // 29 // GEANT4 Class header file 30 // 31 // 32 // File name: G4eeToPGammaModel 33 // 34 // Author: Vladimir Ivanchenko 35 // 36 // Creation date: 25.10.2003 37 // 38 // Modifications: 39 // 40 // 41 // ------------------------------------------------------------------- 42 // 43 44 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 47 48 #include "G4eeToPGammaModel.hh" 49 #include "Randomize.hh" 50 #include "G4PhysicalConstants.hh" 51 #include "G4SystemOfUnits.hh" 52 #include "G4PionZero.hh" 53 #include "G4Eta.hh" 54 #include "G4Gamma.hh" 55 #include "G4DynamicParticle.hh" 56 #include "G4PhysicsVector.hh" 57 #include "G4PhysicsLinearVector.hh" 58 #include "G4eeCrossSections.hh" 59 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 61 62 using namespace std; 63 64 G4eeToPGammaModel::G4eeToPGammaModel(G4eeCrossSections* cr, 65 const G4String& npart, 66 G4double maxkinEnergy, 67 G4double binWidth) 68 : G4Vee2hadrons(cr, 69 npart=="pi0" ? 782.62*MeV:1019.46*MeV, 70 maxkinEnergy, 71 binWidth) 72 { 73 G4cout << "####G4eeToPGammaModel & particle:" << npart 74 << "####" << G4endl; 75 76 pi0 = G4PionZero::PionZero(); 77 if(npart == "pi0") { 78 massR = 782.62*MeV; 79 particle = pi0; 80 } else { 81 massR = 1019.46*MeV; 82 particle = G4Eta::Eta(); 83 } 84 massP = particle->GetPDGMass(); 85 } 86 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 88 89 G4eeToPGammaModel::~G4eeToPGammaModel() 90 {} 91 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 93 94 G4double G4eeToPGammaModel::PeakEnergy() const 95 { 96 return massR; 97 } 98 99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 100 101 G4double G4eeToPGammaModel::ComputeCrossSection(G4double e) const 102 { 103 G4double xs; 104 if(particle == pi0) xs = cross->CrossSectionPi0G(e); 105 else xs = cross->CrossSectionEtaG(e); 106 return xs; 107 } 108 109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 110 111 void G4eeToPGammaModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp, 112 G4double e, const G4ThreeVector& direction) 113 { 114 G4double egam = 0.5*e*(1.0 - massP*massP/(massR*massR)); 115 G4double tkin = e - egam - massP; 116 if(tkin < 0.0) tkin = 0.0; 117 G4double cost; 118 do { 119 cost = 2.0*G4UniformRand() - 1.0; 120 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko 121 } while( 2.0*G4UniformRand() > 1.0 + cost*cost ); 122 123 G4double sint = sqrt(1.0 - cost*cost); 124 G4double phi = twopi * G4UniformRand(); 125 126 G4ThreeVector dir(sint*cos(phi),sint*sin(phi), cost); 127 dir.rotateUz(direction); 128 129 // create G4DynamicParticle objects 130 G4DynamicParticle* p1 = 131 new G4DynamicParticle(particle,dir,tkin); 132 G4DynamicParticle* p2 = 133 new G4DynamicParticle(G4Gamma::Gamma(),-dir,egam); 134 newp->push_back(p1); 135 newp->push_back(p2); 136 } 137 138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 139 140