G4SauterGavrilaAngularDistribution.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4SauterGavrilaAngularDistribution
00034 //
00035 // Author:     Vladimir Ivanchenko using Michel Maire algorithm
00036 //             developed for Geant3
00037 // 
00038 // Creation date: 23 July 2012
00039 //
00040 //
00041 // -------------------------------------------------------------------
00042 //
00043 
00044 #include "G4SauterGavrilaAngularDistribution.hh"
00045 #include "G4PhysicalConstants.hh"
00046 #include "Randomize.hh"
00047 
00048 G4SauterGavrilaAngularDistribution::G4SauterGavrilaAngularDistribution()
00049   : G4VEmAngularDistribution("AngularGenSauterGavrila")
00050 {}    
00051 
00052 G4SauterGavrilaAngularDistribution::~G4SauterGavrilaAngularDistribution() 
00053 {}
00054 
00055 G4ThreeVector& 
00056 G4SauterGavrilaAngularDistribution::SampleDirection(
00057        const G4DynamicParticle* dp, G4double, G4int, const G4Material*)
00058 {
00059   G4double tau = dp->GetKineticEnergy()/electron_mass_c2;
00060   const G4double taulimit = 30.0;
00061 
00062   if (tau > taulimit) {
00063     fLocalDirection = dp->GetMomentumDirection(); 
00064     // Bugzilla 1120
00065     // SI on 05/09/2010 as suggested by JG 04/09/10 
00066   } else {
00067  
00068     G4double invgamma  = 1.0/(tau + 1.0);
00069     G4double beta      = std::sqrt(tau*(tau + 2.0))*invgamma;
00070     G4double b         = 0.5*tau*(tau*tau - 1.0);
00071     G4double invgamma2 = invgamma*invgamma;
00072    
00073     G4double rndm,term,greject,grejsup,costeta,sint2;
00074     if (tau < 1.) { grejsup = (1.+b-beta*b)/invgamma2; }
00075     else          { grejsup = (1.+b+beta*b)/invgamma2; }
00076 
00077     do { 
00078       rndm    = 1 - 2*G4UniformRand();
00079       costeta = (rndm + beta)/(rndm*beta + 1);
00080       term    = invgamma2/(1 + beta*rndm);
00081       sint2   = (1 - costeta)*(1 + costeta);
00082       greject = sint2*(1 + b*term)/(term*term);
00083 
00084     } while(greject < G4UniformRand()*grejsup);
00085        
00086     G4double sinteta = std::sqrt(sint2);
00087     G4double phi  = CLHEP::twopi*G4UniformRand(); 
00088 
00089     fLocalDirection.set(sinteta*std::cos(phi), sinteta*std::sin(phi), costeta);
00090     fLocalDirection.rotateUz(dp->GetMomentumDirection());
00091   }
00092   return fLocalDirection;
00093 }
00094 
00095 void G4SauterGavrilaAngularDistribution::PrintGeneratorInformation() const
00096 {
00097   G4cout << "\n" << G4endl;
00098   G4cout << "Non-polarized photoelectric effect angular generator." << G4endl;
00099   G4cout << "The Sauter-Gavrila distribution for the K-shell is used."<<G4endl;
00100   G4cout << "Originally developed by M.Maire for Geant3" 
00101          << G4endl;
00102 } 

Generated on Mon May 27 17:49:47 2013 for Geant4 by  doxygen 1.4.7