Geant4  10.00.p02
G4eeToPGammaModel.cc
Go to the documentation of this file.
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 // $Id: G4eeToPGammaModel.cc 66241 2012-12-13 18:34:42Z gunter $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4eeToPGammaModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 25.10.2003
38 //
39 // Modifications:
40 //
41 //
42 // -------------------------------------------------------------------
43 //
44 
45 
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 
49 #include "G4eeToPGammaModel.hh"
50 #include "Randomize.hh"
51 #include "G4PhysicalConstants.hh"
52 #include "G4SystemOfUnits.hh"
53 #include "G4PionZero.hh"
54 #include "G4Eta.hh"
55 #include "G4Gamma.hh"
56 #include "G4DynamicParticle.hh"
57 #include "G4PhysicsVector.hh"
58 #include "G4PhysicsLinearVector.hh"
59 #include "G4eeCrossSections.hh"
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 
63 using namespace std;
64 
66  cross(cr)
67 {
69  if(npart == "pi0") {
70  massR = 782.62*MeV;
71  particle = pi0;
72  } else {
73  massR = 1019.46*MeV;
74  particle = G4Eta::Eta();
75  }
77 }
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82 {}
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {
88  return LowEnergy();
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94 {
95  return massR;
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101 {
102  G4double ee = std::min(HighEnergy(),e);
103  G4double xs;
104  if(particle == pi0) xs = cross->CrossSectionPi0G(ee);
105  else xs = cross->CrossSectionEtaG(ee);
106  return xs;
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110 
112  G4double emax) const
113 {
114  G4double tmin = std::max(emin, ThresholdEnergy());
115  G4double tmax = std::max(tmin, emax);
116  G4int nbins = (G4int)((tmax - tmin)/(5.*MeV));
117  G4PhysicsVector* v = new G4PhysicsLinearVector(emin,emax,nbins);
118  v->SetSpline(true);
119  return v;
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
124 void G4eeToPGammaModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
125  G4double e, const G4ThreeVector& direction)
126 {
127  G4double egam = 0.5*e*(1.0 - massP*massP/(massR*massR));
128  G4double tkin = e - egam - massP;
129  if(tkin < 0.0) tkin = 0.0;
130  G4double cost;
131  do {
132  cost = 2.0*G4UniformRand() - 1.0;
133  } while( 2.0*G4UniformRand() > 1.0 + cost*cost );
134 
135  G4double sint = sqrt(1.0 - cost*cost);
136  G4double phi = twopi * G4UniformRand();
137 
138  G4ThreeVector dir(sint*cos(phi),sint*sin(phi), cost);
139  dir.rotateUz(direction);
140 
141  // create G4DynamicParticle objects
142  G4DynamicParticle* p1 =
143  new G4DynamicParticle(particle,dir,tkin);
144  G4DynamicParticle* p2 =
145  new G4DynamicParticle(G4Gamma::Gamma(),-dir,egam);
146  newp->push_back(p1);
147  newp->push_back(p2);
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, G4double, const G4ThreeVector &)
G4double LowEnergy() const
G4ParticleDefinition * particle
static const double MeV
Definition: G4SIunits.hh:193
virtual G4PhysicsVector * PhysicsVector(G4double, G4double) const
virtual ~G4eeToPGammaModel()
CLHEP::Hep3Vector G4ThreeVector
G4eeToPGammaModel(G4eeCrossSections *, const G4String &)
int G4int
Definition: G4Types.hh:78
static G4Eta * Eta()
Definition: G4Eta.cc:109
void SetSpline(G4bool)
G4double CrossSectionEtaG(G4double)
G4double CrossSectionPi0G(G4double)
#define G4UniformRand()
Definition: Randomize.hh:87
virtual G4double ComputeCrossSection(G4double) const
G4eeCrossSections * cross
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4double PeakEnergy() const
G4double HighEnergy() const
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * pi0
virtual G4double ThresholdEnergy() const