Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaTransition Class Reference

#include <G4GammaTransition.hh>

Public Member Functions

 G4GammaTransition ()
 
virtual ~G4GammaTransition ()
 
virtual G4FragmentSampleTransition (G4Fragment *nucleus, G4double newExcEnergy, G4double mpRatio, G4int JP1, G4int JP2, G4int MP, G4int shell, G4bool isDiscrete, G4bool isGamma)
 
virtual void SampleDirection (G4Fragment *nuc, G4double ratio, G4int twoJ1, G4int twoJ2, G4int mp)
 
void SetPolarizationFlag (G4bool val)
 
void SetTwoJMAX (G4int val)
 
void SetVerbose (G4int val)
 

Protected Attributes

G4ThreeVector fDirection
 
G4PolarizationTransition fPolTrans
 
G4int fTwoJMAX
 
G4int fVerbose
 

Detailed Description

Definition at line 52 of file G4GammaTransition.hh.

Constructor & Destructor Documentation

◆ G4GammaTransition()

G4GammaTransition::G4GammaTransition ( )
explicit

Definition at line 48 of file G4GammaTransition.cc.

49 : polarFlag(false), fDirection(0.,0.,0.), fTwoJMAX(10), fVerbose(0)
50{}
G4ThreeVector fDirection

◆ ~G4GammaTransition()

G4GammaTransition::~G4GammaTransition ( )
virtual

Definition at line 52 of file G4GammaTransition.cc.

53{}

Member Function Documentation

◆ SampleDirection()

void G4GammaTransition::SampleDirection ( G4Fragment * nuc,
G4double ratio,
G4int twoJ1,
G4int twoJ2,
G4int mp )
virtual

Definition at line 151 of file G4GammaTransition.cc.

153{
154 G4double cosTheta, phi;
155 G4NuclearPolarization* np = nuc->GetNuclearPolarization();
156 if(fVerbose > 2) {
157 G4cout << "G4GammaTransition::SampleDirection : 2J1= " << twoJ1
158 << " 2J2= " << twoJ2 << " ratio= " << ratio
159 << " mp= " << mp << G4endl;
160 G4cout << " Nucleus: " << *nuc << G4endl;
161 }
162 if(nullptr == np) {
163 cosTheta = 2*G4UniformRand() - 1.0;
164 phi = CLHEP::twopi*G4UniformRand();
165
166 } else {
167 // PhotonEvaporation dataset:
168 // The multipolarity number with 1,2,3,4,5,6,7 representing E0,E1,M1,E2,M2,E3,M3
169 // monopole transition and 100*Nx+Ny representing multipolarity transition with
170 // Ny and Ny taking the value 1,2,3,4,5,6,7 referring to E0,E1,M1,E2,M2,E3,M3,..
171 // For example a M1+E2 transition would be written 304.
172 // M1 is the primary transition (L) and E2 is the secondary (L')
173
174 G4double mpRatio = ratio;
175
176 G4int L0 = 0, Lp = 0;
177 if (mp > 99) {
178 L0 = mp/200;
179 Lp = (mp%100)/2;
180 } else {
181 L0 = mp/2;
182 Lp = 0;
183 mpRatio = 0.;
184 }
185 fPolTrans.SampleGammaTransition(np, twoJ1, twoJ2, L0, Lp, mpRatio, cosTheta, phi);
186 }
187
188 G4double sinTheta = std::sqrt((1.-cosTheta)*(1.+cosTheta));
189 fDirection.set(sinTheta*std::cos(phi),sinTheta*std::sin(phi),cosTheta);
190 if(fVerbose > 3) {
191 G4cout << "G4GammaTransition::SampleDirection done: " << fDirection << G4endl;
192 if(np) { G4cout << *np << G4endl; }
193 }
194}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
void set(double x, double y, double z)
G4PolarizationTransition fPolTrans
void SampleGammaTransition(G4NuclearPolarization *np, G4int twoJ1, G4int twoJ2, G4int L0, G4int Lp, G4double mpRatio, G4double &cosTheta, G4double &phi)

Referenced by SampleTransition().

◆ SampleTransition()

G4Fragment * G4GammaTransition::SampleTransition ( G4Fragment * nucleus,
G4double newExcEnergy,
G4double mpRatio,
G4int JP1,
G4int JP2,
G4int MP,
G4int shell,
G4bool isDiscrete,
G4bool isGamma )
virtual

Definition at line 56 of file G4GammaTransition.cc.

65{
66 G4Fragment* result = nullptr;
67 G4double bond_energy = 0.0;
68
69 if (!isGamma) {
70 if(0 <= shell) {
71 G4int Z = nucleus->GetZ_asInt();
72 if(Z <= 104) {
73 G4int idx = std::min(shell, G4AtomicShells::GetNumberOfShells(Z)-1);
74 bond_energy = G4AtomicShells::GetBindingEnergy(Z, idx);
75 }
76 }
77 }
78 G4double etrans = nucleus->GetExcitationEnergy() - newExcEnergy
79 - bond_energy;
80 if(fVerbose > 2) {
81 G4cout << "G4GammaTransition::GenerateGamma - Etrans(MeV)= "
82 << etrans << " Eexnew= " << newExcEnergy
83 << " Ebond= " << bond_energy << G4endl;
84 }
85 if(etrans <= 0.0) {
86 etrans += bond_energy;
87 bond_energy = 0.0;
88 }
89
90 // Do complete Lorentz computation
91 G4LorentzVector lv = nucleus->GetMomentum();
92 G4double mass = nucleus->GetGroundStateMass() + newExcEnergy;
93
94 // select secondary
96
97 if(isGamma) { part = G4Gamma::Gamma(); }
98 else {
99 part = G4Electron::Electron();
100 G4int ne = std::max(nucleus->GetNumberOfElectrons() - 1, 0);
101 nucleus->SetNumberOfElectrons(ne);
102 }
103
104 if(polarFlag && isDiscrete && JP1 <= fTwoJMAX) {
105 SampleDirection(nucleus, mpRatio, JP1, JP2, MP);
106 } else {
108 }
109
110 G4double emass = part->GetPDGMass();
111
112 // 2-body decay in rest frame
113 G4double ecm = lv.mag();
114 G4ThreeVector bst = lv.boostVector();
115 if(!isGamma) { ecm += (CLHEP::electron_mass_c2 - bond_energy); }
116
117 //G4cout << "Ecm= " << ecm << " mass= " << mass << " emass= " << emass << G4endl;
118
119 ecm = std::max(ecm, mass + emass);
120 G4double energy = 0.5*((ecm - mass)*(ecm + mass) + emass*emass)/ecm;
121 G4double mom = (emass > 0.0) ? std::sqrt((energy - emass)*(energy + emass))
122 : energy;
123
124 // emitted gamma or e-
125 G4LorentzVector res4mom(mom * fDirection.x(),
126 mom * fDirection.y(),
127 mom * fDirection.z(), energy);
128 // residual
129 energy = std::max(ecm - energy, mass);
130 lv.set(-mom*fDirection.x(), -mom*fDirection.y(), -mom*fDirection.z(), energy);
131
132 // Lab system transform for short lived level
133 lv.boost(bst);
134
135 // modified primary fragment
136 nucleus->SetExcEnergyAndMomentum(newExcEnergy, lv);
137
138 // gamma or e- are produced
139 res4mom.boost(bst);
140 result = new G4Fragment(res4mom, part);
141
142 //G4cout << " DeltaE= " << e0 - lv.e() - res4mom.e() + emass
143 // << " Emass= " << emass << G4endl;
144 if(fVerbose > 2) {
145 G4cout << "G4GammaTransition::SampleTransition : " << *result << G4endl;
146 G4cout << " Left nucleus: " << *nucleus << G4endl;
147 }
148 return result;
149}
G4ThreeVector G4RandomDirection()
double z() const
double x() const
double y() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void set(double x, double y, double z, double t)
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
static G4int GetNumberOfShells(G4int Z)
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double GetGroundStateMass() const
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
G4int GetZ_asInt() const
void SetNumberOfElectrons(G4int value)
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
G4int GetNumberOfElectrons() const
virtual void SampleDirection(G4Fragment *nuc, G4double ratio, G4int twoJ1, G4int twoJ2, G4int mp)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4double energy(const ThreeVector &p, const G4double m)

◆ SetPolarizationFlag()

void G4GammaTransition::SetPolarizationFlag ( G4bool val)
inline

Definition at line 73 of file G4GammaTransition.hh.

73{ polarFlag = val; };

Referenced by G4PhotonEvaporation::Initialise().

◆ SetTwoJMAX()

void G4GammaTransition::SetTwoJMAX ( G4int val)
inline

Definition at line 75 of file G4GammaTransition.hh.

75{ fTwoJMAX = val; };

Referenced by G4PhotonEvaporation::Initialise().

◆ SetVerbose()

void G4GammaTransition::SetVerbose ( G4int val)
inline

Definition at line 77 of file G4GammaTransition.hh.

77{ fVerbose = val; fPolTrans.SetVerbose(val); };

Referenced by G4PhotonEvaporation::Initialise().

Member Data Documentation

◆ fDirection

G4ThreeVector G4GammaTransition::fDirection
protected

Definition at line 90 of file G4GammaTransition.hh.

Referenced by SampleDirection(), and SampleTransition().

◆ fPolTrans

G4PolarizationTransition G4GammaTransition::fPolTrans
protected

Definition at line 91 of file G4GammaTransition.hh.

Referenced by SampleDirection(), and SetVerbose().

◆ fTwoJMAX

G4int G4GammaTransition::fTwoJMAX
protected

Definition at line 92 of file G4GammaTransition.hh.

Referenced by SampleTransition(), and SetTwoJMAX().

◆ fVerbose

G4int G4GammaTransition::fVerbose
protected

Definition at line 93 of file G4GammaTransition.hh.

Referenced by SampleDirection(), SampleTransition(), and SetVerbose().


The documentation for this class was generated from the following files: