source: trunk/source/processes/electromagnetic/standard/src/G4ModifiedTsai.cc @ 1350

Last change on this file since 1350 was 1350, checked in by garnier, 13 years ago

update to last version 4.9.4

File size: 3.9 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// $Id: G4ModifiedTsai.cc,v 1.1 2010/10/14 15:17:48 vnivanch Exp $
27// GEANT4 tag $Name: emstand-V09-03-24 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4ModifiedTsai
35//
36// Author:        Andreia Trindade (andreia@lip.pt)
37//                Pedro Rodrigues  (psilva@lip.pt)
38//                Luis Peralta     (luis@lip.pt)
39//
40// Creation date: 21 March 2003
41//
42// Modifications:
43// 21 Mar 2003 A.Trindade First implementation acording with new design
44// 24 Mar 2003 A.Trindade Fix in Tsai generator in order to prevent theta
45//                        generation above pi
46// 13 Oct 2010  V.Ivanchenko  Moved to standard and improved comment
47//
48// Class Description:
49//
50// Bremsstrahlung Angular Distribution Generation
51// suggested by L.Urban (Geant3 manual (1993) Phys211)
52// Derived from Tsai distribution (Rev Mod Phys 49,421(1977))
53//
54// Class Description: End
55//
56// -------------------------------------------------------------------
57//
58
59#include "G4ModifiedTsai.hh"
60#include "Randomize.hh"
61
62G4ModifiedTsai::G4ModifiedTsai(const G4String&)
63  : G4VBremAngularDistribution("AngularGenUrban")
64{}   
65
66G4ModifiedTsai::~G4ModifiedTsai() 
67{}
68
69G4double G4ModifiedTsai::PolarAngle(const G4double initial_energy,
70                                    const G4double, // final_energy
71                                    const G4int ) // Z
72{
73  // Sample gamma angle (Z - axis along the parent particle).
74  // Universal distribution suggested by L. Urban (Geant3 manual (1993)
75  // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
76
77  G4double totalEnergy = initial_energy + electron_mass_c2;   
78
79  const G4double a1 = 0.625, a2 = 3.*a1, d = 27.;
80  G4double u, theta = 0;
81
82  do{
83    u = - std::log(G4UniformRand()*G4UniformRand());
84
85    if (9./(9.+d) > G4UniformRand()) { u /= a1; }
86    else                             { u /= a2; }
87   
88    theta = u*electron_mass_c2/totalEnergy;
89  } while(u*electron_mass_c2 > totalEnergy*pi);
90
91  return theta;
92}
93
94void G4ModifiedTsai::PrintGeneratorInformation() const
95{
96  G4cout << "\n" << G4endl;
97  G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
98  G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)" 
99         << G4endl;
100  G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n" 
101         << G4endl;
102} 
Note: See TracBrowser for help on using the repository browser.