Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/src/G4ModifiedTsai.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:     G4ModifiedTsai
 33 //
 34 // Author:        Andreia Trindade (andreia@lip.pt)
 35 //                Pedro Rodrigues  (psilva@lip.pt)
 36 //                Luis Peralta     (luis@lip.pt)
 37 // 
 38 // Creation date: 21 March 2003
 39 //
 40 // Modifications: 
 41 // 21 Mar 2003 A.Trindade First implementation acording with new design
 42 // 24 Mar 2003 A.Trindade Fix in Tsai generator in order to prevent theta 
 43 //                        generation above pi
 44 // 13 Oct 2010  V.Ivanchenko  Moved to standard and improved comment
 45 // 26.04.2011   V.Grichine    Clean-up of PolarAngle method 
 46 //
 47 // Class Description: 
 48 //
 49 // Bremsstrahlung Angular Distribution Generation 
 50 // suggested by L.Urban (Geant3 manual (1993) Phys211)
 51 // Derived from Tsai distribution (Rev Mod Phys 49,421(1977))
 52 //
 53 // Class Description: End 
 54 //
 55 // -------------------------------------------------------------------
 56 //
 57 
 58 #include "G4ModifiedTsai.hh"
 59 #include "Randomize.hh"
 60 #include "G4Log.hh"
 61 #include <CLHEP/Units/PhysicalConstants.h>
 62 
 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 64 
 65 G4ModifiedTsai::G4ModifiedTsai(const G4String&)
 66   : G4VEmAngularDistribution("ModifiedTsai")
 67 {}    
 68 
 69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 70 
 71 G4ModifiedTsai::~G4ModifiedTsai() = default;
 72 
 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 74 
 75 G4ThreeVector& 
 76 G4ModifiedTsai::SampleDirection(const G4DynamicParticle* dp,
 77                                 G4double, G4int, const G4Material*)
 78 {
 79   // Sample gamma angle (Z - axis along the parent particle).
 80   G4double cost = SampleCosTheta(dp->GetKineticEnergy());
 81   G4double sint = std::sqrt((1 - cost)*(1 + cost));
 82   G4double phi  = CLHEP::twopi*G4UniformRand(); 
 83 
 84   fLocalDirection.set(sint*std::cos(phi), sint*std::sin(phi), cost);
 85   fLocalDirection.rotateUz(dp->GetMomentumDirection());
 86 
 87   return fLocalDirection;
 88 }
 89 
 90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 91 
 92 G4double G4ModifiedTsai::SampleCosTheta(G4double kinEnergy)
 93 {
 94   // Universal distribution suggested by L. Urban (Geant3 manual (1993) 
 95   // Phys211) derived from Tsai distribution (Rev Mod Phys 49,421(1977))
 96 
 97   G4double uMax = 2*(1. + kinEnergy/CLHEP::electron_mass_c2);   
 98 
 99   static const G4double a1     = 1.6;
100   static const G4double a2     = a1/3.;
101   static const G4double border = 0.25;
102   G4double u;
103   CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
104 
105   do {
106     G4double uu = -G4Log(rndmEngine->flat()*rndmEngine->flat());
107     u = (border > rndmEngine->flat()) ? uu*a1 : uu*a2;
108     
109     // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
110   } while(u > uMax);
111 
112   return 1.0 - 2.0*u*u/(uMax*uMax);
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116 
117 void G4ModifiedTsai::SamplePairDirections(const G4DynamicParticle* dp,
118             G4double elecKinEnergy,
119             G4double posiKinEnergy,
120             G4ThreeVector& dirElectron,
121             G4ThreeVector& dirPositron,
122             G4int, const G4Material*)
123 {
124   G4double phi  = CLHEP::twopi * G4UniformRand();
125   G4double sinp = std::sin(phi);
126   G4double cosp = std::cos(phi);
127 
128   G4double cost = SampleCosTheta(elecKinEnergy);
129   G4double sint = std::sqrt((1. - cost)*(1. + cost));
130 
131   dirElectron.set(sint*cosp, sint*sinp, cost);
132   dirElectron.rotateUz(dp->GetMomentumDirection());
133 
134   cost = SampleCosTheta(posiKinEnergy);
135   sint = std::sqrt((1. - cost)*(1. + cost));
136 
137   dirPositron.set(-sint*cosp, -sint*sinp, cost);
138   dirPositron.rotateUz(dp->GetMomentumDirection());
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
143 void G4ModifiedTsai::PrintGeneratorInformation() const
144 {
145   G4cout << "\n" << G4endl;
146   G4cout << "Bremsstrahlung Angular Generator is Modified Tsai" << G4endl;
147   G4cout << "Distribution suggested by L.Urban (Geant3 manual (1993) Phys211)" 
148          << G4endl;
149   G4cout << "Derived from Tsai distribution (Rev Mod Phys 49,421(1977)) \n" 
150          << G4endl;
151 } 
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154