Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4V3DNucleus.hh
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//
28// 20110805 M. Kelsey -- Change nucleon vector to use objects, not pointers
29
30#ifndef G4V3DNucleus_h
31#define G4V3DNucleus_h 1
32
33class G4Nucleon;
35#include "G4DynamicParticle.hh"
36#include "Randomize.hh"
37#include <utility>
38#include <vector>
39
41{
42
43 public:
45 virtual ~G4V3DNucleus();
46
47 private:
48 G4V3DNucleus(const G4V3DNucleus &right);
49 const G4V3DNucleus & operator=(const G4V3DNucleus &right);
50 G4bool operator==(const G4V3DNucleus &right) const;
51 G4bool operator!=(const G4V3DNucleus &right) const;
52
53 public:
54 virtual void Init(G4int theA, G4int theZ) = 0;
55 virtual G4bool StartLoop() = 0;
56 virtual G4Nucleon * GetNextNucleon() = 0;
57 virtual const std::vector<G4Nucleon> & GetNucleons() = 0;
58 virtual G4int GetMassNumber() = 0;
59 virtual G4double GetMass() = 0;
60 virtual G4int GetCharge() = 0;
62 virtual G4double GetNuclearRadius(const G4double maxRelativeDensity) = 0;
63 virtual G4double GetOuterRadius() = 0;
64 virtual G4double CoulombBarrier() = 0;
65 virtual void DoLorentzBoost(const G4LorentzVector & theBoost) = 0;
66 virtual void DoLorentzBoost(const G4ThreeVector & theBeta) = 0;
67 virtual void DoLorentzContraction(const G4LorentzVector & theBoost) = 0;
68 virtual void DoLorentzContraction(const G4ThreeVector & theBeta) = 0;
69 virtual void DoTranslation(const G4ThreeVector & theShift) = 0;
70 virtual const G4VNuclearDensity * GetNuclearDensity() const = 0;
71 virtual void SortNucleonsIncZ() = 0;
72 virtual void SortNucleonsDecZ() = 0;
73
74 public:
75 std::pair<G4double, G4double> ChooseImpactXandY(G4double maxImpact);
76 std::pair<G4double, G4double> RefetchImpactXandY(){return theImpactParameter;}
77
78 private:
79
80 std::pair<G4double, G4double> theImpactParameter;
81
82};
83
84inline
85std::pair<G4double, G4double> G4V3DNucleus::
87{
88 G4double x,y;
89 do
90 {
91 x = 2*G4UniformRand() - 1;
92 y = 2*G4UniformRand() - 1;
93 }
94 while(x*x + y*y > 1); /* Loop checking, 30-Oct-2015, G.Folger */
95
96 G4double impactX = x*(maxImpact);
97 G4double impactY = y*(maxImpact);
98 theImpactParameter.first = impactX;
99 theImpactParameter.second = impactY;
100 return theImpactParameter;
101}
102
103
104#endif
105
106
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
virtual void SortNucleonsIncZ()=0
virtual ~G4V3DNucleus()
Definition: G4V3DNucleus.cc:40
virtual void DoLorentzContraction(const G4ThreeVector &theBeta)=0
virtual G4double CoulombBarrier()=0
std::pair< G4double, G4double > RefetchImpactXandY()
Definition: G4V3DNucleus.hh:76
virtual G4double GetOuterRadius()=0
virtual const G4VNuclearDensity * GetNuclearDensity() const =0
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4bool StartLoop()=0
virtual void DoLorentzBoost(const G4ThreeVector &theBeta)=0
virtual void Init(G4int theA, G4int theZ)=0
virtual G4double GetMass()=0
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0
virtual void DoTranslation(const G4ThreeVector &theShift)=0
virtual G4double GetNuclearRadius()=0
virtual void DoLorentzContraction(const G4LorentzVector &theBoost)=0
virtual const std::vector< G4Nucleon > & GetNucleons()=0
virtual G4double GetNuclearRadius(const G4double maxRelativeDensity)=0
virtual void SortNucleonsDecZ()=0
std::pair< G4double, G4double > ChooseImpactXandY(G4double maxImpact)
Definition: G4V3DNucleus.hh:86
virtual G4int GetMassNumber()=0