source: trunk/source/processes/electromagnetic/highenergy/src/G4hBremsstrahlungModel.cc

Last change on this file was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 3.9 KB
Line 
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// $Id: G4hBremsstrahlungModel.cc,v 1.3 2008/07/22 16:15:16 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4hBremsstrahlungModel
35//
36// Author:        Vladimir Ivanchenko on base of G4MuBremsstrahlungModel
37//
38// Creation date: 28.02.2008
39//
40// Modifications:
41//
42
43//
44// Class Description:
45//
46//
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51
52#include "G4hBremsstrahlungModel.hh"
53
54using namespace std;
55
56G4hBremsstrahlungModel::G4hBremsstrahlungModel(const G4ParticleDefinition* p,
57                                               const G4String& nam)
58  : G4MuBremsstrahlungModel(p, nam)
59{}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
63G4hBremsstrahlungModel::~G4hBremsstrahlungModel()
64{}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68G4double G4hBremsstrahlungModel::ComputeDMicroscopicCrossSection(
69                                           G4double tkin,
70                                           G4double Z,
71                                           G4double gammaEnergy)
72//  differential cross section
73{
74  G4double dxsection = 0.;
75
76  if( gammaEnergy > tkin) return dxsection ;
77  //  G4cout << "G4hBremsstrahlungModel m= " << mass
78  //     << "  " << particle->GetParticleName() << G4endl;
79  G4double E = tkin + mass ;
80  G4double v = gammaEnergy/E ;
81  G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
82  G4double rab0=delta*sqrte ;
83
84  G4int iz = G4int(Z);
85  if(iz < 1) iz = 1;
86
87  G4double z13 = 1.0/nist->GetZ13(iz);
88  G4double dn  = mass*nist->GetA27(iz)/(70.*MeV);
89
90  G4double    b = btf;
91  if(1 == iz) b = bh;
92
93  // nucleus contribution logarithm
94  G4double rab1=b*z13;
95  G4double fn=log(rab1/(dn*(electron_mass_c2+rab0*rab1))*
96              (mass+delta*(dn*sqrte-2.))) ;
97  if(fn <0.) fn = 0. ;
98
99  G4double x = 1.0 - v;
100  if(particle->GetPDGSpin() != 0) x += 0.75*v*v;
101
102  dxsection = coeff*x*Z*Z*fn/gammaEnergy;
103
104  return dxsection;
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.