Geant4  10.02.p02
G4PolarizedGammaTransition.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 // $Id: G4PolarizedGammaTransition.cc 85659 2014-11-03 10:59:10Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 // GEANT 4 class file
30 //
31 // CERN, Geneva, Switzerland
32 //
33 // File name: G4PolarizedGammaTransition
34 //
35 // Author V.Ivanchenko 6 November 2015
36 //
37 // -------------------------------------------------------------------
38 
40 #include "G4AtomicShells.hh"
41 #include "Randomize.hh"
42 #include "G4Gamma.hh"
43 #include "G4Electron.hh"
44 #include "G4LorentzVector.hh"
45 #include "G4SystemOfUnits.hh"
46 #include "G4PhysicalConstants.hh"
48 
51 {
53 }
54 
56 {
57  delete fPolarization;
58 }
59 
60 G4Fragment*
62  G4double newExcEnergy,
63  G4int, size_t shell,
64  G4bool isGamma,
65  G4bool isLongLived)
66 {
67  G4Fragment* result = 0;
68  G4double bond_energy = 0.;
69 
70  if (!isGamma) {
71  G4int Z = nucleus->GetZ_asInt();
72  if(Z <= 100) {
73  G4int idx = (G4int)shell;
75  bond_energy = G4AtomicShells::GetBindingEnergy(Z, idx);
76  }
77  }
78  G4double etrans = nucleus->GetExcitationEnergy() - newExcEnergy
79  - bond_energy;
80  //G4cout << "G4PolarizedGammaTransition::GenerateGamma - Etrans(MeV)= "
81  // << etrans << " Eexnew= " << newExcEnergy << G4endl;
82  if(etrans <= 0.0) { return result; }
83 
84  // Do complete Lorentz computation
85  G4LorentzVector lv = nucleus->GetMomentum();
86  G4double mass = nucleus->GetGroundStateMass() + newExcEnergy;
87 
88  //G4double e0 = lv.e();
89 
90  // select secondary
92 
93  if(isGamma) { part = G4Gamma::Gamma(); }
94  else {
95  part = G4Electron::Electron();
96  G4int ne = nucleus->GetNumberOfElectrons() - 1;
97  if(ne < 0) { ne = 0; }
98  nucleus->SetNumberOfElectrons(ne);
99  lv += G4LorentzVector(0.0,0.0,0.0,
100  CLHEP::electron_mass_c2 - bond_energy);
101  }
102 
103  G4double cosTheta = 1. - 2. * G4UniformRand();
104  G4double sinTheta = std::sqrt(1. - cosTheta * cosTheta);
105  G4double phi = twopi * G4UniformRand();
106 
107  G4double emass = part->GetPDGMass();
108 
109  // 2-body decay in rest frame
110  G4double ecm = lv.mag();
111  G4ThreeVector bst = lv.boostVector();
112 
113  //G4cout << "Ecm= " << ecm << " mass= " << mass << " emass= " << emass
114  // << " isLongLived: " << isLongLived << G4endl;
115 
116  G4double energy = 0.5*((ecm - mass)*(ecm + mass) + emass*emass)/ecm;
117  energy = std::max(energy, emass);
118 
119  G4double mom = std::sqrt((energy - emass)*(energy + emass));
120 
121  G4LorentzVector res4mom(mom * sinTheta * std::cos(phi),
122  mom * sinTheta * std::sin(phi),
123  mom * cosTheta, energy);
124 
125  // Lab system transform for short lived level
126  if(!isLongLived) {
127  res4mom.boost(bst);
128  lv -= res4mom;
129  } else {
130  // In exceptional case sample decay at rest at not correct position
131  // of stopping ion, 4-momentum balance is breaked but gamma energy
132  // is correct
133  lv -= res4mom;
134  G4double E = lv.e();
135  G4double P2= (E - mass)*(E + mass);
136  G4ThreeVector v = lv.vect().unit();
137  G4double p = 0.0;
138  if(P2 > 0.0) { p = std::sqrt(P2); }
139  else { E = mass; }
140  lv.set(v.x()*p, v.y()*p, v.z()*p, E);
141  }
142 
143  // modified primary fragment
144  nucleus->SetMomentum(lv);
145 
146  // gamma or e- are produced
147  result = new G4Fragment(res4mom, part);
148 
149  // G4cout << " DeltaE= " << e0 - lv.e() - res4mom.e() << G4endl;
150  //G4cout << "G4PolarizedGammaTransition::GenerateGamma : " << thePhoton << G4endl;
151  //G4cout << " Left nucleus: " << aNucleus << G4endl;
152  return result;
153 }
CLHEP::Hep3Vector G4ThreeVector
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:379
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:384
G4PolarizationTransition * fPolarization
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
bool G4bool
Definition: G4Types.hh:79
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:289
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:294
static const double twopi
Definition: G4SIunits.hh:75
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:278
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4int GetZ_asInt() const
Definition: G4Fragment.hh:261
static const G4double * P2[nN]
static G4Electron * Electron()
Definition: G4Electron.cc:94
double G4double
Definition: G4Types.hh:76
virtual G4Fragment * SampleTransition(G4Fragment *nucleus, G4double newExcEnergy, G4int deltaS, size_t shell, G4bool isGamma, G4bool isLongLived)
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:273
CLHEP::HepLorentzVector G4LorentzVector
static G4int GetNumberOfShells(G4int Z)