Geant4  10.02.p02
G4Fancy3DNucleus.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 #ifndef G4Fancy3DNucleus_h
27 #define G4Fancy3DNucleus_h 1
28 
29 // ------------------------------------------------------------
30 // GEANT 4 class header file
31 //
32 // ---------------- G4Fancy3DNucleus ----------------
33 // by Gunter Folger, May 1998.
34 // class for a 3D nucleus, arranging nucleons in space and momentum.
35 // ------------------------------------------------------------
36 // 20110805 M. Kelsey -- Remove C-style array (pointer) of G4Nucleons,
37 // make vector a container of objects. Move testSums,
38 // places, momentum and fermiM to class data members for
39 // reuse. Remove args from ReduceSum(), use data members.
40 
41 #include "globals.hh"
42 #include "G4DynamicParticle.hh"
43 #include "G4Nucleon.hh" /* FIXME: This should be forward decl! */
44 #include "G4V3DNucleus.hh"
45 #include "G4VNuclearDensity.hh"
46 #include "G4FermiMomentum.hh"
47 #include <vector>
48 
50 
51 // to test if we can drop old interface for (A,Z), comment next line..
52 //#define NON_INTEGER_A_Z 1
53 
55 {
56 
57  public:
60 
61  private:
63  const G4Fancy3DNucleus & operator=(const G4Fancy3DNucleus &right);
64  int operator==(const G4Fancy3DNucleus &right) const;
65  int operator!=(const G4Fancy3DNucleus &right) const;
66 
67 
68 // Implementation
69  void ChooseNucleons();
70  void ChoosePositions();
71  void ChooseFermiMomenta();
73  G4bool ReduceSum();
74 
75  public:
76 #if defined(NON_INTEGER_A_Z)
77  void Init(G4double theA, G4double theZ);
78 #endif
79  void Init(G4int theA, G4int theZ);
80  G4bool StartLoop();
82  const std::vector<G4Nucleon> & GetNucleons();
84  G4double GetMass();
85  G4int GetCharge();
87  G4double GetNuclearRadius(const G4double maxRelativeDensity);
92  void DoLorentzBoost(const G4LorentzVector & theBoost);
93  void DoLorentzBoost(const G4ThreeVector & theBeta);
94  void DoLorentzContraction(const G4LorentzVector & theBoost);
95  void DoLorentzContraction(const G4ThreeVector & theBeta);
96  void CenterNucleons();
97  void DoTranslation(const G4ThreeVector & theShift);
98  const G4VNuclearDensity * GetNuclearDensity() const;
99  void SortNucleonsIncZ(); // on increased Z-coordinates Uzhi 29.08.08
100  void SortNucleonsDecZ(); // on decreased Z-coordinates Uzhi 29.08.08
101 
102  private:
103 
106  std::vector<G4Nucleon> theNucleons;
107 
113 
114  std::vector<G4ThreeVector> places; // For selecting locations
115  std::vector<G4ThreeVector> momentum; // For selecting nucleon motion
116  std::vector<G4double> fermiM;
117  std::vector<G4Fancy3DNucleusHelper> testSums; // For sorting nucleon configs
118 };
119 
120 
122 {
123  return myZ;
124 }
125 
127 {
128  return myA;
129 }
131 {
132  excitationEnergy +=anE;
133  return excitationEnergy;
134 }
135 
137 {
138  return excitationEnergy;
139 }
140 
141 #endif
std::vector< G4ThreeVector > momentum
CLHEP::Hep3Vector G4ThreeVector
G4double BindingEnergy()
G4double CoulombBarrier()
void DoTranslation(const G4ThreeVector &theShift)
G4double GetExcitationEnergy()
std::vector< G4double > fermiM
void DoLorentzBoost(const G4LorentzVector &theBoost)
std::vector< G4Fancy3DNucleusHelper > testSums
int G4int
Definition: G4Types.hh:78
int operator==(const G4Fancy3DNucleus &right) const
std::vector< G4ThreeVector > places
int operator!=(const G4Fancy3DNucleus &right) const
G4double AddExcitationEnergy(G4double)
std::vector< G4Nucleon > theNucleons
bool G4bool
Definition: G4Types.hh:79
void Init(G4int theA, G4int theZ)
G4double GetNuclearRadius()
const std::vector< G4Nucleon > & GetNucleons()
const G4Fancy3DNucleus & operator=(const G4Fancy3DNucleus &right)
void DoLorentzContraction(const G4LorentzVector &theBoost)
G4FermiMomentum theFermi
G4double GetOuterRadius()
const G4VNuclearDensity * GetNuclearDensity() const
double G4double
Definition: G4Types.hh:76
const G4double nucleondistance
G4Nucleon * GetNextNucleon()
G4VNuclearDensity * theDensity
CLHEP::HepLorentzVector G4LorentzVector