Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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),
63 theEnergy(0.0),
66 theMomentum(ThreeVector(0.,0.,0.)),
69 thePosition(ThreeVector(0.,0.,0.)),
70 nCollisions(0),
71 nDecays(0),
72 nSrcPair(0),
74 rpCorrelated(false),
77 theNKaon(0),
81#endif
82 theHelicity(0.0),
83 emissionTime(0.0),
84 outOfWell(false),
85 theSrcPartner(false),
86 theMass(0.)
87 {
88 ID = nextID;
89 nextID++;
90 }
91
93 ThreeVector const &momentum, ThreeVector const &position)
94 : theEnergy(energy),
97 theMomentum(momentum),
101 nCollisions(0), nDecays(0), nSrcPair(0),
103 rpCorrelated(false),
105 theParticleBias(1.),
106 theNKaon(0),
110#endif
111 theHelicity(0.0),
112 emissionTime(0.0), outOfWell(false), theSrcPartner(false)
113 {
115 ID = nextID;
116 nextID++;
117 if(theEnergy <= 0.0) {
118 INCL_WARN("Particle with energy " << theEnergy << " created." << '\n');
119 }
120 setType(t);
122 }
123
125 ThreeVector const &momentum, ThreeVector const &position)
127 theMomentum(momentum),
131 nCollisions(0), nDecays(0),
132 nSrcPair(0),
134 rpCorrelated(false),
136 theParticleBias(1.),
137 theNKaon(0),
141#endif
142 theHelicity(0.0),
143 emissionTime(0.0), outOfWell(false), theSrcPartner(false)
144 {
146 ID = nextID;
147 nextID++;
148 setType(t);
149 if( isResonance() ) {
150 INCL_ERROR("Cannot create resonance without specifying its momentum four-vector." << '\n');
151 }
152 G4double energy = std::sqrt(theMomentum.mag2() + theMass*theMass);
153 theEnergy = energy;
155 }
156
158 const G4double p2 = theMomentum.mag2();
159 G4double newp2 = theEnergy*theEnergy - theMass*theMass;
160 if( newp2<0.0 ) {
161 INCL_ERROR("Particle has E^2 < m^2." << '\n' << print());
162 newp2 = 0.0;
163 theEnergy = theMass;
164 }
165
166 theMomentum *= std::sqrt(newp2/p2);
167 return theMomentum;
168 }
169
171 theEnergy = std::sqrt(theMomentum.mag2() + theMass*theMass);
172 return theEnergy;
173 }
174
176 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
177 (*i)->rotatePositionAndMomentum(angle, axis);
178 }
179 }
180
181 void ParticleList::rotatePosition(const G4double angle, const ThreeVector &axis) const {
182 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
183 (*i)->rotatePosition(angle, axis);
184 }
185 }
186
187 void ParticleList::rotateMomentum(const G4double angle, const ThreeVector &axis) const {
188 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
189 (*i)->rotateMomentum(angle, axis);
190 }
191 }
192
193 void ParticleList::boost(const ThreeVector &b) const {
194 for(const_iterator i=begin(), e=end(); i!=e; ++i) {
195 (*i)->boost(b);
196 }
197 }
198
200 if(G4int((*this).size())==0) return 1.;
201 std::vector<G4int> MergedVector;
202 for(ParticleIter i = (*this).begin(), e = (*this).end(); i!=e; ++i){
203 MergedVector = Particle::MergeVectorBias(MergedVector,(*i));
204 }
205 return Particle::getBiasFromVector(std::move(MergedVector));
206 }
207
208 std::vector<G4int> ParticleList::getParticleListBiasVector() const {
209 std::vector<G4int> MergedVector;
210 if(G4int((*this).size())==0) return MergedVector;
211 for(ParticleIter i = (*this).begin(), e = (*this).end(); i!=e; ++i){
212 MergedVector = Particle::MergeVectorBias(MergedVector,(*i));
213 }
214 return MergedVector;
215 }
216
218// assert(G4int(Particle::INCLBiasVector.size())==nextBiasedCollisionID);
219 //assert(G4int(Particle::INCLBiasVector.Size())==nextBiasedCollisionID);
220// assert(std::fabs(newBias - 1.) > 1E-6);
221 Particle::INCLBiasVector.push_back(newBias);
222 //Particle::INCLBiasVector.Push_back(newBias);
224 }
225
226 G4double Particle::getBiasFromVector(std::vector<G4int> VectorBias) {
227 if(VectorBias.empty()) return 1.;
228
229 G4double ParticleBias = 1.;
230
231 for(G4int i=0; i<G4int(VectorBias.size()); i++){
232 ParticleBias *= Particle::INCLBiasVector[G4int(VectorBias[i])];
233 }
234
235 return ParticleBias;
236 }
237
238 std::vector<G4int> Particle::MergeVectorBias(Particle const * const p1, Particle const * const p2){
239 std::vector<G4int> MergedVectorBias;
240 std::vector<G4int> VectorBias1 = p1->getBiasCollisionVector();
241 std::vector<G4int> VectorBias2 = p2->getBiasCollisionVector();
242 G4int i = 0;
243 G4int j = 0;
244 if(VectorBias1.size()==0 && VectorBias2.size()==0) return MergedVectorBias;
245 else if(VectorBias1.size()==0) return VectorBias2;
246 else if(VectorBias2.size()==0) return VectorBias1;
247
248 while(i < G4int(VectorBias1.size()) || j < G4int(VectorBias2.size())){
249 if(VectorBias1[i]==VectorBias2[j]){
250 MergedVectorBias.push_back(VectorBias1[i]);
251 i++;
252 j++;
253 if(i == G4int(VectorBias1.size())){
254 for(;j<G4int(VectorBias2.size());j++) MergedVectorBias.push_back(VectorBias2[j]);
255 }
256 else if(j == G4int(VectorBias2.size())){
257 for(;i<G4int(VectorBias1.size());i++) MergedVectorBias.push_back(VectorBias1[i]);
258 }
259 } else if(VectorBias1[i]<VectorBias2[j]){
260 MergedVectorBias.push_back(VectorBias1[i]);
261 i++;
262 if(i == G4int(VectorBias1.size())){
263 for(;j<G4int(VectorBias2.size());j++) MergedVectorBias.push_back(VectorBias2[j]);
264 }
265 }
266 else {
267 MergedVectorBias.push_back(VectorBias2[j]);
268 j++;
269 if(j == G4int(VectorBias2.size())){
270 for(;i<G4int(VectorBias1.size());i++) MergedVectorBias.push_back(VectorBias1[i]);
271 }
272 }
273 }
274 return MergedVectorBias;
275 }
276
277 std::vector<G4int> Particle::MergeVectorBias(std::vector<G4int> p1, Particle const * const p2){
278 std::vector<G4int> MergedVectorBias;
279 std::vector<G4int> VectorBias = p2->getBiasCollisionVector();
280 G4int i = 0;
281 G4int j = 0;
282 if(p1.size()==0 && VectorBias.size()==0) return MergedVectorBias;
283 else if(p1.size()==0) return VectorBias;
284 else if(VectorBias.size()==0) return p1;
285
286 while(i < G4int(p1.size()) || j < G4int(VectorBias.size())){
287 if(p1[i]==VectorBias[j]){
288 MergedVectorBias.push_back(p1[i]);
289 i++;
290 j++;
291 if(i == G4int(p1.size())){
292 for(;j<G4int(VectorBias.size());j++) MergedVectorBias.push_back(VectorBias[j]);
293 }
294 else if(j == G4int(VectorBias.size())){
295 for(;i<G4int(p1.size());i++) MergedVectorBias.push_back(p1[i]);
296 }
297 } else if(p1[i]<VectorBias[j]){
298 MergedVectorBias.push_back(p1[i]);
299 i++;
300 if(i == G4int(p1.size())){
301 for(;j<G4int(VectorBias.size());j++) MergedVectorBias.push_back(VectorBias[j]);
302 }
303 }
304 else {
305 MergedVectorBias.push_back(VectorBias[j]);
306 j++;
307 if(j == G4int(VectorBias.size())){
308 for(;i<G4int(p1.size());i++) MergedVectorBias.push_back(p1[i]);
309 }
310 }
311 }
312 return MergedVectorBias;
313 }
314
316 G4double TotalBias = 1.;
317 for(G4int i=0; i<G4int(INCLBiasVector.size());i++) TotalBias *= Particle::INCLBiasVector[i];
318 return TotalBias;
319 }
320
321 void Particle::setINCLBiasVector(std::vector<G4double> NewVector) {
322 Particle::INCLBiasVector = std::move(NewVector);
323 }
324}
#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 * thePropagationMomentum
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)
G4INCL::ParticleType theType
static std::vector< G4int > MergeVectorBias(Particle const *const p1, Particle const *const p2)
G4int theNKaon
The number of Kaons inside the nucleus (update during the cascade).
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 uncorrelatedMomentum
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
G4double * thePropagationEnergy
static G4double getBiasFromVector(std::vector< G4int > VectorBias)
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theFrozenMomentum
ParticleList::const_iterator ParticleIter
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
#define G4ThreadLocal
Definition tls.hh:77