Geant4 11.1.1
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"
46#include "G4Nucleon.hh"
47#include "G4SystemOfUnits.hh"
48#include "Randomize.hh"
49#include "G4ios.hh"
50#include "G4Pow.hh"
52
53#include "Randomize.hh"
54#include "G4ThreeVector.hh"
55#include "G4RandomDirection.hh"
56#include "G4LorentzRotation.hh"
57#include "G4RotationMatrix.hh"
59
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}
66
68{
69 if(theDensity) delete theDensity;
70}
71
72#if defined(NON_INTEGER_A_Z)
73void G4Fancy3DNucleus::Init(G4double theA, G4double theZ, G4int numberOfLambdas)
74{
75 G4int intZ = G4int(theZ);
76 G4int intA= ( G4UniformRand()>theA-G4int(theA) ) ? G4int(theA) : G4int(theA)+1;
77 // forward to integer Init()
78 Init(intA, intZ, std::max(numberOfLambdas, 0));
79
80}
81#endif
82
83void G4Fancy3DNucleus::Init(G4int theA, G4int theZ, G4int numberOfLambdas)
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}
130
132{
133 currentNucleon=0;
134 return (theNucleons.size()>0);
135}
136
137// Returns by pointer; null pointer indicates end of loop
139{
140 return ( (currentNucleon>=0 && currentNucleon<myA) ?
141 &theNucleons[currentNucleon++] : 0 );
142}
143
144const std::vector<G4Nucleon> & G4Fancy3DNucleus::GetNucleons()
145{
146 return theNucleons;
147}
148
149
150// Class-scope function to sort nucleons by Z coordinate
152{
153 return nuc1.GetPosition().z() < nuc2.GetPosition().z();
154}
155
157{
158 if (theNucleons.size() < 2 ) return; // Avoid unnecesary work
159
160 std::sort(theNucleons.begin(), theNucleons.end(),
162}
163
165{
166 if (theNucleons.size() < 2 ) return; // Avoid unnecessary work
168
169 std::reverse(theNucleons.begin(), theNucleons.end());
170}
171
172
173G4double G4Fancy3DNucleus::BindingEnergy()
174{
176}
177
178
180{
181 return GetNuclearRadius(0.5);
182}
183
185{
186 return theDensity->GetRadius(maxRelativeDensity);
187}
188
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}
202
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}
210
211
212
214{
215 for (G4int i=0; i<myA; i++){
216 theNucleons[i].Boost(theBoost);
217 }
218}
219
221{
222 for (G4int i=0; i<myA; i++){
223 theNucleons[i].Boost(theBeta);
224 }
225}
226
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}
240
242{
243 if (theBoost.e() !=0 ) {
244 G4ThreeVector beta = theBoost.vect()/theBoost.e();
246 }
247}
248
249
250
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}
262
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}
272
274{
275 return theDensity;
276}
277
278//----------------------- private Implementation Methods-------------
279
280void G4Fancy3DNucleus::ChooseNucleons()
281{
282 G4int protons=0, nucleons=0, lambdas=0;
283 G4double probProton = ( G4double(myZ) )/( G4double(myA) );
284 G4double probLambda = myL > 0 ? ( G4double(myL) )/( G4double(myA) ) : 0.0;
285 while ( nucleons < myA ) { /* Loop checking, 30-Oct-2015, G.Folger */
286 G4double rnd = G4UniformRand();
287 if ( rnd < probProton ) {
288 if ( protons < myZ ) {
289 protons++;
290 theNucleons[nucleons++].SetParticleType(G4Proton::Proton());
291 }
292 } else if ( rnd < probProton + probLambda ) {
293 if ( lambdas < myL ) {
294 lambdas++;
295 theNucleons[nucleons++].SetParticleType(G4Lambda::Lambda());
296 }
297 } else {
298 if ( (nucleons - protons - lambdas) < (myA - myZ - myL) ) {
299 theNucleons[nucleons++].SetParticleType(G4Neutron::Neutron());
300 }
301 }
302 }
303 return;
304}
305
306void G4Fancy3DNucleus::ChoosePositions()
307{
308 if( myA != 12) {
309
310 G4int i=0;
311 G4ThreeVector aPos, delta;
312 G4bool freeplace;
313 const G4double nd2=sqr(nucleondistance);
314 G4double maxR=GetNuclearRadius(0.001); // there are no nucleons at a
315 // relative Density of 0.01
316 G4int jr=0;
317 G4int jx,jy;
318 G4double arand[600];
319 G4double *prand=arand;
320 places.clear(); // Reset data buffer
321 G4int interationsLeft=1000*myA;
322 while ( (i < myA) && (--interationsLeft>0)) /* Loop checking, 30-Oct-2015, G.Folger */
323 {
324 do
325 {
326 if ( jr < 3 )
327 {
328 jr=std::min(600,9*(myA - i));
329 G4RandFlat::shootArray(jr,prand);
330 //CLHEP::RandFlat::shootArray(jr, prand );
331 }
332 jx=--jr;
333 jy=--jr;
334 aPos.set((2*arand[jx]-1.), (2*arand[jy]-1.), (2*arand[--jr]-1.));
335 } while (aPos.mag2() > 1. ); /* Loop checking, 30-Oct-2015, G.Folger */
336 aPos *=maxR;
337 G4double density=theDensity->GetRelativeDensity(aPos);
338 if (G4UniformRand() < density)
339 {
340 freeplace= true;
341 std::vector<G4ThreeVector>::iterator iplace;
342 for( iplace=places.begin(); iplace!=places.end() && freeplace;++iplace)
343 {
344 delta = *iplace - aPos;
345 freeplace= delta.mag2() > nd2;
346 }
347 if ( freeplace ) {
348 G4double pFermi=theFermi.GetFermiMomentum(theDensity->GetDensity(aPos));
349 // protons must at least have binding energy of CoulombBarrier, so
350 // assuming the Fermi energy corresponds to a potential, we must place these such
351 // that the Fermi Energy > CoulombBarrier
352 if (theNucleons[i].GetDefinition() == G4Proton::Proton())
353 {
354 G4double nucMass = theNucleons[i].GetDefinition()->GetPDGMass();
355 G4double eFermi= std::sqrt( sqr(pFermi) + sqr(nucMass) ) - nucMass;
356 if (eFermi <= CoulombBarrier() ) freeplace=false;
357 }
358 }
359 if ( freeplace ) {
360 theNucleons[i].SetPosition(aPos);
361 places.push_back(aPos);
362 ++i;
363 }
364 }
365 }
366 if (interationsLeft<=0) {
367 G4Exception("model/util/G4Fancy3DNucleus.cc", "mod_util001", FatalException,
368 "Problem to place nucleons");
369 }
370
371 } else {
372 // Start insertion
373 // Alpha cluster structure of carbon nuclei, C-12, is implemented according to
374 // P. Bozek, W. Broniowski, E.R. Arriola and M. Rybczynski
375 // Phys. Rev. C90, 064902 (2014)
376 const G4double Lbase=3.05*fermi;
377 const G4double Disp=0.552; // 0.91^2*2/3 fermi^2
378 const G4double nd2=sqr(nucleondistance);
379 const G4ThreeVector Corner1=G4ThreeVector( Lbase/2., 0., 0.);
380 const G4ThreeVector Corner2=G4ThreeVector(-Lbase/2., 0., 0.);
381 const G4ThreeVector Corner3=G4ThreeVector( 0.,Lbase*0.866, 0.); // 0.866=sqrt(3)/2
382 G4ThreeVector R1;
383 R1=G4ThreeVector(G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp))*fermi + Corner1;
384 theNucleons[0].SetPosition(R1); // First nucleon of the first He-4
385 G4int loopCounterLeft = 10000;
386 for(G4int ii=1; ii<4; ii++) // 2 - 4 nucleons of the first He-4
387 {
388 G4bool Continue;
389 do
390 {
391 R1=G4ThreeVector(G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp))*fermi + Corner1;
392 theNucleons[ii].SetPosition(R1);
393 Continue=false;
394 for(G4int jj=0; jj < ii; jj++)
395 {
396 if( (theNucleons[ii].GetPosition() - theNucleons[jj].GetPosition()).mag2() <= nd2 ) {Continue = true; break;}
397 }
398 } while( Continue && --loopCounterLeft > 0 ); /* Loop checking, 12-Dec-2017, A.Ribon */
399 }
400 if ( loopCounterLeft <= 0 ) {
401 G4Exception("model/util/G4Fancy3DNucleus.cc", "mod_util002", FatalException,
402 "Unable to find a good position for the first alpha cluster");
403 }
404 loopCounterLeft = 10000;
405 for(G4int ii=4; ii<8; ii++) // 5 - 8 nucleons of the second He-4
406 {
407 G4bool Continue;
408 do
409 {
410 R1=G4ThreeVector(G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp))*fermi + Corner2;
411 theNucleons[ii].SetPosition(R1);
412 Continue=false;
413 for(G4int jj=0; jj < ii; jj++)
414 {
415 if( (theNucleons[ii].GetPosition() - theNucleons[jj].GetPosition()).mag2() <= nd2 ) {Continue = true; break;}
416 }
417 } while( Continue && --loopCounterLeft > 0 ); /* Loop checking, 12-Dec-2017, A.Ribon */
418 }
419 if ( loopCounterLeft <= 0 ) {
420 G4Exception("model/util/G4Fancy3DNucleus.cc", "mod_util003", FatalException,
421 "Unable to find a good position for the second alpha cluster");
422 }
423 loopCounterLeft = 10000;
424 for(G4int ii=8; ii<12; ii++) // 9 - 12 nucleons of the third He-4
425 {
426 G4bool Continue;
427 do
428 {
429 R1=G4ThreeVector(G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp), G4RandGauss::shoot(0.,Disp))*fermi + Corner3;
430 theNucleons[ii].SetPosition(R1);
431 Continue=false;
432 for(G4int jj=0; jj < ii; jj++)
433 {
434 if( (theNucleons[ii].GetPosition() - theNucleons[jj].GetPosition()).mag2() <= nd2 ) {Continue = true; break;}
435 }
436 } while( Continue && --loopCounterLeft > 0 ); /* Loop checking, 12-Dec-2017, A.Ribon */
437 }
438 if ( loopCounterLeft <= 0 ) {
439 G4Exception("model/util/G4Fancy3DNucleus.cc", "mod_util004", FatalException,
440 "Unable to find a good position for the third alpha cluster");
441 }
442 G4LorentzRotation RandomRotation;
443 RandomRotation.rotateZ(2.*pi*G4UniformRand());
444 RandomRotation.rotateY(std::acos(2.*G4UniformRand()-1.));
445 // Randomly rotation of the created nucleus
447 for(G4int ii=0; ii<myA; ii++ )
448 {
449 Pos=G4LorentzVector(theNucleons[ii].GetPosition(),0.); Pos *=RandomRotation;
450 G4ThreeVector NewPos = Pos.vect();
451 theNucleons[ii].SetPosition(NewPos);
452 }
453
454 }
455}
456
457void G4Fancy3DNucleus::ChooseFermiMomenta()
458{
459 G4int i;
460 G4double density;
461
462 // Pre-allocate buffers for filling by index
463 momentum.resize(myA, G4ThreeVector(0.,0.,0.));
464 fermiM.resize(myA, 0.*GeV);
465
466 for (G4int ntry=0; ntry<1 ; ntry ++ )
467 {
468 for (i=0; i < myA; i++ ) // momenta for all, including last, in case we swap nucleons
469 {
470 density = theDensity->GetDensity(theNucleons[i].GetPosition());
471 fermiM[i] = theFermi.GetFermiMomentum(density);
472 G4ThreeVector mom=theFermi.GetMomentum(density);
473 if (theNucleons[i].GetDefinition() == G4Proton::Proton())
474 {
475 G4double eMax = std::sqrt(sqr(fermiM[i]) +sqr(theNucleons[i].GetDefinition()->GetPDGMass()) )
476 - CoulombBarrier();
477 if ( eMax > theNucleons[i].GetDefinition()->GetPDGMass() )
478 {
479 G4double pmax2= sqr(eMax) - sqr(theNucleons[i].GetDefinition()->GetPDGMass());
480 fermiM[i] = std::sqrt(pmax2);
481 while ( mom.mag2() > pmax2 ) /* Loop checking, 30-Oct-2015, G.Folger */
482 {
483 mom=theFermi.GetMomentum(density, fermiM[i]);
484 }
485 } else
486 {
487 //AR-21Dec2017 : emit a "JustWarning" exception instead of writing on the error stream.
488 //G4cerr << "G4Fancy3DNucleus: difficulty finding proton momentum" << G4endl;
490 ed << "Nucleus Z A " << myZ << " " << myA << G4endl;
491 ed << "proton with eMax=" << eMax << G4endl;
492 G4Exception( "G4Fancy3DNucleus::ChooseFermiMomenta(): difficulty finding proton momentum, set it to (0,0,0)",
493 "HAD_FANCY3DNUCLEUS_001", JustWarning, ed );
494 mom=G4ThreeVector(0,0,0);
495 }
496
497 }
498 momentum[i]= mom;
499 }
500
501 if ( ReduceSum() ) break;
502// G4cout <<" G4FancyNucleus: iterating to find momenta: "<< ntry<< G4endl;
503 }
504
505// G4ThreeVector sum;
506// for (G4int index=0; index<myA;sum+=momentum[index++])
507// ;
508// G4cout << "final sum / mag() " << sum << " / " << sum.mag() << G4endl;
509
511 for ( i=0; i< myA ; i++ )
512 {
513 energy = theNucleons[i].GetParticleType()->GetPDGMass()
514 - BindingEnergy()/myA;
515 G4LorentzVector tempV(momentum[i],energy);
516 theNucleons[i].SetMomentum(tempV);
517 // GF 11-05-2011: set BindingEnergy to be T of Nucleon with p , ~ p**2/2m
518 //theNucleons[i].SetBindingEnergy(
519 // 0.5*sqr(fermiM[i])/theNucleons[i].GetParticleType()->GetPDGMass());
520 }
521}
522
523
524G4bool G4Fancy3DNucleus::ReduceSum()
525{
526 G4ThreeVector sum;
527 G4double PFermi=fermiM[myA-1];
528
529 for (G4int i=0; i < myA-1 ; i++ )
530 { sum+=momentum[i]; }
531
532// check if have to do anything at all..
533 if ( sum.mag() <= PFermi )
534 {
535 momentum[myA-1]=-sum;
536 return true;
537 }
538
539// find all possible changes in momentum, changing only the component parallel to sum
540 G4ThreeVector testDir=sum.unit();
541 testSums.clear();
542 testSums.resize(myA-1); // Allocate block for filling below
543
544 G4ThreeVector delta;
545 for (G4int aNucleon=0; aNucleon < myA-1; ++aNucleon) {
546 delta = 2.*((momentum[aNucleon]*testDir)*testDir);
547
548 testSums[aNucleon].Fill(delta, delta.mag(), aNucleon);
549 }
550
551 std::sort(testSums.begin(), testSums.end());
552
553// reduce Momentum Sum until the next would be allowed.
554 G4int index=(G4int)testSums.size();
555 while ( (sum-testSums[--index].Vector).mag()>PFermi && index>0) /* Loop checking, 30-Oct-2015, G.Folger */
556 {
557 // Only take one which improve, ie. don't change sign and overshoot...
558 if ( sum.mag() > (sum-testSums[index].Vector).mag() ) {
559 momentum[testSums[index].Index]-=testSums[index].Vector;
560 sum-=testSums[index].Vector;
561 }
562 }
563
564 if ( (sum-testSums[index].Vector).mag() <= PFermi )
565 {
566 G4int best=-1;
567 G4double pBest=2*PFermi; // anything larger than PFermi
568 for ( G4int aNucleon=0; aNucleon<=index; aNucleon++)
569 {
570 // find the momentum closest to choosen momentum for last Nucleon.
571 G4double pTry=(testSums[aNucleon].Vector-sum).mag();
572 if ( pTry < PFermi
573 && std::abs(momentum[myA-1].mag() - pTry ) < pBest )
574 {
575 pBest=std::abs(momentum[myA-1].mag() - pTry );
576 best=aNucleon;
577 }
578 }
579 if ( best < 0 )
580 {
581 G4String text = "G4Fancy3DNucleus.cc: Logic error in ReduceSum()";
582 throw G4HadronicException(__FILE__, __LINE__, text);
583 }
584 momentum[testSums[best].Index]-=testSums[best].Vector;
585 momentum[myA-1]=testSums[best].Vector-sum;
586
587 return true;
588
589 }
590
591 // try to compensate momentum using another Nucleon....
592 G4int swapit=-1;
593 while (swapit<myA-1) /* Loop checking, 30-Oct-2015, G.Folger */
594 {
595 if ( fermiM[++swapit] > PFermi ) break;
596 }
597 if (swapit == myA-1 ) return false;
598
599 // Now we have a nucleon with a bigger Fermi Momentum.
600 // Exchange with last nucleon.. and iterate.
601 std::swap(theNucleons[swapit], theNucleons[myA-1]);
602 std::swap(momentum[swapit], momentum[myA-1]);
603 std::swap(fermiM[swapit], fermiM[myA-1]);
604 return ReduceSum();
605}
606
608{
609 static const G4double cfactor = (1.44/1.14) * MeV;
610 return cfactor*myZ/(1.0 + G4Pow::GetInstance()->Z13(myA));
611}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
bool G4Fancy3DNucleusHelperForSortInZ(const G4Nucleon &nuc1, const G4Nucleon &nuc2)
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
double mag2() const
double mag() const
void set(double x, double y, double z)
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
Hep3Vector vect() const
const std::vector< G4Nucleon > & GetNucleons()
G4double CoulombBarrier()
G4Nucleon * GetNextNucleon()
const G4VNuclearDensity * GetNuclearDensity() const
void DoLorentzContraction(const G4LorentzVector &theBoost)
void Init(G4int theA, G4int theZ, G4int numberOfLambdas=0)
G4double GetNuclearRadius()
void DoTranslation(const G4ThreeVector &theShift)
G4double GetOuterRadius()
void DoLorentzBoost(const G4LorentzVector &theBoost)
G4ThreeVector GetMomentum(G4double density, G4double maxMomentum=-1.)
G4double GetFermiMomentum(G4double density)
void Init(G4int anA, G4int aZ)
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4double GetBindingEnergy(const G4int A, const G4int Z)
const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:140
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4double GetDensity(const G4ThreeVector &aPosition) const
virtual G4double GetRelativeDensity(const G4ThreeVector &aPosition) const =0
virtual G4double GetRadius(const G4double maxRelativeDenisty) const =0
ush Pos
Definition: deflate.h:91
G4double energy(const ThreeVector &p, const G4double m)
T sqr(const T &x)
Definition: templates.hh:128