Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4INCLParticle.cc
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// INCL++ intra-nuclear cascade model
27// Alain Boudard, CEA-Saclay, France
28// Joseph Cugnon, University of Liege, Belgium
29// Jean-Christophe David, CEA-Saclay, France
30// Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31// Sylvie Leray, CEA-Saclay, France
32// Davide Mancusi, CEA-Saclay, France
33//
34#define INCLXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
38/*
39 * Particle.cc
40 *
41 * \date Jun 5, 2009
42 * \author Pekka Kaitaniemi
43 */
44
45#include "G4INCLParticle.hh"
47
48namespace G4INCL {
49
50#ifdef INCLXX_IN_GEANT4_MODE
51 std::vector<G4double> Particle::INCLBiasVector;
52#else
53 G4ThreadLocal std::vector<G4double> Particle::INCLBiasVector;
54 //G4VectorCache<G4double> Particle::INCLBiasVector;
55#endif
56 G4ThreadLocal long Particle::nextID = 1;
58
60 : theZ(0), theA(0), theS(0),
61 theParticipantType(TargetSpectator),
62 theType(UnknownParticle),
63 theEnergy(0.0),
64 thePropagationEnergy(&theEnergy),
65 theFrozenEnergy(theEnergy),
66 theMomentum(ThreeVector(0.,0.,0.)),
67 thePropagationMomentum(&theMomentum),
68 theFrozenMomentum(theMomentum),
69 thePosition(ThreeVector(0.,0.,0.)),
70 nCollisions(0),
71 nDecays(0),
72 thePotentialEnergy(0.0),
73 rpCorrelated(false),
74 uncorrelatedMomentum(0.),
75 theParticleBias(1.),
76 theNKaon(0),
77 theParentResonancePDGCode(0),
78 theParentResonanceID(0),
79 theHelicity(0.0),
80 emissionTime(0.0),
81 outOfWell(false),
82 theMass(0.)
83 {
84 ID = nextID;
85 nextID++;
86 }
87
89 ThreeVector const &momentum, ThreeVector const &position)
90 : theEnergy(energy),
91 thePropagationEnergy(&theEnergy),
92 theFrozenEnergy(theEnergy),
93 theMomentum(momentum),
94 thePropagationMomentum(&theMomentum),
95 theFrozenMomentum(theMomentum),
96 thePosition(position),
97 nCollisions(0), nDecays(0),
98 thePotentialEnergy(0.),
99 rpCorrelated(false),
100 uncorrelatedMomentum(theMomentum.mag()),
101 theParticleBias(1.),
102 theNKaon(0),
103 theParentResonancePDGCode(0),
104 theParentResonanceID(0),
105 theHelicity(0.0),
106 emissionTime(0.0), outOfWell(false)
107 {
109 ID = nextID;
110 nextID++;
111 if(theEnergy <= 0.0) {
112 INCL_WARN("Particle with energy " << theEnergy << " created." << '\n');
113 }
114 setType(t);
116 }
117
119 ThreeVector const &momentum, ThreeVector const &position)
120 : thePropagationEnergy(&theEnergy),
121 theMomentum(momentum),
122 thePropagationMomentum(&theMomentum),
123 theFrozenMomentum(theMomentum),
124 thePosition(position),
125 nCollisions(0), nDecays(0),
126 thePotentialEnergy(0.),
127 rpCorrelated(false),
128 uncorrelatedMomentum(theMomentum.mag()),
129 theParticleBias(1.),
130 theNKaon(0),
131 theParentResonancePDGCode(0),
132 theParentResonanceID(0),
133 theHelicity(0.0),
134 emissionTime(0.0), outOfWell(false)
135 {
137 ID = nextID;
138 nextID++;
139 setType(t);
140 if( isResonance() ) {
141 INCL_ERROR("Cannot create resonance without specifying its momentum four-vector." << '\n');
142 }
143 G4double energy = std::sqrt(theMomentum.mag2() + theMass*theMass);
144 theEnergy = energy;
146 }
147
149 const G4double p2 = theMomentum.mag2();
150 G4double newp2 = theEnergy*theEnergy - theMass*theMass;
151 if( newp2<0.0 ) {
152 INCL_ERROR("Particle has E^2 < m^2." << '\n' << print());
153 newp2 = 0.0;
154 theEnergy = theMass;
155 }
156
157 theMomentum *= std::sqrt(newp2/p2);
158 return theMomentum;
159 }
160
162 theEnergy = std::sqrt(theMomentum.mag2() + theMass*theMass);
163 return theEnergy;
164 }
165
166 void ParticleList::rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) const {
167 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
168 (*i)->rotatePositionAndMomentum(angle, axis);
169 }
170 }
171
172 void ParticleList::rotatePosition(const G4double angle, const ThreeVector &axis) const {
173 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
174 (*i)->rotatePosition(angle, axis);
175 }
176 }
177
178 void ParticleList::rotateMomentum(const G4double angle, const ThreeVector &axis) const {
179 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
180 (*i)->rotateMomentum(angle, axis);
181 }
182 }
183
184 void ParticleList::boost(const ThreeVector &b) const {
185 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
186 (*i)->boost(b);
187 }
188 }
189
191 if(G4int((*this).size())==0) return 1.;
192 std::vector<G4int> MergedVector;
193 for(ParticleIter i = (*this).begin(), e = (*this).end(); i!=e; ++i){
194 MergedVector = Particle::MergeVectorBias(MergedVector,(*i));
195 }
196 return Particle::getBiasFromVector(MergedVector);
197 }
198
199 std::vector<G4int> ParticleList::getParticleListBiasVector() const {
200 std::vector<G4int> MergedVector;
201 if(G4int((*this).size())==0) return MergedVector;
202 for(ParticleIter i = (*this).begin(), e = (*this).end(); i!=e; ++i){
203 MergedVector = Particle::MergeVectorBias(MergedVector,(*i));
204 }
205 return MergedVector;
206 }
207
209// assert(G4int(Particle::INCLBiasVector.size())==nextBiasedCollisionID);
210 //assert(G4int(Particle::INCLBiasVector.Size())==nextBiasedCollisionID);
211// assert(std::fabs(newBias - 1.) > 1E-6);
212 Particle::INCLBiasVector.push_back(newBias);
213 //Particle::INCLBiasVector.Push_back(newBias);
215 }
216
217 G4double Particle::getBiasFromVector(std::vector<G4int> VectorBias) {
218 if(VectorBias.empty()) return 1.;
219
220 G4double ParticleBias = 1.;
221
222 for(G4int i=0; i<G4int(VectorBias.size()); i++){
223 ParticleBias *= Particle::INCLBiasVector[G4int(VectorBias[i])];
224 }
225
226 return ParticleBias;
227 }
228
229 std::vector<G4int> Particle::MergeVectorBias(Particle const * const p1, Particle const * const p2){
230 std::vector<G4int> MergedVectorBias;
231 std::vector<G4int> VectorBias1 = p1->getBiasCollisionVector();
232 std::vector<G4int> VectorBias2 = p2->getBiasCollisionVector();
233 G4int i = 0;
234 G4int j = 0;
235 if(VectorBias1.size()==0 && VectorBias2.size()==0) return MergedVectorBias;
236 else if(VectorBias1.size()==0) return VectorBias2;
237 else if(VectorBias2.size()==0) return VectorBias1;
238
239 while(i < G4int(VectorBias1.size()) || j < G4int(VectorBias2.size())){
240 if(VectorBias1[i]==VectorBias2[j]){
241 MergedVectorBias.push_back(VectorBias1[i]);
242 i++;
243 j++;
244 if(i == G4int(VectorBias1.size())){
245 for(;j<G4int(VectorBias2.size());j++) MergedVectorBias.push_back(VectorBias2[j]);
246 }
247 else if(j == G4int(VectorBias2.size())){
248 for(;i<G4int(VectorBias1.size());i++) MergedVectorBias.push_back(VectorBias1[i]);
249 }
250 } else if(VectorBias1[i]<VectorBias2[j]){
251 MergedVectorBias.push_back(VectorBias1[i]);
252 i++;
253 if(i == G4int(VectorBias1.size())){
254 for(;j<G4int(VectorBias2.size());j++) MergedVectorBias.push_back(VectorBias2[j]);
255 }
256 }
257 else {
258 MergedVectorBias.push_back(VectorBias2[j]);
259 j++;
260 if(j == G4int(VectorBias2.size())){
261 for(;i<G4int(VectorBias1.size());i++) MergedVectorBias.push_back(VectorBias1[i]);
262 }
263 }
264 }
265 return MergedVectorBias;
266 }
267
268 std::vector<G4int> Particle::MergeVectorBias(std::vector<G4int> p1, Particle const * const p2){
269 std::vector<G4int> MergedVectorBias;
270 std::vector<G4int> VectorBias = p2->getBiasCollisionVector();
271 G4int i = 0;
272 G4int j = 0;
273 if(p1.size()==0 && VectorBias.size()==0) return MergedVectorBias;
274 else if(p1.size()==0) return VectorBias;
275 else if(VectorBias.size()==0) return p1;
276
277 while(i < G4int(p1.size()) || j < G4int(VectorBias.size())){
278 if(p1[i]==VectorBias[j]){
279 MergedVectorBias.push_back(p1[i]);
280 i++;
281 j++;
282 if(i == G4int(p1.size())){
283 for(;j<G4int(VectorBias.size());j++) MergedVectorBias.push_back(VectorBias[j]);
284 }
285 else if(j == G4int(VectorBias.size())){
286 for(;i<G4int(p1.size());i++) MergedVectorBias.push_back(p1[i]);
287 }
288 } else if(p1[i]<VectorBias[j]){
289 MergedVectorBias.push_back(p1[i]);
290 i++;
291 if(i == G4int(p1.size())){
292 for(;j<G4int(VectorBias.size());j++) MergedVectorBias.push_back(VectorBias[j]);
293 }
294 }
295 else {
296 MergedVectorBias.push_back(VectorBias[j]);
297 j++;
298 if(j == G4int(VectorBias.size())){
299 for(;i<G4int(p1.size());i++) MergedVectorBias.push_back(p1[i]);
300 }
301 }
302 }
303 return MergedVectorBias;
304 }
305
307 G4double TotalBias = 1.;
308 for(G4int i=0; i<G4int(INCLBiasVector.size());i++) TotalBias *= Particle::INCLBiasVector[i];
309 return TotalBias;
310 }
311
312 void Particle::setINCLBiasVector(std::vector<G4double> NewVector) {
313 Particle::INCLBiasVector = NewVector;
314 }
315}
#define INCL_ERROR(x)
#define INCL_WARN(x)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double getParticleListBias() const
std::vector< G4int > getParticleListBiasVector() const
void rotateMomentum(const G4double angle, const ThreeVector &axis) const
void boost(const ThreeVector &b) const
void rotatePosition(const G4double angle, const ThreeVector &axis) const
void rotatePositionAndMomentum(const G4double angle, const ThreeVector &axis) const
G4INCL::ThreeVector theMomentum
static std::vector< G4double > INCLBiasVector
Time ordered vector of all bias applied.
void setMass(G4double mass)
ParticipantType theParticipantType
static void FillINCLBiasVector(G4double newBias)
static std::vector< G4int > MergeVectorBias(Particle const *const p1, Particle const *const p2)
std::vector< G4int > getBiasCollisionVector() const
Get the vector list of biased vertices on the particle path.
G4double adjustEnergyFromMomentum()
Recompute the energy to match the momentum.
static void setINCLBiasVector(std::vector< G4double > NewVector)
static G4double getTotalBias()
General bias vector function.
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
G4double getInvariantMass() const
Get the the particle invariant mass.
G4bool isResonance() const
Is it a resonance?
void setType(ParticleType t)
std::string print() const
static G4ThreadLocal G4int nextBiasedCollisionID
static G4double getBiasFromVector(std::vector< G4int > VectorBias)
G4double mag2() const
ParticleList::const_iterator ParticleIter
#define G4ThreadLocal
Definition: tls.hh:77