Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GammaTransition.cc
이 파일의 문서화 페이지로 가기
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.), fTwoJMAX(10), 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 
96  // select secondary
98 
99  if(isGamma) { part = G4Gamma::Gamma(); }
100  else {
101  part = G4Electron::Electron();
102  G4int ne = std::max(nucleus->GetNumberOfElectrons() - 1, 0);
103  nucleus->SetNumberOfElectrons(ne);
104  }
105 
106  if(polarFlag && isDiscrete && JP1 <= fTwoJMAX) {
107  SampleDirection(nucleus, mpRatio, JP1, JP2, MP);
108  } else {
110  }
111 
112  G4double emass = part->GetPDGMass();
113 
114  // 2-body decay in rest frame
115  G4double ecm = lv.mag();
116  G4ThreeVector bst = lv.boostVector();
117  if(!isGamma) { ecm += (CLHEP::electron_mass_c2 - bond_energy); }
118 
119  //G4cout << "Ecm= " << ecm << " mass= " << mass << " emass= " << emass << G4endl;
120 
121  ecm = std::max(ecm, mass + emass);
122  G4double energy = 0.5*((ecm - mass)*(ecm + mass) + emass*emass)/ecm;
123  G4double mom = (emass > 0.0) ? std::sqrt((energy - emass)*(energy + emass))
124  : energy;
125 
126  // emitted gamma or e-
127  G4LorentzVector res4mom(mom * fDirection.x(),
128  mom * fDirection.y(),
129  mom * fDirection.z(), energy);
130  // residual
131  energy = std::max(ecm - energy, mass);
132  lv.set(-mom*fDirection.x(), -mom*fDirection.y(), -mom*fDirection.z(), energy);
133 
134  // Lab system transform for short lived level
135  lv.boost(bst);
136 
137  // modified primary fragment
138  nucleus->SetExcEnergyAndMomentum(newExcEnergy, lv);
139 
140  // gamma or e- are produced
141  res4mom.boost(bst);
142  result = new G4Fragment(res4mom, part);
143 
144  //G4cout << " DeltaE= " << e0 - lv.e() - res4mom.e() + emass
145  // << " Emass= " << emass << G4endl;
146  if(fVerbose > 1) {
147  G4cout << "G4GammaTransition::SampleTransition : " << *result << G4endl;
148  G4cout << " Left nucleus: " << *nucleus << G4endl;
149  }
150  return result;
151 }
152 
154  G4int twoJ1, G4int twoJ2, G4int mp)
155 {
156  G4double cosTheta, phi;
158  if(fVerbose > 1) {
159  G4cout << "G4GammaTransition::SampleDirection : 2J1= " << twoJ1
160  << " 2J2= " << twoJ2 << " ratio= " << ratio
161  << " mp= " << mp << G4endl;
162  G4cout << " Nucleus: " << *nuc << G4endl;
163  }
164  if(nullptr == np) {
165  cosTheta = 2*G4UniformRand() - 1.0;
166  phi = CLHEP::twopi*G4UniformRand();
167 
168  } else {
169  // PhotonEvaporation dataset:
170  // The multipolarity number with 1,2,3,4,5,6,7 representing E0,E1,M1,E2,M2,E3,M3
171  // monopole transition and 100*Nx+Ny representing multipolarity transition with
172  // Ny and Ny taking the value 1,2,3,4,5,6,7 referring to E0,E1,M1,E2,M2,E3,M3,..
173  // For example a M1+E2 transition would be written 304.
174  // M1 is the primary transition (L) and E2 is the secondary (L')
175 
176  G4double mpRatio = ratio;
177 
178  G4int L0 = 0, Lp = 0;
179  if (mp > 99) {
180  L0 = mp/200;
181  Lp = (mp%100)/2;
182  } else {
183  L0 = mp/2;
184  Lp = 0;
185  mpRatio = 0.;
186  }
187  fPolTrans.SampleGammaTransition(np, twoJ1, twoJ2, L0, Lp, mpRatio, cosTheta, phi);
188  }
189 
190  G4double sinTheta = std::sqrt((1.-cosTheta)*(1.+cosTheta));
191  fDirection.set(sinTheta*std::cos(phi),sinTheta*std::sin(phi),cosTheta);
192  if(fVerbose > 1) {
193  G4cout << "G4GammaTransition::SampleDirection done: " << fDirection << G4endl;
194  if(np) { G4cout << *np << G4endl; }
195  }
196 }
void set(double x, double y, double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:300
G4ThreeVector fDirection
void SampleGammaTransition(G4NuclearPolarization *np, G4int twoJ1, G4int twoJ2, G4int L0, G4int Lp, G4double mpRatio, G4double &cosTheta, G4double &phi)
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:456
#define G4endl
Definition: G4ios.hh:61
double z() const
G4ThreeVector G4RandomDirection()
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:281
G4double GetPDGMass() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Float_t Z
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.hh:286
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
static constexpr double electron_mass_c2
G4PolarizationTransition fPolTrans
#define G4UniformRand()
Definition: Randomize.hh:53
double energy
Definition: plottest35.C:25
static G4int GetNumberOfShells(G4int Z)
void set(double x, double y, double z, double t)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4int GetZ_asInt() const
Definition: G4Fragment.hh:264
G4double G4ParticleHPJENDLHEData::G4double result
static constexpr double twopi
Definition: SystemOfUnits.h:55
virtual G4Fragment * SampleTransition(G4Fragment *nucleus, G4double newExcEnergy, G4double mpRatio, G4int JP1, G4int JP2, G4int MP, G4int shell, G4bool isDiscrete, G4bool isGamma)
double mag() const
virtual ~G4GammaTransition()
int G4int
Definition: G4Types.hh:78
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:390
virtual void SampleDirection(G4Fragment *nuc, G4double ratio, G4int twoJ1, G4int twoJ2, G4int mp)
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
G4GLOB_DLL std::ostream G4cout
double x() const
Hep3Vector boostVector() const
double y() const
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:276
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:395
HepLorentzVector & boost(double, double, double)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments