Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Fancy3DNucleus.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// ------------------------------------------------------------
27// GEANT 4 class implementation file
28//
29// ---------------- G4Fancy3DNucleus ----------------
30// by Gunter Folger, May 1998.
31// class for a 3D nucleus, arranging nucleons in space and momentum.
32// ------------------------------------------------------------
33// 20110805 M. Kelsey -- Remove C-style array (pointer) of G4Nucleons,
34// make vector a container of objects. Move Helper class
35// to .hh. Move testSums, places, momentum and fermiM to
36// class data members for reuse.
37
38#include <algorithm>
39
40#include "G4Fancy3DNucleus.hh"
44#include "G4NucleiProperties.hh"
45#include "G4Nucleon.hh"
46#include "G4SystemOfUnits.hh"
47#include "Randomize.hh"
48#include "G4ios.hh"
50
51
53 : myA(0), myZ(0), theNucleons(250), currentNucleon(-1), theDensity(0),
54 nucleondistance(0.8*fermi), places(250), momentum(250), fermiM(250),
55 testSums(250)
56{
57//G4cout <<"G4Fancy3DNucleus::G4Fancy3DNucleus()"<<G4endl;
58}
59
61{
62 if(theDensity) delete theDensity;
63}
64
65#if defined(NON_INTEGER_A_Z)
67{
68 G4int intZ = G4int(theZ);
69 G4int intA= ( G4UniformRand()>theA-G4int(theA) ) ? G4int(theA) : G4int(theA)+1;
70 // forward to integer Init()
71 Init(intA, intZ);
72
73}
74#endif
75
77{
78// G4cout << "G4Fancy3DNucleus::Init(theA, theZ) called"<<G4endl;
79 currentNucleon=-1;
80 theNucleons.clear();
81
82 myZ = theZ;
83 myA= theA;
84
85 theNucleons.resize(myA); // Pre-loads vector with empty elements
86
87// G4cout << "myA, myZ" << myA << ", " << myZ << G4endl;
88
89 if(theDensity) delete theDensity;
90 if ( myA < 17 ) {
91 theDensity = new G4NuclearShellModelDensity(myA, myZ);
92 } else {
93 theDensity = new G4NuclearFermiDensity(myA, myZ);
94 }
95
96 theFermi.Init(myA, myZ);
97
98 ChooseNucleons();
99
100 ChoosePositions();
101
102// CenterNucleons(); // This would introduce a bias
103
104 ChooseFermiMomenta();
105
106 G4double Ebinding= BindingEnergy()/myA;
107
108 for (G4int aNucleon=0; aNucleon < myA; aNucleon++)
109 {
110 theNucleons[aNucleon].SetBindingEnergy(Ebinding);
111 }
112
113
114 return;
115}
116
118{
119 currentNucleon=0;
120 return (theNucleons.size()>0);
121}
122
123// Returns by pointer; null pointer indicates end of loop
125{
126 return ( (currentNucleon>=0 && currentNucleon<myA) ?
127 &theNucleons[currentNucleon++] : 0 );
128}
129
130const std::vector<G4Nucleon> & G4Fancy3DNucleus::GetNucleons()
131{
132 return theNucleons;
133}
134
135
136// Class-scope function to sort nucleons by Z coordinate
138{
139 return nuc1.GetPosition().z() < nuc2.GetPosition().z();
140}
141
142void G4Fancy3DNucleus::SortNucleonsIncZ() // on increased Z-coordinates Uzhi 29.08.08
143{
144 if (theNucleons.size() < 2 ) return; // Avoid unnecesary work
145
146 std::sort(theNucleons.begin(), theNucleons.end(),
148}
149
150void G4Fancy3DNucleus::SortNucleonsDecZ() // on decreased Z-coordinates Uzhi 29.08.08
151{
152 if (theNucleons.size() < 2 ) return; // Avoid unnecessary work
154
155 std::reverse(theNucleons.begin(), theNucleons.end());
156}
157
158
159G4double G4Fancy3DNucleus::BindingEnergy()
160{
162}
163
164
166{
167 return GetNuclearRadius(0.5);
168}
169
171{
172 return theDensity->GetRadius(maxRelativeDensity);
173}
174
176{
177 G4double maxradius2=0;
178
179 for (int i=0; i<myA; i++)
180 {
181 if ( theNucleons[i].GetPosition().mag2() > maxradius2 )
182 {
183 maxradius2=theNucleons[i].GetPosition().mag2();
184 }
185 }
186 return std::sqrt(maxradius2)+nucleondistance;
187}
188
190{
191 return myZ*G4Proton::Proton()->GetPDGMass() +
192 (myA-myZ)*G4Neutron::Neutron()->GetPDGMass() -
193 BindingEnergy();
194}
195
196
197
199{
200 for (G4int i=0; i<myA; i++){
201 theNucleons[i].Boost(theBoost);
202 }
203}
204
206{
207 for (G4int i=0; i<myA; i++){
208 theNucleons[i].Boost(theBeta);
209 }
210}
211
213{
214 G4double factor=(1-std::sqrt(1-theBeta.mag2()))/theBeta.mag2(); // (gamma-1)/gamma/beta**2
215 G4ThreeVector rprime;
216 for (G4int i=0; i< myA; i++) {
217 rprime = theNucleons[i].GetPosition() -
218 factor * (theBeta*theNucleons[i].GetPosition()) * theBeta;
219 theNucleons[i].SetPosition(rprime);
220 }
221}
222
224{
225 G4ThreeVector beta = theBoost.vect()/theBoost.e();
226 // DoLorentzBoost(beta);
228}
229
230
231
233{
234 G4ThreeVector center;
235
236 for (G4int i=0; i<myA; i++ )
237 {
238 center+=theNucleons[i].GetPosition();
239 }
240 center /= -myA;
241 DoTranslation(center);
242}
243
245{
246 G4ThreeVector tempV;
247 for (G4int i=0; i<myA; i++ )
248 {
249 tempV = theNucleons[i].GetPosition() + theShift;
250 theNucleons[i].SetPosition(tempV);
251 }
252}
253
255{
256 return theDensity;
257}
258
259//----------------------- private Implementation Methods-------------
260
261void G4Fancy3DNucleus::ChooseNucleons()
262{
263 G4int protons=0,nucleons=0;
264
265 while (nucleons < myA )
266 {
267 if ( protons < myZ && G4UniformRand() < (G4double)(myZ-protons)/(G4double)(myA-nucleons) )
268 {
269 protons++;
270 theNucleons[nucleons++].SetParticleType(G4Proton::Proton());
271 }
272 else if ( (nucleons-protons) < (myA-myZ) )
273 {
274 theNucleons[nucleons++].SetParticleType(G4Neutron::Neutron());
275 }
276 else G4cout << "G4Fancy3DNucleus::ChooseNucleons not efficient" << G4endl;
277 }
278 return;
279}
280
281void G4Fancy3DNucleus::ChoosePositions()
282{
283 G4int i=0;
284 G4ThreeVector aPos, delta;
285 G4bool freeplace;
286 static G4double nd2 = sqr(nucleondistance);
287 G4double maxR=GetNuclearRadius(0.001); // there are no nucleons at a
288 // relative Density of 0.01
289 G4int jr=0;
290 G4int jx,jy;
291 G4double arand[600];
292 G4double *prand=arand;
293
294 places.clear(); // Reset data buffer
295
296 while ( i < myA )
297 {
298 do
299 {
300 if ( jr < 3 )
301 {
302 jr=std::min(600,9*(myA - i));
303 CLHEP::RandFlat::shootArray(jr, prand );
304 }
305 jx=--jr;
306 jy=--jr;
307 aPos.set((2*arand[jx]-1.), (2*arand[jy]-1.), (2*arand[--jr]-1.));
308 } while (aPos.mag2() > 1. );
309 aPos *=maxR;
310 G4double density=theDensity->GetRelativeDensity(aPos);
311 if (G4UniformRand() < density)
312 {
313 freeplace= true;
314 std::vector<G4ThreeVector>::iterator iplace;
315 for( iplace=places.begin(); iplace!=places.end() && freeplace;++iplace)
316 {
317 delta = *iplace - aPos;
318 freeplace= delta.mag2() > nd2;
319 }
320
321 if ( freeplace )
322 {
323 G4double pFermi=theFermi.GetFermiMomentum(theDensity->GetDensity(aPos));
324 // protons must at least have binding energy of CoulombBarrier, so
325 // assuming the Fermi energy corresponds to a potential, we must place these such
326 // that the Fermi Energy > CoulombBarrier
327 if (theNucleons[i].GetDefinition() == G4Proton::Proton())
328 {
329 G4double nucMass = theNucleons[i].GetDefinition()->GetPDGMass();
330 G4double eFermi= std::sqrt( sqr(pFermi) + sqr(nucMass) )
331 - nucMass;
332 if (eFermi <= CoulombBarrier() ) freeplace=false;
333 }
334 }
335 if ( freeplace )
336 {
337 theNucleons[i].SetPosition(aPos);
338 places.push_back(aPos);
339 ++i;
340 }
341 }
342 }
343
344}
345
346void G4Fancy3DNucleus::ChooseFermiMomenta()
347{
348 G4int i;
349 G4double density;
350
351 // Pre-allocate buffers for filling by index
352 momentum.resize(myA, G4ThreeVector(0.,0.,0.));
353 fermiM.resize(myA, 0.*GeV);
354
355 for (G4int ntry=0; ntry<1 ; ntry ++ )
356 {
357 for (i=0; i < myA; i++ ) // momenta for all, including last, in case we swap nucleons
358 {
359 density = theDensity->GetDensity(theNucleons[i].GetPosition());
360 fermiM[i] = theFermi.GetFermiMomentum(density);
361 G4ThreeVector mom=theFermi.GetMomentum(density);
362 if (theNucleons[i].GetDefinition() == G4Proton::Proton())
363 {
364 G4double eMax = std::sqrt(sqr(fermiM[i]) +sqr(theNucleons[i].GetDefinition()->GetPDGMass()) )
365 - CoulombBarrier();
366 if ( eMax > theNucleons[i].GetDefinition()->GetPDGMass() )
367 {
368 G4double pmax2= sqr(eMax) - sqr(theNucleons[i].GetDefinition()->GetPDGMass());
369 fermiM[i] = std::sqrt(pmax2);
370 while ( mom.mag2() > pmax2 )
371 {
372 mom=theFermi.GetMomentum(density, fermiM[i]);
373 }
374 } else
375 {
376 G4cerr << "G4Fancy3DNucleus: difficulty finding proton momentum" << G4endl;
377 mom=G4ThreeVector(0,0,0);
378 }
379
380 }
381 momentum[i]= mom;
382 }
383
384 if ( ReduceSum() ) break;
385// G4cout <<" G4FancyNucleus: iterating to find momenta: "<< ntry<< G4endl;
386 }
387
388// G4ThreeVector sum;
389// for (G4int index=0; index<myA;sum+=momentum[index++])
390// ;
391// G4cout << "final sum / mag() " << sum << " / " << sum.mag() << G4endl;
392
393 G4double energy;
394 for ( i=0; i< myA ; i++ )
395 {
396 energy = theNucleons[i].GetParticleType()->GetPDGMass()
397 - BindingEnergy()/myA;
398 G4LorentzVector tempV(momentum[i],energy);
399 theNucleons[i].SetMomentum(tempV);
400 // GF 11-05-2011: set BindingEnergy to be T of Nucleon with p , ~ p**2/2m
401 //theNucleons[i].SetBindingEnergy(
402 // 0.5*sqr(fermiM[i])/theNucleons[i].GetParticleType()->GetPDGMass());
403 }
404}
405
406
407G4bool G4Fancy3DNucleus::ReduceSum()
408{
409 G4ThreeVector sum;
410 G4double PFermi=fermiM[myA-1];
411
412 for (G4int i=0; i < myA-1 ; i++ )
413 { sum+=momentum[i]; }
414
415// check if have to do anything at all..
416 if ( sum.mag() <= PFermi )
417 {
418 momentum[myA-1]=-sum;
419 return true;
420 }
421
422// find all possible changes in momentum, changing only the component parallel to sum
423 G4ThreeVector testDir=sum.unit();
424 testSums.clear();
425 testSums.resize(myA-1); // Allocate block for filling below
426
427 G4ThreeVector delta;
428 for (G4int aNucleon=0; aNucleon < myA-1; aNucleon++) {
429 delta = 2.*((momentum[aNucleon]*testDir)*testDir);
430
431 testSums[aNucleon].Fill(delta, delta.mag(), aNucleon);
432 }
433
434 std::sort(testSums.begin(), testSums.end());
435
436// reduce Momentum Sum until the next would be allowed.
437 G4int index=testSums.size();
438 while ( (sum-testSums[--index].Vector).mag()>PFermi && index>0)
439 {
440 // Only take one which improve, ie. don't change sign and overshoot...
441 if ( sum.mag() > (sum-testSums[index].Vector).mag() ) {
442 momentum[testSums[index].Index]-=testSums[index].Vector;
443 sum-=testSums[index].Vector;
444 }
445 }
446
447 if ( (sum-testSums[index].Vector).mag() <= PFermi )
448 {
449 G4int best=-1;
450 G4double pBest=2*PFermi; // anything larger than PFermi
451 for ( G4int aNucleon=0; aNucleon<=index; aNucleon++)
452 {
453 // find the momentum closest to choosen momentum for last Nucleon.
454 G4double pTry=(testSums[aNucleon].Vector-sum).mag();
455 if ( pTry < PFermi
456 && std::abs(momentum[myA-1].mag() - pTry ) < pBest )
457 {
458 pBest=std::abs(momentum[myA-1].mag() - pTry );
459 best=aNucleon;
460 }
461 }
462 if ( best < 0 )
463 {
464 G4String text = "G4Fancy3DNucleus.cc: Logic error in ReduceSum()";
465 throw G4HadronicException(__FILE__, __LINE__, text);
466 }
467 momentum[testSums[best].Index]-=testSums[best].Vector;
468 momentum[myA-1]=testSums[best].Vector-sum;
469
470 return true;
471
472 }
473
474 // try to compensate momentum using another Nucleon....
475 G4int swapit=-1;
476 while (swapit<myA-1)
477 {
478 if ( fermiM[++swapit] > PFermi ) break;
479 }
480 if (swapit == myA-1 ) return false;
481
482 // Now we have a nucleon with a bigger Fermi Momentum.
483 // Exchange with last nucleon.. and iterate.
484 std::swap(theNucleons[swapit], theNucleons[myA-1]);
485 std::swap(momentum[swapit], momentum[myA-1]);
486 std::swap(fermiM[swapit], fermiM[myA-1]);
487 return ReduceSum();
488}
489
491{
492 G4double coulombBarrier = (1.44/1.14) * MeV * myZ / (1.0 + std::pow(G4double(myA),1./3.));
493 return coulombBarrier;
494}
bool G4Fancy3DNucleusHelperForSortInZ(const G4Nucleon &nuc1, const G4Nucleon &nuc2)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
Hep3Vector unit() const
double mag2() const
double mag() const
void set(double x, double y, double z)
Hep3Vector vect() const
static void shootArray(const int size, double *vect)
Definition: RandFlat.cc:63
const std::vector< G4Nucleon > & GetNucleons()
G4double CoulombBarrier()
G4Nucleon * GetNextNucleon()
const G4VNuclearDensity * GetNuclearDensity() const
void DoLorentzContraction(const G4LorentzVector &theBoost)
G4double GetNuclearRadius()
void DoTranslation(const G4ThreeVector &theShift)
G4double GetOuterRadius()
void DoLorentzBoost(const G4LorentzVector &theBoost)
void Init(G4int theA, G4int theZ)
G4ThreeVector GetMomentum(G4double density, G4double maxMomentum=-1.)
G4double GetFermiMomentum(G4double density)
void Init(G4int anA, G4int aZ)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetBindingEnergy(const G4int A, const G4int Z)
virtual const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:68
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetDensity(const G4ThreeVector &aPosition) const
virtual G4double GetRelativeDensity(const G4ThreeVector &aPosition) const =0
virtual G4double GetRadius(const G4double maxRelativeDenisty) const =0
T sqr(const T &x)
Definition: templates.hh:145