G4DipBustGenerator.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:     G4DipBustGenerator
00034 //
00035 // Author:  Vladimir Grichine     
00036 // 
00037 // Creation date: 17 May 2011
00038 //
00039 // Modifications: 
00040 // 
00041 //
00042 // Class Description: 
00043 //
00044 // Bremsstrahlung Angular Distribution Generation 
00045 // suggested  the dipole approximation in the rest frame of electron 
00046 // busted in the laboratory frame.
00047 //
00048 // Class Description: End 
00049 //
00050 // -------------------------------------------------------------------
00051 //
00052 
00053 #include "G4DipBustGenerator.hh"
00054 #include "G4PhysicalConstants.hh"
00055 #include "Randomize.hh"
00056 
00057 G4DipBustGenerator::G4DipBustGenerator(const G4String&)
00058   : G4VEmAngularDistribution("DipBustGen")
00059 {}    
00060 
00061 G4DipBustGenerator::~G4DipBustGenerator() 
00062 {}
00063 
00064 G4ThreeVector& 
00065 G4DipBustGenerator::SampleDirection(const G4DynamicParticle* dp,
00066                                     G4double, G4int, const G4Material*)
00067 {
00068   G4double c, cosTheta, delta, cofA, signc = 1., a, power = 1./3.;
00069 
00070   G4double eTkin = dp->GetKineticEnergy();
00071 
00072   c = 4. - 8.*G4UniformRand();
00073   a = c;
00074  
00075   if( c < 0. )
00076   {
00077     signc = -1.;
00078     a     = -c;
00079   }
00080   delta  = std::sqrt(a*a+4.);
00081   delta += a;
00082   delta *= 0.5; 
00083 
00084   cofA = -signc*std::pow(delta, power);
00085 
00086   cosTheta = cofA - 1./cofA;
00087 
00088   G4double tau = eTkin/electron_mass_c2;
00089   G4double beta = std::sqrt(tau*(tau + 2.))/(tau + 1.);
00090 
00091   cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
00092 
00093   G4double sinTheta = std::sqrt((1 - cosTheta)*(1 + cosTheta));
00094   G4double phi  = twopi*G4UniformRand(); 
00095 
00096   fLocalDirection.set(sinTheta*std::cos(phi), sinTheta*std::sin(phi),cosTheta);
00097   fLocalDirection.rotateUz(dp->GetMomentumDirection());
00098 
00099   return fLocalDirection;
00100 
00101 }
00102 
00103 G4double G4DipBustGenerator::PolarAngle(const G4double eTkin,
00104                                     const G4double, // final_energy
00105                                     const G4int ) // Z
00106 {
00107   G4double c, cosTheta, delta, cofA, signc = 1., a, power = 1./3.;
00108   G4double gamma, beta, theta;
00109 
00110   c = 4. - 8.*G4UniformRand();
00111   a = c;
00112  
00113   if( c < 0. )
00114   {
00115     signc = -1.;
00116     a     = -c;
00117   }
00118   delta  = std::sqrt(a*a+4.);
00119   delta += a;
00120   delta *= 0.5; 
00121 
00122   cofA = -signc*std::pow(delta, power);
00123 
00124   cosTheta = cofA - 1./cofA;
00125 
00126   gamma = 1. + eTkin/electron_mass_c2;
00127   beta = std::sqrt(1. - 1./gamma/gamma);
00128 
00129   cosTheta = (cosTheta + beta)/(1 + cosTheta*beta);
00130 
00131   theta = std::acos(cosTheta);
00132 
00133   if( theta < 0. )  theta = 0.;
00134   if( theta > pi )  theta = pi;
00135   // G4cout <<"theta = "<<theta<<"; ";
00136 
00137   return theta;
00138 }
00139 
00140 void G4DipBustGenerator::PrintGeneratorInformation() const
00141 {
00142   G4cout << "\n" << G4endl;
00143   G4cout << "Angular Generator based on classical formula from" << G4endl;
00144   G4cout << "J.D. Jackson, Classical Electrodynamics, Wiley, New York 1975" 
00145          << G4endl;
00146 } 

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