Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/eRosita/physics/src/G4RDPhotoElectricAngularGeneratorSauterGavrila.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:     G4RDPhotoElectricAngularGeneratorSauterGavrila
 33 //
 34 // Creation date: 10 May 2004
 35 //
 36 // Modifications: 
 37 // 10 May 2003     P. Rodrigues    First implementation acording with new design
 38 //
 39 // Class Description: 
 40 //
 41 // Concrete class for PhotoElectric Electron Angular Distribution Generation 
 42 // This model is a re-implementation of the Photolectric angular distribution
 43 // developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect 
 44 //
 45 // Class Description: End 
 46 //
 47 // -------------------------------------------------------------------
 48 //
 49 //    
 50 
 51 #include "G4RDPhotoElectricAngularGeneratorSauterGavrila.hh"
 52 #include "Randomize.hh"
 53 #include "G4PhysicalConstants.hh"
 54 
 55 //    
 56 
 57 G4RDPhotoElectricAngularGeneratorSauterGavrila::G4RDPhotoElectricAngularGeneratorSauterGavrila(const G4String& name):G4RDVPhotoElectricAngularDistribution(name)
 58 {;}
 59 
 60 //    
 61 
 62 G4RDPhotoElectricAngularGeneratorSauterGavrila::~G4RDPhotoElectricAngularGeneratorSauterGavrila() 
 63 {;}
 64 
 65 //
 66 
 67 G4ThreeVector G4RDPhotoElectricAngularGeneratorSauterGavrila::GetPhotoElectronDirection(const G4ThreeVector& direction, const G4double eKineticEnergy, const G4ThreeVector&, const G4int) const
 68 {
 69 
 70   // Compute Theta distribution of the emitted electron, with respect to the
 71   // incident Gamma.
 72   // The Sauter-Gavrila distribution for the K-shell is used. (adapted from G4PhotoElectricEffect)
 73 
 74   G4double costeta = 1.;
 75   G4double Phi     = twopi * G4UniformRand();
 76   G4double cosphi = std::cos(Phi);
 77   G4double sinphi = std::sin(Phi);
 78   G4double sinteta = 0;
 79   G4double gamma   = 1. + eKineticEnergy/electron_mass_c2;
 80 
 81   if (gamma > 5.) {
 82     G4ThreeVector direction (sinteta*cosphi, sinteta*sinphi, costeta);
 83     return direction;
 84   }
 85 
 86   G4double beta  = std::sqrt(gamma*gamma-1.)/gamma;
 87   G4double b     = 0.5*gamma*(gamma-1.)*(gamma-2);
 88     
 89   G4double rndm,term,greject,grejsup;
 90   if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b);
 91   else            grejsup = gamma*gamma*(1.+b+beta*b);
 92   
 93   do { rndm = 1.-2*G4UniformRand();
 94        costeta = (rndm+beta)/(rndm*beta+1.);
 95        term = 1.-beta*costeta;
 96        greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
 97   } while(greject < G4UniformRand()*grejsup);
 98        
 99 
100   sinteta = std::sqrt(1.-costeta*costeta);
101   G4ThreeVector photoelectrondirection (sinteta*cosphi, sinteta*sinphi, costeta);
102   photoelectrondirection.rotateUz(direction);
103   return photoelectrondirection;
104 }
105 
106 //
107 
108 void G4RDPhotoElectricAngularGeneratorSauterGavrila::PrintGeneratorInformation() const
109 {
110   G4cout << "\n" << G4endl;
111   G4cout << "" << G4endl;
112   G4cout << "Re-implementation of the photolectric angular distribution" << G4endl;
113   G4cout << "developed my M. Maire for the Standard EM Physics G4PhotoElectricEffect" << G4endl;
114   G4cout << "It computes the theta distribution of the emitted electron, with respect to the" << G4endl;
115   G4cout << "incident Gamma, using the Sauter-Gavrila distribution for the K-shell\n" << G4endl;
116 } 
117