Geant4  10.00.p01
G4NumIntTwoBodyAngDst.icc
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: G4NumIntTwoBodyAngDst.icc 71719 2013-06-21 00:01:54Z mkelsey $
27 // Author: Dennis Wright (SLAC)
28 // Date: 28 January 2013
29 //
30 // Description: implementation base for numerically integrated two-body
31 // final state angular distributions in Bertini-style cascade
32 //
33 // 20130219 M. Kelsey: Move to .icc file, required for templated classes
34 // 20130620 Fix Coverity #50295, recursive #include.
35 // 20130924 M. Kelsey -- Use G4Log, G4Exp for CPU speedup
36 
37 #ifndef G4NumIntTwoBodyAngDst_icc
38 #define G4NumIntTwoBodyAngDst_icc 1
39 
40 #include "G4Log.hh"
41 #include "G4Exp.hh"
42 #include "Randomize.hh"
43 
44 
45 template <G4int NKEBINS, G4int NANGLES> G4double
46 G4NumIntTwoBodyAngDst<NKEBINS,NANGLES>::
47 GetCosTheta(const G4double& ekin, const G4double& pcm) const
48 {
49  G4double costh=1.0, slope=1.0; // Buffers for final result
50  G4double xrand = G4UniformRand(); // Random value to sample CDF
51 
52  if (ekin < labKE[nDists-1]) {
53  Interpolate(ekin); // Get intermediate CDF @ energy
54  for (G4int i=1; i<nAngles; i++) {
55  if (xrand < angDist[i]) {
56  slope = (cosBins[i] - cosBins[i-1])/(angDist[i] - angDist[i-1]);
57  costh = (xrand - angDist[i-1])*slope + cosBins[i-1];
58  break;
59  }
60  }
61  } else {
62  slope = 2.*tcoeff*pcm*pcm;
63  costh = G4Log(1.0 - xrand*(1. - G4Exp(2.*slope) ) )/slope - 1.0;
64  }
65  return costh;
66 }
67 
68 
69 // Find normalized indefinite integral of angular distribution at interpolated
70 // kinetic energy
71 
72 template <G4int NKEBINS, G4int NANGLES>
73 void G4NumIntTwoBodyAngDst<NKEBINS,NANGLES>::Interpolate(const G4double& ekin) const
74 {
75  for (G4int j=1; j<nDists; j++) {
76  if (ekin >= labKE[j]) continue; // Find desired pair of cos(theta) CDFs
77 
78  G4double frac = (ekin - labKE[j-1])/(labKE[j] - labKE[j-1]);
79  for (G4int i=0; i < nAngles; i++) {
80  angDist[i] = (1.0 - frac)*angDists[j-1][i] + frac*angDists[j][i];
81  }
82 
83  break;
84  }
85 }
86 
87 #endif /* G4NumIntTwoBodyAngDst_icc */