Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/src/G4Generator2BS.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // 12 Oct 2010  V.Ivanchenko Moved RejectionFunction inline, use G4Pow to speadup
 45 // 09 May 2011  L.Pandola    Initialize private members, to avoid Coverity warning
 46 //
 47 // Class Description: 
 48 //
 49 // Concrete base class for Bremsstrahlung Angular Distribution Generation 
 50 // 2BS Distribution
 51 //
 52 // Class Description: End 
 53 //
 54 // -------------------------------------------------------------------
 55 //
 56 
 57 #include "G4Generator2BS.hh"
 58 #include "Randomize.hh"   
 59 #include "G4PhysicalConstants.hh"
 60 #include "G4SystemOfUnits.hh"
 61 #include "G4Pow.hh"   
 62 
 63 G4Generator2BS::G4Generator2BS(const G4String&)
 64   : G4VEmAngularDistribution("AngularGen2BS"),fz(1.),ratio(1.),
 65     ratio1(1.),ratio2(1.),delta(0.)
 66 {
 67   g4pow = G4Pow::GetInstance();
 68   nwarn = 0;
 69 }
 70 
 71 G4Generator2BS::~G4Generator2BS() 
 72 {}
 73 
 74 G4ThreeVector& G4Generator2BS::SampleDirection(const G4DynamicParticle* dp,
 75                  G4double final_energy,
 76                  G4int Z,
 77                  const G4Material*)
 78 {
 79   // Adapted from "Improved bremsstrahlung photon angular sampling in the EGS4 code system"
 80   // by Alex F. Bielajew, Rahde Mohan anc Chen-Shou Chui, PIRS-0203
 81   // Ionizing Radiation Standards
 82   // Institute for National Measurement Standards 
 83   // National Research Council of Canada
 84   // Departement of Medical Physics, Memorial Sloan-Kettering Cancer Center, New York
 85 
 86   G4double energy = dp->GetTotalEnergy();
 87   ratio = final_energy/energy;
 88   ratio1 = (1 + ratio)*(1 + ratio);
 89   ratio2 = 1 + ratio*ratio;
 90 
 91   G4double gamma = energy/electron_mass_c2;
 92   G4double beta  = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
 93 
 94   // VI speadup
 95   fz = 0.00008116224*g4pow->Z13(Z)*g4pow->Z13(Z+1);
 96 
 97   // majoranta
 98   G4double ymax = 2*beta*(1 + beta)*gamma*gamma;
 99   G4double gMax = RejectionFunction(0.0);
100   gMax = std::max(gMax,RejectionFunction(ymax));
101 
102   G4double y, gfun;
103 
104   do{
105     G4double q = G4UniformRand();
106     y = q*ymax/(1 + ymax*(1 - q));
107     gfun = RejectionFunction(y);
108 
109     // violation point
110     if(gfun > gMax && nwarn >= 20) {
111       ++nwarn;
112       G4cout << "### WARNING in G4Generator2BS: Etot(MeV)= " << energy/MeV 
113        << "  Egamma(MeV)" << (energy - final_energy)/MeV
114        << " gMax= " << gMax << "  < " << gfun
115        << "  results are not reliable!" 
116        << G4endl;
117       if(20 == nwarn) { 
118   G4cout << "   WARNING in G4Generator2BS is closed" << G4endl; 
119       }
120     }
121 
122   } while(G4UniformRand()*gMax > gfun || y > ymax);
123   
124   G4double cost = 1 - 2*y/ymax;
125   G4double sint = std::sqrt((1 - cost)*(1 + cost));
126   G4double phi  = twopi*G4UniformRand(); 
127 
128   fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi),cost);
129   fLocalDirection.rotateUz(dp->GetMomentumDirection());
130 
131   return fLocalDirection;
132 }
133 
134 void G4Generator2BS::PrintGeneratorInformation() const
135 {
136   G4cout << "\n" << G4endl;
137   G4cout << "Bremsstrahlung Angular Generator is 2BS Generator "
138    << "from 2BS Koch & Motz distribution (Rev Mod Phys 31(4), 920 (1959))" << G4endl;
139   G4cout << "Sampling algorithm adapted from PIRS-0203" << G4endl;
140   G4cout << "\n" << G4endl;
141 } 
142 
143