source: trunk/source/processes/hadronic/models/util/src/G4Fancy3DNucleus.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 15.6 KB
Line 
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//
28// ------------------------------------------------------------
29// GEANT 4 class implementation file
30//
31// ---------------- G4Fancy3DNucleus ----------------
32// by Gunter Folger, May 1998.
33// class for a 3D nucleus, arranging nucleons in space and momentum.
34// ------------------------------------------------------------
35
36#include "G4Fancy3DNucleus.hh"
37#include "G4NuclearFermiDensity.hh"
38#include "G4NuclearShellModelDensity.hh"
39#include "G4NucleiProperties.hh"
40#include "Randomize.hh"
41#include "G4ios.hh"
42#include <algorithm>
43#include "G4HadronicException.hh"
44
45
46G4Fancy3DNucleus::G4Fancy3DNucleus()
47 : nucleondistance(0.8*fermi)
48{
49 theDensity=0;
50 theNucleons=0;
51 currentNucleon=-1;
52 myA=0;
53 myZ=0;
54//G4cout <<"G4Fancy3DNucleus::G4Fancy3DNucleus()"<<G4endl;
55}
56
57G4Fancy3DNucleus::~G4Fancy3DNucleus()
58{
59 if(theNucleons) delete [] theNucleons;
60 if(theDensity) delete theDensity;
61}
62
63
64void G4Fancy3DNucleus::Init(G4double theA, G4double theZ)
65{
66// G4cout << "G4Fancy3DNucleus::Init(theA, theZ) called"<<G4endl;
67 currentNucleon=-1;
68 if(theNucleons) delete [] theNucleons;
69
70 theRWNucleons.clear();
71
72 myZ = G4int(theZ);
73 myA= ( G4UniformRand()>theA-G4int(theA) ) ? G4int(theA) : G4int(theA)+1;
74
75 theNucleons = new G4Nucleon[myA];
76
77// G4cout << "myA, myZ" << myA << ", " << myZ << G4endl;
78
79 if(theDensity) delete theDensity;
80 if ( myA < 17 ) {
81 theDensity = new G4NuclearShellModelDensity(myA, myZ);
82 } else {
83 theDensity = new G4NuclearFermiDensity(myA, myZ);
84 }
85
86 theFermi.Init(myA, myZ);
87
88 ChooseNucleons();
89
90 ChoosePositions();
91
92// CenterNucleons(); // This would introduce a bias
93
94 ChooseFermiMomenta();
95
96 G4double Ebinding= BindingEnergy()/myA;
97
98 for (G4int aNucleon=0; aNucleon < myA; aNucleon++)
99 {
100 theNucleons[aNucleon].SetBindingEnergy(Ebinding);
101 }
102
103
104 return;
105}
106
107G4bool G4Fancy3DNucleus::StartLoop()
108{
109 currentNucleon=0;
110 return theNucleons;
111}
112
113G4Nucleon * G4Fancy3DNucleus::GetNextNucleon()
114{
115 return ( currentNucleon>=0 && currentNucleon<myA ) ?
116 theNucleons+currentNucleon++ : 0;
117}
118
119const std::vector<G4Nucleon *> & G4Fancy3DNucleus::GetNucleons()
120{
121 if ( theRWNucleons.size()==0 )
122 {
123 for (G4int i=0; i< myA; i++)
124 {
125 theRWNucleons.push_back(theNucleons+i);
126 }
127 }
128 return theRWNucleons;
129}
130
131//void G4Fancy3DNucleus::SortNucleonsIncZ() // on increased Z-coordinates Uzhi 29.08.08
132
133 bool G4Fancy3DNucleusHelperForSortInZ(const G4Nucleon* nuc1, const G4Nucleon* nuc2)
134{
135 return nuc1->GetPosition().z() < nuc2->GetPosition().z();
136}
137
138//void G4Fancy3DNucleus::SortNucleonsInZ()
139void G4Fancy3DNucleus::SortNucleonsIncZ() // on increased Z-coordinates Uzhi 29.08.08
140{
141
142 GetNucleons(); // make sure theRWNucleons is initialised
143
144 if (theRWNucleons.size() < 2 ) return;
145
146 sort( theRWNucleons.begin(),theRWNucleons.end(),G4Fancy3DNucleusHelperForSortInZ);
147
148// now copy sorted nucleons to theNucleons array. TheRWNucleons are pointers in theNucleons
149// so we need to copy to new, and then swap.
150 G4Nucleon * sortedNucleons = new G4Nucleon[myA];
151 for ( unsigned int i=0; i<theRWNucleons.size(); i++ )
152 {
153 sortedNucleons[i]= *(theRWNucleons[i]);
154 }
155
156 theRWNucleons.clear(); // about to delete array these point to....
157 delete [] theNucleons;
158
159 theNucleons=sortedNucleons;
160
161 return;
162}
163
164void G4Fancy3DNucleus::SortNucleonsDecZ() // on decreased Z-coordinates Uzhi 29.08.08
165{
166 G4Nucleon * sortedNucleons = new G4Nucleon[myA];
167
168 GetNucleons(); // make sure theRWNucleons is initialised
169
170 if (theRWNucleons.size() < 2 ) return;
171 sort( theRWNucleons.begin(),theRWNucleons.end(),G4Fancy3DNucleusHelperForSortInZ);
172
173// now copy sorted nucleons to theNucleons array. TheRWNucleons are pointers in theNucleons
174// so we need to copy to new, and then swap.
175 for ( unsigned int i=0; i<theRWNucleons.size(); i++ )
176 {
177 sortedNucleons[i]= *(theRWNucleons[myA-1-i]); // Uzhi 29.08.08
178 }
179 theRWNucleons.clear(); // about to delete array elements these point to....
180 delete [] theNucleons;
181 theNucleons=sortedNucleons;
182
183 return;
184}
185
186G4double G4Fancy3DNucleus::BindingEnergy()
187{
188 return G4NucleiProperties::GetBindingEnergy(myA,myZ);
189}
190
191
192G4double G4Fancy3DNucleus::GetNuclearRadius()
193{
194 return GetNuclearRadius(0.5);
195}
196
197G4double G4Fancy3DNucleus::GetNuclearRadius(const G4double maxRelativeDensity)
198{
199 return theDensity->GetRadius(maxRelativeDensity);
200}
201
202G4double G4Fancy3DNucleus::GetOuterRadius()
203{
204 G4double maxradius2=0;
205
206 for (int i=0; i<myA; i++)
207 {
208 if ( theNucleons[i].GetPosition().mag2() > maxradius2 )
209 {
210 maxradius2=theNucleons[i].GetPosition().mag2();
211 }
212 }
213 return std::sqrt(maxradius2)+nucleondistance;
214}
215
216G4double G4Fancy3DNucleus::GetMass()
217{
218 return myZ*G4Proton::Proton()->GetPDGMass() +
219 (myA-myZ)*G4Neutron::Neutron()->GetPDGMass() -
220 BindingEnergy();
221}
222
223
224
225void G4Fancy3DNucleus::DoLorentzBoost(const G4LorentzVector & theBoost)
226{
227 for (G4int i=0; i<myA; i++){
228 theNucleons[i].Boost(theBoost);
229 }
230}
231
232void G4Fancy3DNucleus::DoLorentzBoost(const G4ThreeVector & theBeta)
233{
234 for (G4int i=0; i<myA; i++){
235 theNucleons[i].Boost(theBeta);
236 }
237}
238
239void G4Fancy3DNucleus::DoLorentzContraction(const G4ThreeVector & theBeta)
240{
241 G4double factor=(1-std::sqrt(1-theBeta.mag2()))/theBeta.mag2(); // (gamma-1)/gamma/beta**2
242 for (G4int i=0; i< myA; i++)
243 {
244 G4ThreeVector rprime=theNucleons[i].GetPosition() -
245 factor * (theBeta*theNucleons[i].GetPosition()) *
246 // theNucleons[i].GetPosition();
247 theBeta;
248 theNucleons[i].SetPosition(rprime);
249 }
250}
251
252void G4Fancy3DNucleus::DoLorentzContraction(const G4LorentzVector & theBoost)
253{
254 G4ThreeVector beta= 1/theBoost.e() * theBoost.vect();
255 // DoLorentzBoost(beta);
256 DoLorentzContraction(beta);
257}
258
259
260
261void G4Fancy3DNucleus::CenterNucleons()
262{
263 G4ThreeVector center;
264
265 for (G4int i=0; i<myA; i++ )
266 {
267 center+=theNucleons[i].GetPosition();
268 }
269 center *= -1./myA;
270 DoTranslation(center);
271}
272
273void G4Fancy3DNucleus::DoTranslation(const G4ThreeVector & theShift)
274{
275 for (G4int i=0; i<myA; i++ )
276 {
277 G4ThreeVector tempV = theNucleons[i].GetPosition() + theShift;
278 theNucleons[i].SetPosition(tempV);
279 }
280}
281
282const G4VNuclearDensity * G4Fancy3DNucleus::GetNuclearDensity() const
283{
284 return theDensity;
285}
286
287//----------------------- private Implementation Methods-------------
288
289void G4Fancy3DNucleus::ChooseNucleons()
290{
291 G4int protons=0,nucleons=0;
292
293 while (nucleons < myA )
294 {
295 if ( protons < myZ && G4UniformRand() < (G4double)(myZ-protons)/(G4double)(myA-nucleons) )
296 {
297 protons++;
298 theNucleons[nucleons++].SetParticleType(G4Proton::Proton());
299 }
300 else if ( (nucleons-protons) < (myA-myZ) )
301 {
302 theNucleons[nucleons++].SetParticleType(G4Neutron::Neutron());
303 }
304 else G4cout << "G4Fancy3DNucleus::ChooseNucleons not efficient" << G4endl;
305 }
306 return;
307}
308
309void G4Fancy3DNucleus::ChoosePositions()
310{
311 G4int i=0;
312 G4ThreeVector aPos, delta;
313 std::vector<G4ThreeVector> places;
314 places.reserve(myA);
315 G4bool freeplace;
316 static G4double nd2 = sqr(nucleondistance);
317 G4double maxR=GetNuclearRadius(0.001); // there are no nucleons at a
318 // relative Density of 0.01
319 G4int jr=0;
320 G4int jx,jy;
321 G4double arand[600];
322 G4double *prand=arand;
323 while ( i < myA )
324 {
325 do
326 {
327 if ( jr < 3 )
328 {
329 jr=std::min(600,9*(myA - i));
330 CLHEP::RandFlat::shootArray(jr, prand );
331 }
332 jx=--jr;
333 jy=--jr;
334 aPos=G4ThreeVector( (2*arand[jx]-1.),
335 (2*arand[jy]-1.),
336 (2*arand[--jr]-1.));
337 } while (aPos.mag2() > 1. );
338 aPos *=maxR;
339 G4double density=theDensity->GetRelativeDensity(aPos);
340 if (G4UniformRand() < density)
341 {
342 freeplace= true;
343 std::vector<G4ThreeVector>::iterator iplace;
344 for( iplace=places.begin(); iplace!=places.end() && freeplace;++iplace)
345 {
346 delta = *iplace - aPos;
347 freeplace= delta.mag2() > nd2;
348 }
349
350 if ( freeplace )
351 {
352 G4double pFermi=theFermi.GetFermiMomentum(theDensity->GetDensity(aPos));
353 // protons must at least have binding energy of CoulombBarrier, so
354 // assuming the Fermi energy corresponds to a potential, we must place these such
355 // that the Fermi Energy > CoulombBarrier
356 if (theNucleons[i].GetDefinition() == G4Proton::Proton())
357 {
358 G4double eFermi= std::sqrt( sqr(pFermi) + sqr(theNucleons[i].GetDefinition()->GetPDGMass()) )
359 - theNucleons[i].GetDefinition()->GetPDGMass();
360 if (eFermi <= CoulombBarrier() ) freeplace=false;
361 }
362 }
363 if ( freeplace )
364 {
365 theNucleons[i].SetPosition(aPos);
366 places.push_back(aPos);
367 ++i;
368 }
369 }
370 }
371
372}
373
374void G4Fancy3DNucleus::ChooseFermiMomenta()
375{
376 G4int i;
377 G4double density;
378 G4ThreeVector * momentum=new G4ThreeVector[myA];
379
380 G4double * fermiM=new G4double[myA];
381
382 for (G4int ntry=0; ntry<1 ; ntry ++ )
383 {
384 for (i=0; i < myA; i++ ) // momenta for all, including last, in case we swap nucleons
385 {
386 density = theDensity->GetDensity(theNucleons[i].GetPosition());
387 fermiM[i] = theFermi.GetFermiMomentum(density);
388 G4ThreeVector mom=theFermi.GetMomentum(density);
389 if (theNucleons[i].GetDefinition() == G4Proton::Proton())
390 {
391 G4double eMax = std::sqrt(sqr(fermiM[i]) +sqr(theNucleons[i].GetDefinition()->GetPDGMass()) )
392 - CoulombBarrier();
393 if ( eMax > theNucleons[i].GetDefinition()->GetPDGMass() )
394 {
395 G4double pmax2= sqr(eMax) - sqr(theNucleons[i].GetDefinition()->GetPDGMass());
396 fermiM[i] = std::sqrt(pmax2);
397 while ( mom.mag2() > pmax2 )
398 {
399 mom=theFermi.GetMomentum(density, fermiM[i]);
400 }
401 } else
402 {
403 G4cerr << "G4Fancy3DNucleus: difficulty finding proton momentum" << G4endl;
404 mom=G4ThreeVector(0,0,0);
405 }
406
407 }
408 momentum[i]= mom;
409 }
410
411 if (ReduceSum(momentum,fermiM) )
412 break;
413// G4cout <<" G4FancyNucleus: iterating to find momenta: "<< ntry<< G4endl;
414 }
415
416// G4ThreeVector sum;
417// for (G4int index=0; index<myA;sum+=momentum[index++])
418// ;
419// G4cout << "final sum / mag() " << sum << " / " << sum.mag() << G4endl;
420
421 G4double energy;
422 for ( i=0; i< myA ; i++ )
423 {
424 energy = theNucleons[i].GetParticleType()->GetPDGMass()
425 - BindingEnergy()/myA;
426 G4LorentzVector tempV(momentum[i],energy);
427 theNucleons[i].SetMomentum(tempV);
428 }
429
430 delete [] momentum;
431 delete [] fermiM;
432}
433
434 class G4Fancy3DNucleusHelper // Helper class
435 {
436 public:
437 G4Fancy3DNucleusHelper(const G4ThreeVector &vec,const G4double size,const G4int index)
438 : Vector(vec), Size(size), anInt(index) {}
439 int operator ==(const G4Fancy3DNucleusHelper &right) const
440 {
441 return this==&right;
442 }
443 int operator < (const G4Fancy3DNucleusHelper &right) const
444 {
445 return size()<right.size();
446 }
447 const G4ThreeVector& vector() const
448 {
449 return Vector;
450 }
451 G4double size() const
452 {
453 return Size;
454 }
455 G4int index() const
456 {
457 return anInt;
458 }
459 G4Fancy3DNucleusHelper operator =(const G4Fancy3DNucleusHelper &right)
460 {
461 Vector = right.Vector;
462 Size = right.Size;
463 anInt = right.anInt;
464 return *this;
465 }
466
467 private:
468 G4Fancy3DNucleusHelper(): Vector(0), Size(0), anInt(0) {G4cout << "def ctor for MixMasch" << G4endl;}
469 G4ThreeVector Vector;
470 G4double Size;
471 G4int anInt;
472 };
473
474G4bool G4Fancy3DNucleus::ReduceSum(G4ThreeVector * momentum, G4double *pFermiM)
475{
476 G4ThreeVector sum;
477 G4double PFermi=pFermiM[myA-1];
478
479 for (G4int i=0; i < myA-1 ; i++ )
480 { sum+=momentum[i]; }
481
482// check if have to do anything at all..
483 if ( sum.mag() <= PFermi )
484 {
485 momentum[myA-1]=-sum;
486 return true;
487 }
488
489// find all possible changes in momentum, changing only the component parallel to sum
490 G4ThreeVector testDir=sum.unit();
491 std::vector<G4Fancy3DNucleusHelper> testSums; // Sorted on delta.mag()
492
493 for ( G4int aNucleon=0; aNucleon < myA-1; aNucleon++){
494 G4ThreeVector delta=2*((momentum[aNucleon]*testDir)* testDir);
495 testSums.push_back(G4Fancy3DNucleusHelper(delta,delta.mag(),aNucleon));
496 }
497 std::sort(testSums.begin(), testSums.end());
498
499// reduce Momentum Sum until the next would be allowed.
500 G4int index=testSums.size();
501 while ( (sum-testSums[--index].vector()).mag()>PFermi && index>0)
502 {
503 // Only take one which improve, ie. don't change sign and overshoot...
504 if ( sum.mag() > (sum-testSums[index].vector()).mag() ) {
505 momentum[testSums[index].index()]-=testSums[index].vector();
506 sum-=testSums[index].vector();
507 }
508 }
509
510 if ( (sum-testSums[index].vector()).mag() <= PFermi )
511 {
512 G4int best=-1;
513 G4double pBest=2*PFermi; // anything larger than PFermi
514 for ( G4int aNucleon=0; aNucleon<=index; aNucleon++)
515 {
516 // find the momentum closest to choosen momentum for last Nucleon.
517 G4double pTry=(testSums[aNucleon].vector()-sum).mag();
518 if ( pTry < PFermi
519 && std::abs(momentum[myA-1].mag() - pTry ) < pBest )
520 {
521 pBest=std::abs(momentum[myA-1].mag() - pTry );
522 best=aNucleon;
523 }
524 }
525 if ( best < 0 )
526 {
527 G4String text = "G4Fancy3DNucleus.cc: Logic error in ReduceSum()";
528 throw G4HadronicException(__FILE__, __LINE__, text);
529 }
530 momentum[testSums[best].index()]-=testSums[best].vector();
531 momentum[myA-1]=testSums[best].vector()-sum;
532
533 testSums.clear();
534 return true;
535
536 }
537 testSums.clear();
538
539 // try to compensate momentum using another Nucleon....
540
541 G4int swapit=-1;
542 while (swapit<myA-1)
543 {
544 if ( pFermiM[++swapit] > PFermi ) break;
545 }
546 if (swapit == myA-1 ) return false;
547
548 // Now we have a nucleon with a bigger Fermi Momentum.
549 // Exchange with last nucleon.. and iterate.
550 G4Nucleon swap= theNucleons[swapit];
551 G4ThreeVector mom_swap=momentum[swapit];
552 G4double pf=pFermiM[swapit];
553 theNucleons[swapit]=theNucleons[myA-1];
554 momentum[swapit]=momentum[myA-1];
555 pFermiM[swapit]=pFermiM[myA-1];
556 theNucleons[myA-1]=swap;
557 momentum[myA-1]=mom_swap;
558 pFermiM[myA-1]=pf;
559 return ReduceSum(momentum,pFermiM);
560}
561
562G4double G4Fancy3DNucleus::CoulombBarrier()
563{
564 G4double coulombBarrier = (1.44/1.14) * MeV * myZ / (1.0 + std::pow(G4double(myA),1./3.));
565 return coulombBarrier;
566}
Note: See TracBrowser for help on using the repository browser.