G4ModifiedTsai.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:     G4ModifiedTsai
00034 //
00035 // Author:        Andreia Trindade (andreia@lip.pt)
00036 //                Pedro Rodrigues  (psilva@lip.pt)
00037 //                Luis Peralta     (luis@lip.pt)
00038 // 
00039 // Creation date: 21 March 2003
00040 //
00041 // Modifications: 
00042 // 21 Mar 2003 A.Trindade First implementation acording with new design
00043 // 24 Mar 2003 A.Trindade Fix in Tsai generator in order to prevent theta 
00044 //                        generation above pi
00045 // 13 Oct 2010  V.Ivanchenko  Moved to standard and improved comment
00046 // 26.04.2011   V.Grichine    Clean-up of PolarAngle method 
00047 //
00048 // Class Description: 
00049 //
00050 // Bremsstrahlung Angular Distribution Generation 
00051 // suggested by L.Urban (Geant3 manual (1993) Phys211)
00052 // Derived from Tsai distribution (Rev Mod Phys 49,421(1977))
00053 //
00054 // Class Description: End 
00055 //
00056 // -------------------------------------------------------------------
00057 //
00058 
00059 #include "G4ModifiedTsai.hh"
00060 #include "G4PhysicalConstants.hh"
00061 #include "Randomize.hh"
00062 
00063 G4ModifiedTsai::G4ModifiedTsai(const G4String&)
00064   : G4VEmAngularDistribution("AngularGenUrban")
00065 {}    
00066 
00067 G4ModifiedTsai::~G4ModifiedTsai() 
00068 {}
00069 
00070 G4ThreeVector& 
00071 G4ModifiedTsai::SampleDirection(const G4DynamicParticle* dp,
00072                                 G4double, G4int, const G4Material*)
00073 {
00074   // Sample gamma angle (Z - axis along the parent particle).
00075   // Universal distribution suggested by L. Urban (Geant3 manual (1993) 
00076   // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
00077 
00078   G4double uMax = 2*(1. + dp->GetKineticEnergy()/electron_mass_c2);   
00079 
00080   const G4double a1     = 0.625;
00081   const G4double a2     = 1.875;
00082   const G4double border = 0.25;
00083   G4double u;
00084 
00085   do {
00086     u = - std::log(G4UniformRand()*G4UniformRand());
00087 
00088     if ( border > G4UniformRand() ) { u /= a1; }
00089     else                            { u /= a2; }
00090     
00091   } while(u > uMax);
00092 
00093   G4double cost = 1.0 - 2*u*u/(uMax*uMax);
00094   G4double sint = std::sqrt((1 - cost)*(1 + cost));
00095   G4double phi  = CLHEP::twopi*G4UniformRand(); 
00096 
00097   fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
00098   fLocalDirection.rotateUz(dp->GetMomentumDirection());
00099 
00100   return fLocalDirection;
00101 }
00102 
00103 void G4ModifiedTsai::PrintGeneratorInformation() const
00104 {
00105   G4cout << "\n" << G4endl;
00106   G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
00107   G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)" 
00108          << G4endl;
00109   G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n" 
00110          << G4endl;
00111 } 

Generated on Mon May 27 17:48:52 2013 for Geant4 by  doxygen 1.4.7