Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4InuclParticle.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// $Id$
27//
28// 20100112 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29// 20100409 M. Kelsey -- Drop unused string argument from ctors.
30// 20100519 M. Kelsey -- Add public access to G4DynamicParticle content
31// 20100715 M. Kelsey -- Add setKineticEnergy() function
32// 20100915 M. Kelsey -- Add constructor to copy G4DynamicParticle input
33// 20110214 M. Kelsey -- Replace integer "model" with enum
34// 20110225 M. Kelsey -- Add equality operator (NOT sorting!)
35// 20110721 M. Kelsey -- Add model ID as optional ctor argument (so subclasses
36// don't have to call SetModel()).
37// 20110917 M. Kelsey -- Add coalesence to model ID list
38// 20110919 M. Kelsey -- Move setDefinition code to .cc to allow null argument
39// 20110922 M. Kelsey -- Add stream argument to printParticle() => print(),
40// add operator<<().
41
42#ifndef G4INUCL_PARTICLE_HH
43#define G4INUCL_PARTICLE_HH
44
45#include <iosfwd>
47
48#include "globals.hh"
49#include "G4LorentzVector.hh"
50#include "G4DynamicParticle.hh"
51
53public:
54 // used to indicate model that created instance of G4InuclParticle
55 // 0 default
56 // 1 bullet
57 // 2 target
58 // 3 G4ElementaryParticleCollider
59 // 4 G4IntraNucleiCascader
60 // 5 G4NonEquilibriumEvaporator
61 // 6 G4EquilibriumEvaporator
62 // 7 G4Fissioner
63 // 8 G4BigBanger
64 // 9 G4PreCompound
65 // 10 G4CascadeCoalescence
69
70public:
72
74 : pDP(dynPart), modelId(model) {}
75
77 : modelId(model) { pDP.Set4Momentum(mom*CLHEP::GeV/CLHEP::MeV); } // Bertini to G4 units
78
79 virtual ~G4InuclParticle() {}
80
81 // Copy and assignment constructors for use with std::vector<>
83 : pDP(right.pDP), modelId(right.modelId) {}
84
86
87 // Equality (comparison) operator -- NOT SORTING
88 bool operator==(const G4InuclParticle& right) {
89 return ( (&right == this) || (pDP == right.pDP) ); // Ignore model code
90 }
91
92 bool operator!=(const G4InuclParticle& right) {
93 return !operator==(right);
94 }
95
96 // This is no longer required, as setMomentum() handles mass adjustment
97 void setEnergy() { ; }
98
99 // These are call-throughs to G4DynamicParticle
100 void setMomentum(const G4LorentzVector& mom);
101
102 void setKineticEnergy(G4double ekin) { pDP.SetKineticEnergy(ekin*CLHEP::GeV/CLHEP::MeV); }
103
104 void setMass(G4double mass) { pDP.SetMass(mass*CLHEP::GeV/CLHEP::MeV); }
105
107 return pDP.GetMass()*CLHEP::MeV/CLHEP::GeV; // From G4 to Bertini units
108 }
109
111 return pDP.GetCharge();
112 }
113
115 return pDP.GetKineticEnergy()*CLHEP::MeV/CLHEP::GeV; // From G4 to Bertini units
116 }
117
119 return pDP.GetTotalEnergy()*CLHEP::MeV/CLHEP::GeV; // From G4 to Bertini units
120 }
121
123 return pDP.GetTotalMomentum()*CLHEP::MeV/CLHEP::GeV; // From G4 to Bertini units
124 }
125
127 return pDP.Get4Momentum()*CLHEP::MeV/CLHEP::GeV; // From G4 to Bertini units
128 }
129
130 virtual void print(std::ostream& os) const;
131
133 return pDP.GetDefinition();
134 }
135
137 return pDP;
138 }
139
140public:
141 void setModel(Model model) { modelId = model; }
142 Model getModel() const { return modelId; }
143
144protected:
145 // Special constructors for subclasses to set particle type correctly
147 : modelId(model) { setDefinition(pd); }
148
149 // FIXME: Bertini code doesn't pass valid 4-vectors, so force mass value
150 // from supplied PartDefn, with required unit conversions
152 Model model=DefaultModel);
153
154 // NOTE: Momentum forced along Z direction
156 Model model=DefaultModel)
157 : pDP(pd,G4ThreeVector(0.,0.,1.),ekin*CLHEP::GeV/CLHEP::MeV), modelId(model) {}
158
160
161private:
162 G4DynamicParticle pDP; // Carries all the kinematics and info
163 Model modelId;
164};
165
166// Proper stream output (just calls print())
167
168std::ostream& operator<<(std::ostream& os, const G4InuclParticle& part);
169
170#endif // G4INUCL_PARTICLE_HH
std::ostream & operator<<(std::ostream &os, const G4InuclParticle &part)
double G4double
Definition: G4Types.hh:64
G4double GetMass() const
G4double GetCharge() const
void SetMass(G4double mass)
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetTotalMomentum() const
void SetKineticEnergy(G4double aEnergy)
G4InuclParticle(const G4LorentzVector &mom, Model model=DefaultModel)
bool operator==(const G4InuclParticle &right)
G4InuclParticle(const G4InuclParticle &right)
void setMass(G4double mass)
G4ParticleDefinition * getDefinition() const
virtual void print(std::ostream &os) const
G4double getKineticEnergy() const
G4InuclParticle & operator=(const G4InuclParticle &right)
G4double getMass() const
G4LorentzVector getMomentum() const
G4double getMomModule() const
void setKineticEnergy(G4double ekin)
G4double getEnergy() const
G4InuclParticle(G4ParticleDefinition *pd, Model model=DefaultModel)
void setMomentum(const G4LorentzVector &mom)
bool operator!=(const G4InuclParticle &right)
virtual ~G4InuclParticle()
const G4DynamicParticle & getDynamicParticle() const
G4InuclParticle(const G4DynamicParticle &dynPart, Model model=DefaultModel)
G4double getCharge() const
void setDefinition(G4ParticleDefinition *pd)
Model getModel() const
G4InuclParticle(G4ParticleDefinition *pd, G4double ekin, Model model=DefaultModel)
void setModel(Model model)
Definition: DoubConv.h:17