source: trunk/source/processes/electromagnetic/lowenergy/src/G4Generator2BS.cc @ 1196

Last change on this file since 1196 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 5.1 KB
Line 
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 file
30//
31//
32// File name:     G4Generator2BS
33//
34// Author:        Andreia Trindade (andreia@lip.pt)
35//                Pedro Rodrigues  (psilva@lip.pt)
36//                Luis Peralta     (luis@lip.pt)
37//
38// Creation date: 2 June 2003
39//
40// Modifications:
41// 02 Jun 2003                    First implementation acording with new design
42// 05 Nov 2003  MGP               Fixed std namespace
43// 17 Nov 2003  MGP               Fixed compilation problem on Windows                 
44//
45// Class Description:
46//
47// Concrete base class for Bremsstrahlung Angular Distribution Generation - 2BS Distribution
48//
49// Class Description: End
50//
51// -------------------------------------------------------------------
52//
53//   
54
55#include "G4Generator2BS.hh"
56#include "Randomize.hh"
57//   
58
59G4Generator2BS::G4Generator2BS(const G4String& name):G4VBremAngularDistribution(name)
60{;}
61
62//   
63
64G4Generator2BS::~G4Generator2BS() 
65{;}
66
67//
68
69G4double G4Generator2BS::PolarAngle(const G4double initial_energy,
70                                    const G4double final_energy,
71                                    const G4int Z)
72{
73
74  // Adapted from "Improved bremsstrahlung photon angular sampling in the EGS4 code system"
75  // by Alex F. Bielajew, Rahde Mohan anc Chen-Shou Chui, PIRS-0203
76  // Ionizing Radiation Standards
77  // Institute for National Measurement Standards
78  // National Research Council of Canada
79  // Departement of Medical Physics, Memorial Sloan-Kettering Cancer Center, New York
80
81
82  G4double theta = 0;
83
84  G4double initialTotalEnergy = (initial_energy+electron_mass_c2)/electron_mass_c2;
85  G4double finalTotalEnergy = (final_energy+electron_mass_c2)/electron_mass_c2;
86  EnergyRatio = finalTotalEnergy/initialTotalEnergy;
87  G4double gMaxEnergy = (pi*initialTotalEnergy)*(pi*initialTotalEnergy);
88
89  G4double Zeff = std::sqrt(static_cast<G4double>(Z) * (static_cast<G4double>(Z) + 1.0));
90  z = (0.00008116224*(std::pow(Zeff,0.3333333)));
91
92  // Rejection arguments
93  rejection_argument1 = (1.0+EnergyRatio*EnergyRatio); 
94  rejection_argument2 = - 2.0*EnergyRatio + 3.0*rejection_argument1;
95  rejection_argument3 = ((1-EnergyRatio)/(2.0*initialTotalEnergy*EnergyRatio))*
96    ((1-EnergyRatio)/(2.0*initialTotalEnergy*EnergyRatio));
97
98  // Calculate rejection function at 0, 1 and Emax
99  G4double gfunction0 = RejectionFunction(0);
100  G4double gfunction1 = RejectionFunction(1);
101  G4double gfunctionEmax = RejectionFunction(gMaxEnergy);
102
103
104  // Calculate Maximum value
105  G4double gMaximum = std::max(gfunction0,gfunction1);
106  gMaximum = std::max(gMaximum,gfunctionEmax);
107
108  G4double rand, gfunctionTest, randTest;
109
110  do{
111    rand = G4UniformRand();
112    rand = rand/(1-rand+1.0/gMaxEnergy);
113    gfunctionTest = RejectionFunction(rand);
114    randTest = G4UniformRand();
115
116  }while(randTest > (gfunctionTest/gMaximum));
117
118  theta = std::sqrt(rand)/initialTotalEnergy;
119
120
121  return theta;
122}
123//
124
125G4double G4Generator2BS::RejectionFunction(G4double value) const
126{
127
128  G4double argument = (1+value)*(1+value);
129
130  G4double gfunction = (4+std::log(rejection_argument3+(z/argument)))*
131    ((4*EnergyRatio*value/argument)-rejection_argument1)+rejection_argument2;
132
133  return gfunction;
134
135}
136
137void G4Generator2BS::PrintGeneratorInformation() const
138{
139  G4cout << "\n" << G4endl;
140  G4cout << "Bremsstrahlung Angular Generator is 2BS Generator from 2BS Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
141  G4cout << "Sampling algorithm adapted from PIRS-0203" << G4endl;
142  G4cout << "\n" << G4endl;
143} 
144
Note: See TracBrowser for help on using the repository browser.