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

#include <G4Fancy3DNucleus.hh>

+ Inheritance diagram for G4Fancy3DNucleus:

Public Member Functions

 G4Fancy3DNucleus ()
 
 ~G4Fancy3DNucleus ()
 
void Init (G4int theA, G4int theZ, G4int numberOfLambdas=0)
 
G4bool StartLoop ()
 
G4NucleonGetNextNucleon ()
 
const std::vector< G4Nucleon > & GetNucleons ()
 
G4int GetMassNumber ()
 
G4double GetMass ()
 
G4int GetCharge ()
 
G4int GetNumberOfLambdas ()
 
G4double GetNuclearRadius ()
 
G4double GetNuclearRadius (const G4double maxRelativeDensity)
 
G4double GetOuterRadius ()
 
G4double AddExcitationEnergy (G4double)
 
G4double GetExcitationEnergy ()
 
G4double CoulombBarrier ()
 
void DoLorentzBoost (const G4LorentzVector &theBoost)
 
void DoLorentzBoost (const G4ThreeVector &theBeta)
 
void DoLorentzContraction (const G4LorentzVector &theBoost)
 
void DoLorentzContraction (const G4ThreeVector &theBeta)
 
void CenterNucleons ()
 
void DoTranslation (const G4ThreeVector &theShift)
 
const G4VNuclearDensityGetNuclearDensity () const
 
void SortNucleonsIncZ ()
 
void SortNucleonsDecZ ()
 
- Public Member Functions inherited from G4V3DNucleus
 G4V3DNucleus ()
 
virtual ~G4V3DNucleus ()
 
virtual void Init (G4int theA, G4int theZ, G4int numberOfLambdas=0)=0
 
virtual G4bool StartLoop ()=0
 
virtual G4NucleonGetNextNucleon ()=0
 
virtual const std::vector< G4Nucleon > & GetNucleons ()=0
 
virtual G4int GetMassNumber ()=0
 
virtual G4double GetMass ()=0
 
virtual G4int GetCharge ()=0
 
virtual G4int GetNumberOfLambdas ()=0
 
virtual G4double GetNuclearRadius ()=0
 
virtual G4double GetNuclearRadius (const G4double maxRelativeDensity)=0
 
virtual G4double GetOuterRadius ()=0
 
virtual G4double CoulombBarrier ()=0
 
virtual void DoLorentzBoost (const G4LorentzVector &theBoost)=0
 
virtual void DoLorentzBoost (const G4ThreeVector &theBeta)=0
 
virtual void DoLorentzContraction (const G4LorentzVector &theBoost)=0
 
virtual void DoLorentzContraction (const G4ThreeVector &theBeta)=0
 
virtual void DoTranslation (const G4ThreeVector &theShift)=0
 
virtual const G4VNuclearDensityGetNuclearDensity () const =0
 
virtual void SortNucleonsIncZ ()=0
 
virtual void SortNucleonsDecZ ()=0
 
std::pair< G4double, G4doubleChooseImpactXandY (G4double maxImpact)
 
std::pair< G4double, G4doubleRefetchImpactXandY ()
 

Detailed Description

Definition at line 54 of file G4Fancy3DNucleus.hh.

Constructor & Destructor Documentation

◆ G4Fancy3DNucleus()

G4Fancy3DNucleus::G4Fancy3DNucleus ( )

Definition at line 60 of file G4Fancy3DNucleus.cc.

61 : myA(0), myZ(0), myL(0), theNucleons(250), currentNucleon(-1), theDensity(0),
62 nucleondistance(0.8*fermi),excitationEnergy(0.),
63 places(250), momentum(250), fermiM(250), testSums(250)
64{
65}

◆ ~G4Fancy3DNucleus()

G4Fancy3DNucleus::~G4Fancy3DNucleus ( )

Definition at line 67 of file G4Fancy3DNucleus.cc.

68{
69 if(theDensity) delete theDensity;
70}

Member Function Documentation

◆ AddExcitationEnergy()

G4double G4Fancy3DNucleus::AddExcitationEnergy ( G4double  anE)
inline

Definition at line 138 of file G4Fancy3DNucleus.hh.

139{
140 excitationEnergy +=anE;
141 return excitationEnergy;
142}

◆ CenterNucleons()

void G4Fancy3DNucleus::CenterNucleons ( )

Definition at line 251 of file G4Fancy3DNucleus.cc.

252{
253 G4ThreeVector center;
254
255 for (G4int i=0; i<myA; i++ )
256 {
257 center+=theNucleons[i].GetPosition();
258 }
259 center /= -myA;
260 DoTranslation(center);
261}
int G4int
Definition: G4Types.hh:85
void DoTranslation(const G4ThreeVector &theShift)

Referenced by Init().

◆ CoulombBarrier()

G4double G4Fancy3DNucleus::CoulombBarrier ( )
virtual

Implements G4V3DNucleus.

Definition at line 607 of file G4Fancy3DNucleus.cc.

608{
609 static const G4double cfactor = (1.44/1.14) * MeV;
610 return cfactor*myZ/(1.0 + G4Pow::GetInstance()->Z13(myA));
611}
double G4double
Definition: G4Types.hh:83
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123

◆ DoLorentzBoost() [1/2]

void G4Fancy3DNucleus::DoLorentzBoost ( const G4LorentzVector theBoost)
virtual

Implements G4V3DNucleus.

Definition at line 213 of file G4Fancy3DNucleus.cc.

214{
215 for (G4int i=0; i<myA; i++){
216 theNucleons[i].Boost(theBoost);
217 }
218}

◆ DoLorentzBoost() [2/2]

void G4Fancy3DNucleus::DoLorentzBoost ( const G4ThreeVector theBeta)
virtual

Implements G4V3DNucleus.

Definition at line 220 of file G4Fancy3DNucleus.cc.

221{
222 for (G4int i=0; i<myA; i++){
223 theNucleons[i].Boost(theBeta);
224 }
225}

◆ DoLorentzContraction() [1/2]

void G4Fancy3DNucleus::DoLorentzContraction ( const G4LorentzVector theBoost)
virtual

Implements G4V3DNucleus.

Definition at line 241 of file G4Fancy3DNucleus.cc.

242{
243 if (theBoost.e() !=0 ) {
244 G4ThreeVector beta = theBoost.vect()/theBoost.e();
246 }
247}
Hep3Vector vect() const
void DoLorentzContraction(const G4LorentzVector &theBoost)

Referenced by DoLorentzContraction().

◆ DoLorentzContraction() [2/2]

void G4Fancy3DNucleus::DoLorentzContraction ( const G4ThreeVector theBeta)
virtual

Implements G4V3DNucleus.

Definition at line 227 of file G4Fancy3DNucleus.cc.

228{
229 G4double beta2=theBeta.mag2();
230 if (beta2 > 0) {
231 G4double factor=(1-std::sqrt(1-beta2))/beta2; // (gamma-1)/gamma/beta**2
232 G4ThreeVector rprime;
233 for (G4int i=0; i< myA; i++) {
234 rprime = theNucleons[i].GetPosition() -
235 factor * (theBeta*theNucleons[i].GetPosition()) * theBeta;
236 theNucleons[i].SetPosition(rprime);
237 }
238 }
239}
double mag2() const

◆ DoTranslation()

void G4Fancy3DNucleus::DoTranslation ( const G4ThreeVector theShift)
virtual

Implements G4V3DNucleus.

Definition at line 263 of file G4Fancy3DNucleus.cc.

264{
265 G4ThreeVector tempV;
266 for (G4int i=0; i<myA; i++ )
267 {
268 tempV = theNucleons[i].GetPosition() + theShift;
269 theNucleons[i].SetPosition(tempV);
270 }
271}

Referenced by CenterNucleons().

◆ GetCharge()

G4int G4Fancy3DNucleus::GetCharge ( )
inlinevirtual

Implements G4V3DNucleus.

Definition at line 123 of file G4Fancy3DNucleus.hh.

124{
125 return myZ;
126}

◆ GetExcitationEnergy()

G4double G4Fancy3DNucleus::GetExcitationEnergy ( )
inline

Definition at line 144 of file G4Fancy3DNucleus.hh.

145{
146 return excitationEnergy;
147}

◆ GetMass()

G4double G4Fancy3DNucleus::GetMass ( )
virtual

Implements G4V3DNucleus.

Definition at line 203 of file G4Fancy3DNucleus.cc.

204{
205 if ( myL <= 0 ) return myZ*G4Proton::Proton()->GetPDGMass() +
206 (myA-myZ)*G4Neutron::Neutron()->GetPDGMass() -
207 BindingEnergy();
208 else return G4HyperNucleiProperties::GetNuclearMass(myA, myZ, myL);
209}
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92

◆ GetMassNumber()

G4int G4Fancy3DNucleus::GetMassNumber ( )
inlinevirtual

Implements G4V3DNucleus.

Definition at line 128 of file G4Fancy3DNucleus.hh.

129{
130 return myA;
131}

◆ GetNextNucleon()

G4Nucleon * G4Fancy3DNucleus::GetNextNucleon ( )
virtual

Implements G4V3DNucleus.

Definition at line 138 of file G4Fancy3DNucleus.cc.

139{
140 return ( (currentNucleon>=0 && currentNucleon<myA) ?
141 &theNucleons[currentNucleon++] : 0 );
142}

Referenced by G4LowEIonFragmentation::ApplyYourself().

◆ GetNuclearDensity()

const G4VNuclearDensity * G4Fancy3DNucleus::GetNuclearDensity ( ) const
virtual

Implements G4V3DNucleus.

Definition at line 273 of file G4Fancy3DNucleus.cc.

274{
275 return theDensity;
276}

◆ GetNuclearRadius() [1/2]

G4double G4Fancy3DNucleus::GetNuclearRadius ( )
virtual

Implements G4V3DNucleus.

Definition at line 179 of file G4Fancy3DNucleus.cc.

180{
181 return GetNuclearRadius(0.5);
182}
G4double GetNuclearRadius()

Referenced by GetNuclearRadius().

◆ GetNuclearRadius() [2/2]

G4double G4Fancy3DNucleus::GetNuclearRadius ( const G4double  maxRelativeDensity)
virtual

Implements G4V3DNucleus.

Definition at line 184 of file G4Fancy3DNucleus.cc.

185{
186 return theDensity->GetRadius(maxRelativeDensity);
187}
virtual G4double GetRadius(const G4double maxRelativeDenisty) const =0

◆ GetNucleons()

const std::vector< G4Nucleon > & G4Fancy3DNucleus::GetNucleons ( )
virtual

Implements G4V3DNucleus.

Definition at line 144 of file G4Fancy3DNucleus.cc.

145{
146 return theNucleons;
147}

Referenced by G4MuMinusCapturePrecompound::ApplyYourself().

◆ GetNumberOfLambdas()

G4int G4Fancy3DNucleus::GetNumberOfLambdas ( )
inlinevirtual

Implements G4V3DNucleus.

Definition at line 133 of file G4Fancy3DNucleus.hh.

134{
135 return myL;
136}

◆ GetOuterRadius()

G4double G4Fancy3DNucleus::GetOuterRadius ( )
virtual

Implements G4V3DNucleus.

Definition at line 189 of file G4Fancy3DNucleus.cc.

190{
191 G4double maxradius2=0;
192
193 for (int i=0; i<myA; i++)
194 {
195 if ( theNucleons[i].GetPosition().mag2() > maxradius2 )
196 {
197 maxradius2=theNucleons[i].GetPosition().mag2();
198 }
199 }
200 return std::sqrt(maxradius2)+nucleondistance;
201}

Referenced by G4LowEIonFragmentation::ApplyYourself().

◆ Init()

void G4Fancy3DNucleus::Init ( G4int  theA,
G4int  theZ,
G4int  numberOfLambdas = 0 
)
virtual

Implements G4V3DNucleus.

Definition at line 83 of file G4Fancy3DNucleus.cc.

84{
85 currentNucleon=-1;
86 theNucleons.clear();
87 nucleondistance = 0.8*fermi;
88 places.clear();
89 momentum.clear();
90 fermiM.clear();
91 testSums.clear();
92
93 myZ = theZ;
94 myA = theA;
95 myL = std::max(numberOfLambdas, 0); // Cannot be negative
96 excitationEnergy=0;
97
98 theNucleons.resize(myA); // Pre-loads vector with empty elements
99
100 // For simplicity, we neglect eventual Lambdas in the nucleus as far as the
101 // density of nucler levels and the Fermi level are concerned.
102
103 if(theDensity) delete theDensity;
104 if ( myA < 17 ) {
105 theDensity = new G4NuclearShellModelDensity(myA, myZ);
106 if( myA == 12 ) nucleondistance=0.9*fermi;
107 } else {
108 theDensity = new G4NuclearFermiDensity(myA, myZ);
109 }
110
111 theFermi.Init(myA, myZ);
112
113 ChooseNucleons();
114
115 ChoosePositions();
116
117 if( myA == 12 ) CenterNucleons(); // This would introduce a bias
118
119 ChooseFermiMomenta();
120
121 G4double Ebinding= BindingEnergy()/myA;
122
123 for (G4int aNucleon=0; aNucleon < myA; aNucleon++)
124 {
125 theNucleons[aNucleon].SetBindingEnergy(Ebinding);
126 }
127
128 return;
129}
void Init(G4int anA, G4int aZ)

Referenced by G4MuMinusCapturePrecompound::ApplyYourself(), and G4LowEIonFragmentation::ApplyYourself().

◆ SortNucleonsDecZ()

void G4Fancy3DNucleus::SortNucleonsDecZ ( )
virtual

Implements G4V3DNucleus.

Definition at line 164 of file G4Fancy3DNucleus.cc.

165{
166 if (theNucleons.size() < 2 ) return; // Avoid unnecessary work
168
169 std::reverse(theNucleons.begin(), theNucleons.end());
170}

◆ SortNucleonsIncZ()

void G4Fancy3DNucleus::SortNucleonsIncZ ( )
virtual

Implements G4V3DNucleus.

Definition at line 156 of file G4Fancy3DNucleus.cc.

157{
158 if (theNucleons.size() < 2 ) return; // Avoid unnecesary work
159
160 std::sort(theNucleons.begin(), theNucleons.end(),
162}
bool G4Fancy3DNucleusHelperForSortInZ(const G4Nucleon &nuc1, const G4Nucleon &nuc2)

Referenced by SortNucleonsDecZ().

◆ StartLoop()

G4bool G4Fancy3DNucleus::StartLoop ( )
virtual

Implements G4V3DNucleus.

Definition at line 131 of file G4Fancy3DNucleus.cc.

132{
133 currentNucleon=0;
134 return (theNucleons.size()>0);
135}

Referenced by G4LowEIonFragmentation::ApplyYourself().


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