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

Last change on this file since 1346 was 1340, checked in by garnier, 15 years ago

update ti head

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