Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QSynchRad.cc
Go to the documentation of this file.
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 // $Id$
28 //
29 // Created by Mikhail Kosov 6-Nov-2009
30 //
31 // --------------------------------------------------------------
32 // Short description: Algorithm of Synchrotron Radiation from PDG
33 // gamma>>1: dI/dw=(8pi/9)*alpha*gamma*F(w/wc), wc=3*gamma^3*c/2/R
34 // F(y)=(9*sqrt(3)/8/pi)*y*int{y,inf}(K_(5/3)(x)dx) (approximated)
35 // N_gamma=[5pi/sqrt(3)]*alpha*gamma; <w>=[8/15/sqrt(3)]*wc
36 // for electrons/positrons: wc(keV)=2.22*[E(GeV)]^3/R(m)
37 // dE per revolution = (4pi/3)*e^2*beta^3*gamma/R
38 // at beta=1, dE(MeV)=.o885*[E(GeV)]^4/R(m)
39 //---------------------------------------------------------------
40 
41 //#define debug
42 //#define pdebug
43 
44 #include "G4QSynchRad.hh"
45 #include "G4PhysicalConstants.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4HadronicDeprecate.hh"
48 
49 
50 // Constructor
52  G4VDiscreteProcess (Name, Type), minGamma(227.), Polarization(0.,0.,1.) {
53  G4HadronicDeprecate("G4QSynchRad");
54 }
55 
56 // Calculates MeanFreePath in GEANT4 internal units
58 {
59  static const G4double coef = 0.4*std::sqrt(3.)/fine_structure_const;
60  const G4DynamicParticle* particle = track.GetDynamicParticle();
61  *cond = NotForced ;
62  G4double gamma = particle->GetTotalEnergy() / particle->GetMass();
63 #ifdef debug
64  G4cout<<"G4QSynchRad::MeanFreePath: gamma = "<<gamma<<G4endl;
65 #endif
66  G4double MFP = DBL_MAX;
67  if( gamma > minGamma ) // For smalle gamma neglect the process
68  {
69  G4double R = GetRadius(track);
70 #ifdef debug
71  G4cout<<"G4QSynchRad::MeanFreePath: Radius = "<<R/meter<<" [m]"<<G4endl;
72 #endif
73  if(R > 0.) MFP= coef*R/gamma;
74  }
75 #ifdef debug
76  G4cout<<"G4QSynchRad::MeanFreePath = "<<MFP/centimeter<<" [cm]"<<G4endl;
77 #endif
78  return MFP;
79 }
80 
82 
83 {
84  static const G4double hc = 1.5 * c_light * hbar_Planck; // E_c=h*w_c=1.5*(hc)*(gamma^3)/R
86  const G4DynamicParticle* particle=track.GetDynamicParticle();
87  G4double gamma = particle->GetTotalEnergy() / particle->GetMass();
88  if(gamma <= minGamma )
89  {
90 #ifdef debug
91  G4cout<<"-Warning-G4QSynchRad::PostStepDoIt is called for small gamma="<<gamma<<G4endl;
92 #endif
93  return G4VDiscreteProcess::PostStepDoIt(track,step);
94  }
95  // Photon energy calculation (E < 8.1*Ec restriction)
96  G4double R = GetRadius(track);
97  if(R <= 0.)
98  {
99 #ifdef debug
100  G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: zero or negativ radius ="
101  <<R/meter<<" [m]"<<G4endl;
102 #endif
103  return G4VDiscreteProcess::PostStepDoIt(track, step);
104  }
105  G4double EPhoton = hc * gamma * gamma * gamma / R; // E_c
106  G4double dd=5.e-8;
107  G4double rnd=G4UniformRand()*(1.+dd);
108  if (rnd < 0.5 ) EPhoton *= .65 * rnd * rnd * rnd;
109  else if(rnd > .997) EPhoton *= 15.-1.03*std::log((1.-rnd)/dd+1.);
110  else
111  {
112  G4double r2=rnd*rnd;
113  G4double dr=1.-rnd;
114  EPhoton*=(2806.+28./rnd)/(1.+500./r2/r2+6500.*(std::sqrt(dr)+28.*dr*dr*dr));
115  }
116 #ifdef debug
117  G4cout<<"G4SynchRad::PostStepDoIt: PhotonEnergy = "<<EPhoton/keV<<" [keV]"<<G4endl;
118 #endif
119  if(EPhoton <= 0.)
120  {
121  G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: zero or negativ photon energy="
122  <<EPhoton/keV<<" [keV]"<<G4endl;
123  return G4VDiscreteProcess::PostStepDoIt(track, step);
124  }
125  G4double kinEn = particle->GetKineticEnergy();
126  G4double newEn = kinEn - EPhoton ;
127  if (newEn > 0.)
128  {
131  }
132  else // Very low probable event
133  {
134  G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: PhotonEnergy > TotalKinEnergy"<<G4endl;
135  EPhoton = kinEn;
139  }
140  G4ThreeVector MomDir = particle->GetMomentumDirection();
141  G4DynamicParticle* Photon = new G4DynamicParticle(G4Gamma::Gamma(), MomDir, EPhoton);
142  Photon->SetPolarization(Polarization.x(), Polarization.y(), Polarization.z());
145  return G4VDiscreteProcess::PostStepDoIt(track,step);
146 }
147 
148 // Revolution Radius in independent units for the particle (general member function)
150 {
151  static const G4double unk = meter*tesla/0.3/gigaelectronvolt;
152  const G4DynamicParticle* particle = track.GetDynamicParticle();
153  G4double z = particle->GetDefinition()->GetPDGCharge();
154  if(z == 0.) return 0.; // --> neutral particle
155  if(z < 0.) z=-z;
157  G4PropagatorInField* Field = transMan->GetPropagatorInField();
158  G4FieldManager* fMan = Field->FindAndSetFieldManager(track.GetVolume());
159  if(!fMan || !fMan->GetDetectorField()) return 0.; // --> no field at all
160  const G4Field* pField = fMan->GetDetectorField();
162  G4double PosArray[3]={position.x(), position.y(), position.z()};
163  G4double BArray[3];
164  pField->GetFieldValue(PosArray, BArray);
165  G4ThreeVector B3D(BArray[0], BArray[1], BArray[2]);
166 #ifdef debug
167  G4cout<<"G4QSynchRad::GetRadius: Pos="<<position/meter<<", B(tesla)="<<B3D/tesla<<G4endl;
168 #endif
169  G4ThreeVector MomDir = particle->GetMomentumDirection();
170  G4ThreeVector Ort = B3D.cross(MomDir);
171  G4double OrtB = Ort.mag(); // not negative (independent units)
172  if(OrtB == 0.) return 0.; // --> along the field line
173  Polarization = Ort/OrtB; // Polarization unit vector
174  G4double mom = particle->GetTotalMomentum(); // Momentum of the particle
175 #ifdef debug
176  G4cout<<"G4QSynchRad::GetRadius: P(GeV)="<<mom/GeV<<", B(tesla)="<<OrtB/tesla<<G4endl;
177 #endif
178  // R [m]= mom [GeV]/(0.3 * z * OrtB [tesla])
179  return mom * unk / z / OrtB;
180 }