File indexing completed on 2025-01-18 09:59:18
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030 #ifndef G4V3DNucleus_h
0031 #define G4V3DNucleus_h 1
0032
0033 class G4Nucleon;
0034 class G4VNuclearDensity;
0035 #include "G4DynamicParticle.hh"
0036 #include "Randomize.hh"
0037 #include <utility>
0038 #include <vector>
0039
0040 class G4V3DNucleus
0041 {
0042
0043 public:
0044 G4V3DNucleus();
0045 virtual ~G4V3DNucleus();
0046
0047 private:
0048 G4V3DNucleus(const G4V3DNucleus &right);
0049 const G4V3DNucleus & operator=(const G4V3DNucleus &right);
0050 G4bool operator==(const G4V3DNucleus &right) const;
0051 G4bool operator!=(const G4V3DNucleus &right) const;
0052
0053 public:
0054 virtual void Init(G4int theA, G4int theZ, G4int numberOfLambdas = 0) = 0;
0055 virtual G4bool StartLoop() = 0;
0056 virtual G4Nucleon * GetNextNucleon() = 0;
0057 virtual const std::vector<G4Nucleon> & GetNucleons() = 0;
0058 virtual G4int GetMassNumber() = 0;
0059 virtual G4double GetMass() = 0;
0060 virtual G4int GetCharge() = 0;
0061 virtual G4int GetNumberOfLambdas() = 0;
0062 virtual G4double GetNuclearRadius() = 0;
0063 virtual G4double GetNuclearRadius(const G4double maxRelativeDensity) = 0;
0064 virtual G4double GetOuterRadius() = 0;
0065 virtual G4double CoulombBarrier() = 0;
0066 virtual void DoLorentzBoost(const G4LorentzVector & theBoost) = 0;
0067 virtual void DoLorentzBoost(const G4ThreeVector & theBeta) = 0;
0068 virtual void DoLorentzContraction(const G4LorentzVector & theBoost) = 0;
0069 virtual void DoLorentzContraction(const G4ThreeVector & theBeta) = 0;
0070 virtual void DoTranslation(const G4ThreeVector & theShift) = 0;
0071 virtual const G4VNuclearDensity * GetNuclearDensity() const = 0;
0072 virtual void SortNucleonsIncZ() = 0;
0073 virtual void SortNucleonsDecZ() = 0;
0074
0075 public:
0076 std::pair<G4double, G4double> ChooseImpactXandY(G4double maxImpact);
0077 std::pair<G4double, G4double> RefetchImpactXandY(){return theImpactParameter;}
0078
0079 private:
0080
0081 std::pair<G4double, G4double> theImpactParameter;
0082
0083 };
0084
0085 inline
0086 std::pair<G4double, G4double> G4V3DNucleus::
0087 ChooseImpactXandY(G4double maxImpact)
0088 {
0089 G4double x,y;
0090 do
0091 {
0092 x = 2*G4UniformRand() - 1;
0093 y = 2*G4UniformRand() - 1;
0094 }
0095 while(x*x + y*y > 1);
0096
0097 G4double impactX = x*(maxImpact);
0098 G4double impactY = y*(maxImpact);
0099 theImpactParameter.first = impactX;
0100 theImpactParameter.second = impactY;
0101 return theImpactParameter;
0102 }
0103
0104
0105 #endif
0106
0107