Geant4  10.01.p01
G4Nucleon.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 #ifndef G4Nucleon_h
29 #define G4Nucleon_h 1
30 
31 // ------------------------------------------------------------
32 // GEANT 4 class header file
33 //
34 // ---------------- G4Nucleon ----------------
35 // by Gunter Folger, May 1998.
36 // class for a nucleon (inside a 3D Nucleus)
37 // ------------------------------------------------------------
38 
39 #include "G4ThreeVector.hh"
40 #include "G4LorentzVector.hh"
41 #include "globals.hh"
42 #include "G4ParticleDefinition.hh"
43 #include "G4Proton.hh"
44 #include "G4Neutron.hh"
45 
46 #include "G4AntiProton.hh" // Uzhi Feb. 2011
47 #include "G4AntiNeutron.hh" // Uzhi Feb. 2011
48 
49 #include "G4VKineticNucleon.hh"
50 
51 //#include "G4VSplitableHadron.hh"
52 class G4VSplitableHadron;
53 
55 {
56 
57  public:
58  G4Nucleon();
59  ~G4Nucleon();
60 
61  inline int operator==(const G4Nucleon &right) const;
62  inline int operator!=(const G4Nucleon &right) const;
63  G4Nucleon& operator=(const G4Nucleon& right);
64 
65  public:
66 
67  inline void SetPosition(G4ThreeVector & aPosition) {thePosition = aPosition;}
68  virtual inline const G4ThreeVector & GetPosition() const {return thePosition;}
69 
70  inline void SetMomentum(G4LorentzVector & aMomentum) {theMomentum = aMomentum;}
71  inline const G4LorentzVector& GetMomentum() const {return theMomentum;}
72  virtual inline const G4LorentzVector & Get4Momentum() const {return theMomentum;}
73 
74  inline void SetBindingEnergy(G4double anEnergy) {theBindingE = anEnergy;}
75  inline G4double GetBindingEnergy() const {return theBindingE;}
76 
77  inline void SetParticleType(G4Proton * aProton) {theParticleType = aProton;}
78  inline void SetParticleType(G4Neutron *aNeutron){theParticleType = aNeutron;}
79 
80  inline void SetParticleType(G4AntiProton * aAntiProton) {theParticleType =aAntiProton;} //VU
81  inline void SetParticleType(G4AntiNeutron *aAntiNeutron){theParticleType =aAntiNeutron;}//VU
82 
83 
84  inline const G4ParticleDefinition* GetParticleType() const {return theParticleType;}
85  virtual const G4ParticleDefinition* GetDefinition() const {return theParticleType;}
86 
87  inline void Boost(const G4ThreeVector & beta){ theMomentum.boost(beta); }
88  void Boost(const G4LorentzVector & aMomentum);
89 
90  inline void Hit(G4VSplitableHadron *aHit) { theSplitableHadron=aHit;}
91 // inline void Hit(G4int ) { isHit=true;}
92  inline void Hit(G4int )
93  {
94  theSplitableHadron=reinterpret_cast<G4VSplitableHadron *>(1111);
95  }
97  inline G4bool AreYouHit() const { return theSplitableHadron!=0;}
98 
99  private:
100 
106 
107 
108 };
109 
110 std::ostream & operator << (std::ostream &, const G4Nucleon&);
111 
112 inline int G4Nucleon::operator==(const G4Nucleon &right) const
113 {
114  return this==&right;
115 }
116 inline int G4Nucleon::operator!=(const G4Nucleon &right) const
117 {
118  return this!=&right;
119 }
120 
122 {
123  if (this != &right)
124  {
125  thePosition=right.GetPosition();
126  theMomentum=right.Get4Momentum();
130  }
131  return *this;
132 }
133 
134 #endif
135 
136 
void SetParticleType(G4AntiNeutron *aAntiNeutron)
Definition: G4Nucleon.hh:81
const G4ParticleDefinition * GetParticleType() const
Definition: G4Nucleon.hh:84
void SetPosition(G4ThreeVector &aPosition)
Definition: G4Nucleon.hh:67
void SetParticleType(G4Proton *aProton)
Definition: G4Nucleon.hh:77
const G4LorentzVector & GetMomentum() const
Definition: G4Nucleon.hh:71
void SetMomentum(G4LorentzVector &aMomentum)
Definition: G4Nucleon.hh:70
CLHEP::Hep3Vector G4ThreeVector
const G4ParticleDefinition * theParticleType
Definition: G4Nucleon.hh:104
G4LorentzVector theMomentum
Definition: G4Nucleon.hh:102
int operator!=(const G4Nucleon &right) const
Definition: G4Nucleon.hh:116
void SetParticleType(G4AntiProton *aAntiProton)
Definition: G4Nucleon.hh:80
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
G4ThreeVector thePosition
Definition: G4Nucleon.hh:101
virtual const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:68
int G4int
Definition: G4Types.hh:78
virtual const G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
int operator==(const G4Nucleon &right) const
Definition: G4Nucleon.hh:112
G4VSplitableHadron * GetSplitableHadron() const
Definition: G4Nucleon.hh:96
G4double theBindingE
Definition: G4Nucleon.hh:103
bool G4bool
Definition: G4Types.hh:79
G4bool AreYouHit() const
Definition: G4Nucleon.hh:97
G4Nucleon & operator=(const G4Nucleon &right)
Definition: G4Nucleon.hh:121
void SetBindingEnergy(G4double anEnergy)
Definition: G4Nucleon.hh:74
std::ostream & operator<<(std::ostream &, const G4Nucleon &)
Definition: G4Nucleon.cc:58
G4VSplitableHadron * theSplitableHadron
Definition: G4Nucleon.hh:105
void Hit(G4int)
Definition: G4Nucleon.hh:92
void Boost(const G4ThreeVector &beta)
Definition: G4Nucleon.hh:87
void SetParticleType(G4Neutron *aNeutron)
Definition: G4Nucleon.hh:78
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:90
double G4double
Definition: G4Types.hh:76
~G4Nucleon()
Definition: G4Nucleon.cc:42
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
CLHEP::HepLorentzVector G4LorentzVector