Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GammaTransition.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: G4GammaTransition.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: G4GammaTransition
34 //
35 // Author V.Ivanchenko 27 February 2015
36 //
37 // -------------------------------------------------------------------
38 
39 #include "G4GammaTransition.hh"
40 #include "G4AtomicShells.hh"
41 #include "Randomize.hh"
42 #include "G4RandomDirection.hh"
43 #include "G4Gamma.hh"
44 #include "G4Electron.hh"
45 #include "G4LorentzVector.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4PhysicalConstants.hh"
48 
50  : polarFlag(false), fDirection(0.,0.,0.), fVerbose(0)
51 {}
52 
54 {}
55 
56 G4Fragment*
58  G4double newExcEnergy,
59  G4double mpRatio,
60  G4int JP1,
61  G4int JP2,
62  G4int MP,
63  G4int shell,
64  G4bool isDiscrete,
65  G4bool isGamma)
66 {
67  G4Fragment* result = nullptr;
68  G4double bond_energy = 0.0;
69 
70  if (!isGamma) {
71  if(0 <= shell) {
72  G4int Z = nucleus->GetZ_asInt();
73  if(Z <= 100) {
74  G4int idx = (G4int)shell;
76  bond_energy = G4AtomicShells::GetBindingEnergy(Z, idx);
77  }
78  }
79  }
80  G4double etrans = nucleus->GetExcitationEnergy() - newExcEnergy
81  - bond_energy;
82  if(fVerbose > 1) {
83  G4cout << "G4GammaTransition::GenerateGamma - Etrans(MeV)= "
84  << etrans << " Eexnew= " << newExcEnergy
85  << " Ebond= " << bond_energy << G4endl;
86  }
87  if(etrans <= 0.0) {
88  etrans += bond_energy;
89  bond_energy = 0.0;
90  }
91 
92  // Do complete Lorentz computation
93  G4LorentzVector lv = nucleus->GetMomentum();
94  G4double mass = nucleus->GetGroundStateMass() + newExcEnergy;
95  //G4double e0 = lv.e();
96 
97  // select secondary
99 
100  if(isGamma) { part = G4Gamma::Gamma(); }
101  else {
102  part = G4Electron::Electron();
103  G4int ne = std::max(nucleus->GetNumberOfElectrons() - 1, 0);
104  nucleus->SetNumberOfElectrons(ne);
105  }
106 
107  if(polarFlag && isDiscrete) {
108  SampleDirection(nucleus, mpRatio, JP1, JP2, MP);
109  } else {
111  }
112 
113  G4double emass = part->GetPDGMass();
114 
115  // 2-body decay in rest frame
116  G4double ecm = lv.mag();
117  G4ThreeVector bst = lv.boostVector();
118  if(!isGamma) { ecm += (CLHEP::electron_mass_c2 - bond_energy); }
119 
120  //G4cout << "Ecm= " << ecm << " mass= " << mass << " emass= " << emass << G4endl;
121 
122  ecm = std::max(ecm, mass + emass);
123  G4double energy = 0.5*((ecm - mass)*(ecm + mass) + emass*emass)/ecm;
124  G4double mom = (emass > 0.0) ? std::sqrt((energy - emass)*(energy + emass))
125  : energy;
126 
127  // emitted gamma or e-
128  G4LorentzVector res4mom(mom * fDirection.x(),
129  mom * fDirection.y(),
130  mom * fDirection.z(), energy);
131  // residual
132  energy = std::max(ecm - energy, mass);
133  lv.set(-mom*fDirection.x(), -mom*fDirection.y(), -mom*fDirection.z(), energy);
134 
135  // Lab system transform for short lived level
136  lv.boost(bst);
137 
138  // modified primary fragment
139  nucleus->SetExcEnergyAndMomentum(newExcEnergy, lv);
140 
141  // gamma or e- are produced
142  res4mom.boost(bst);
143  result = new G4Fragment(res4mom, part);
144 
145  //G4cout << " DeltaE= " << e0 - lv.e() - res4mom.e() + emass
146  // << " Emass= " << emass << G4endl;
147  if(fVerbose > 1) {
148  G4cout << "G4GammaTransition::SampleTransition : " << result << G4endl;
149  G4cout << " Left nucleus: " << nucleus << G4endl;
150  }
151  return result;
152 }
153 
155  G4int twoJ1, G4int twoJ2, G4int mp)
156 {
157  // PhotonEvaporation dataset:
158  // The multipolarity number with 1,2,3,4,5,6,7 representing E0,E1,M1,E2,M2,E3,M3
159  // monopole transition and 100*Nx+Ny representing multipolarity transition with
160  // Ny and Ny taking the value 1,2,3,4,5,6,7 referring to E0,E1,M1,E2,M2,E3,M3,..
161  // For example a M1+E2 transition would be written 304.
162  // M1 is the primary transition (L) and E2 is the secondary (L')
163 
164  G4double mpRatio = ratio;
165 
166  G4int L0 = 0, Lp = 0;
167  if (mp > 99) {
168  L0 = mp/200;
169  Lp = (mp%100)/2;
170  } else {
171  L0 = mp/2;
172  Lp = 0;
173  mpRatio = 0.;
174  }
175 
176  fPolTrans.SetGammaTransitionData(twoJ1, twoJ2, L0, mpRatio, Lp);
178 
179  G4double cosTheta, phi;
180  if(!np) {
181  // initial state is non-polarized - create polarization
182  np = new G4NuclearPolarization();
183  nuc->SetNuclearPolarization(np);
184  cosTheta = 2*G4UniformRand() - 1.0;
185  phi = CLHEP::twopi*G4UniformRand();
186 
187  } else {
188 
189  // initial state is polarized - generate correlation
191  phi = fPolTrans.GenerateGammaPhi(cosTheta, np->GetPolarization());
192  }
193  fPolTrans.UpdatePolarizationToFinalState(cosTheta, phi, nuc);
194 
195  G4double sinTheta = std::sqrt((1.-cosTheta)*(1.+cosTheta));
196  fDirection.set(sinTheta*std::cos(phi),sinTheta*std::sin(phi),cosTheta);
197  if(fVerbose > 1) {
198  G4cout << "G4GammaTransition::SampleDirection : " << fDirection << G4endl;
199  G4cout << "Polarisation : " << *np << G4endl;
200  }
201 }
G4double G4ParticleHPJENDLHEData::G4double result
void set(double x, double y, double z)
Hep3Vector boostVector() const
double x() const
virtual void SampleDirection(G4Fragment *nuc, G4double ratio, G4int twoJ1, G4int twoJ2, G4int mp)
G4ThreeVector fDirection
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:397
G4double GenerateGammaPhi(G4double cosTheta, const POLAR &)
G4double GenerateGammaCosTheta(const POLAR &)
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:402
G4ThreeVector G4RandomDirection()
int G4int
Definition: G4Types.hh:78
double z() const
static constexpr double electron_mass_c2
virtual ~G4GammaTransition()
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
void SetNuclearPolarization(G4NuclearPolarization *)
Definition: G4Fragment.hh:220
double mag() const
std::vector< std::vector< G4complex > > & GetPolarization()
bool G4bool
Definition: G4Types.hh:79
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
HepLorentzVector & boost(double, double, double)
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:288
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetGammaTransitionData(G4int twoJ1, G4int twoJ2, G4int Lbar, G4double delta=0, G4int Lprime=1)
virtual G4Fragment * SampleTransition(G4Fragment *nucleus, G4double newExcEnergy, G4double mpRatio, G4int JP1, G4int JP2, G4int MP, G4int shell, G4bool isDiscrete, G4bool isGamma)
void set(double x, double y, double z, double t)
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)
double y() const
G4PolarizationTransition fPolTrans
void UpdatePolarizationToFinalState(G4double cosTheta, G4double phi, G4Fragment *)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.hh:293
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double twopi
Definition: SystemOfUnits.h:55
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:458
static G4int GetNumberOfShells(G4int Z)