Geant4 10.7.0
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 152 of file G4GammaTransition.cc.

154{
155 G4double cosTheta, phi;
156 G4NuclearPolarization* np = nuc->GetNuclearPolarization();
157 if(fVerbose > 2) {
158 G4cout << "G4GammaTransition::SampleDirection : 2J1= " << twoJ1
159 << " 2J2= " << twoJ2 << " ratio= " << ratio
160 << " mp= " << mp << G4endl;
161 G4cout << " Nucleus: " << *nuc << G4endl;
162 }
163 if(nullptr == np) {
164 cosTheta = 2*G4UniformRand() - 1.0;
165 phi = CLHEP::twopi*G4UniformRand();
166
167 } else {
168 // PhotonEvaporation dataset:
169 // The multipolarity number with 1,2,3,4,5,6,7 representing E0,E1,M1,E2,M2,E3,M3
170 // monopole transition and 100*Nx+Ny representing multipolarity transition with
171 // Ny and Ny taking the value 1,2,3,4,5,6,7 referring to E0,E1,M1,E2,M2,E3,M3,..
172 // For example a M1+E2 transition would be written 304.
173 // M1 is the primary transition (L) and E2 is the secondary (L')
174
175 G4double mpRatio = ratio;
176
177 G4int L0 = 0, Lp = 0;
178 if (mp > 99) {
179 L0 = mp/200;
180 Lp = (mp%100)/2;
181 } else {
182 L0 = mp/2;
183 Lp = 0;
184 mpRatio = 0.;
185 }
186 fPolTrans.SampleGammaTransition(np, twoJ1, twoJ2, L0, Lp, mpRatio, cosTheta, phi);
187 }
188
189 G4double sinTheta = std::sqrt((1.-cosTheta)*(1.+cosTheta));
190 fDirection.set(sinTheta*std::cos(phi),sinTheta*std::sin(phi),cosTheta);
191 if(fVerbose > 3) {
192 G4cout << "G4GammaTransition::SampleDirection done: " << fDirection << G4endl;
193 if(np) { G4cout << *np << G4endl; }
194 }
195}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
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 <= 100) {
73 G4int idx = (G4int)shell;
74 idx = std::min(idx, G4AtomicShells::GetNumberOfShells(Z)-1);
75 bond_energy = G4AtomicShells::GetBindingEnergy(Z, idx);
76 }
77 }
78 }
79 G4double etrans = nucleus->GetExcitationEnergy() - newExcEnergy
80 - bond_energy;
81 if(fVerbose > 2) {
82 G4cout << "G4GammaTransition::GenerateGamma - Etrans(MeV)= "
83 << etrans << " Eexnew= " << newExcEnergy
84 << " Ebond= " << bond_energy << G4endl;
85 }
86 if(etrans <= 0.0) {
87 etrans += bond_energy;
88 bond_energy = 0.0;
89 }
90
91 // Do complete Lorentz computation
92 G4LorentzVector lv = nucleus->GetMomentum();
93 G4double mass = nucleus->GetGroundStateMass() + newExcEnergy;
94
95 // select secondary
97
98 if(isGamma) { part = G4Gamma::Gamma(); }
99 else {
100 part = G4Electron::Electron();
101 G4int ne = std::max(nucleus->GetNumberOfElectrons() - 1, 0);
102 nucleus->SetNumberOfElectrons(ne);
103 }
104
105 if(polarFlag && isDiscrete && JP1 <= fTwoJMAX) {
106 SampleDirection(nucleus, mpRatio, JP1, JP2, MP);
107 } else {
109 }
110
111 G4double emass = part->GetPDGMass();
112
113 // 2-body decay in rest frame
114 G4double ecm = lv.mag();
115 G4ThreeVector bst = lv.boostVector();
116 if(!isGamma) { ecm += (CLHEP::electron_mass_c2 - bond_energy); }
117
118 //G4cout << "Ecm= " << ecm << " mass= " << mass << " emass= " << emass << G4endl;
119
120 ecm = std::max(ecm, mass + emass);
121 G4double energy = 0.5*((ecm - mass)*(ecm + mass) + emass*emass)/ecm;
122 G4double mom = (emass > 0.0) ? std::sqrt((energy - emass)*(energy + emass))
123 : energy;
124
125 // emitted gamma or e-
126 G4LorentzVector res4mom(mom * fDirection.x(),
127 mom * fDirection.y(),
128 mom * fDirection.z(), energy);
129 // residual
130 energy = std::max(ecm - energy, mass);
131 lv.set(-mom*fDirection.x(), -mom*fDirection.y(), -mom*fDirection.z(), energy);
132
133 // Lab system transform for short lived level
134 lv.boost(bst);
135
136 // modified primary fragment
137 nucleus->SetExcEnergyAndMomentum(newExcEnergy, lv);
138
139 // gamma or e- are produced
140 res4mom.boost(bst);
141 result = new G4Fragment(res4mom, part);
142
143 //G4cout << " DeltaE= " << e0 - lv.e() - res4mom.e() + emass
144 // << " Emass= " << emass << G4endl;
145 if(fVerbose > 2) {
146 G4cout << "G4GammaTransition::SampleTransition : " << *result << G4endl;
147 G4cout << " Left nucleus: " << *nucleus << G4endl;
148 }
149 return result;
150}
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:93
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:280
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:275
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
G4int GetZ_asInt() const
Definition: G4Fragment.hh:263
void SetNumberOfElectrons(G4int value)
Definition: G4Fragment.hh:394
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.hh:285
G4int GetNumberOfElectrons() const
Definition: G4Fragment.hh:389
virtual void SampleDirection(G4Fragment *nuc, G4double ratio, G4int twoJ1, G4int twoJ2, G4int mp)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
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: