Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCLEventInfo.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/** \file G4INCLEventInfo.cc
39 * \brief Simple container for output of event results.
40 *
41 * Contains the results of an INCL cascade.
42 *
43 * \date 21 January 2011
44 * \author Davide Mancusi
45 */
46
47#include "G4INCLEventInfo.hh"
48#include "G4INCLGlobals.hh"
50#include "G4INCLParticle.hh"
51#include <cmath>
52
53namespace G4INCL {
54
56
58 const Double_t beta = std::sqrt(1.-1./(gamma*gamma));
59 for(Int_t i=0; i<nParticles; ++i) {
60 // determine the particle mass from the kinetic energy and the momentum;
61 // this ensures consistency with the masses uses by the models
62 Double_t mass;
63 if(EKin[i]>0.) {
64 mass = std::max(
65 0.5 * (px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i]-EKin[i]*EKin[i]) / EKin[i],
66 0.0);
67 } else {
68 INCL_WARN("Particle with null kinetic energy in fillInverseKinematics, cannot determine its mass:\n"
69 << " A=" << A[i] << ", Z=" << Z[i] << ", S=" << S[i] << '\n'
70 << " EKin=" << EKin[i] << ", px=" << px[i] << ", py=" << py[i] << ", pz=" << pz[i] << '\n'
71 << " Falling back to the mass from the INCL ParticleTable" << '\n');
72 mass = ParticleTable::getRealMass(A[i], Z[i], S[i]);
73 }
74
75 const Double_t ETot = EKin[i] + mass;
76 const Double_t ETotPrime = gamma*(ETot - beta*pz[i]);
77 EKinPrime[i] = ETotPrime - mass;
78 pzPrime[i] = -gamma*(pz[i] - beta*ETot);
79 const Double_t pPrime = std::sqrt(px[i]*px[i] + py[i]*py[i] + pzPrime[i]*pzPrime[i]);
80 const Double_t cosThetaPrime = (pPrime>0.) ? (pzPrime[i]/pPrime) : 1.;
81 if(cosThetaPrime>=1.)
82 thetaPrime[i] = 0.;
83 else if(cosThetaPrime<=-1.)
84 thetaPrime[i] = 180.;
85 else
86 thetaPrime[i] = Math::toDegrees(Math::arcCos(cosThetaPrime));
87 }
88 }
89
90 void EventInfo::remnantToParticle(const G4int remnantIndex) {
91
92 INCL_DEBUG("remnantToParticle function used\n");
93
94 A[nParticles] = ARem[remnantIndex];
95 Z[nParticles] = ZRem[remnantIndex];
96 S[nParticles] = SRem[remnantIndex];
97
100
103
104 px[nParticles] = pxRem[remnantIndex];
105 py[nParticles] = pyRem[remnantIndex];
106 pz[nParticles] = pzRem[remnantIndex];
107
108 const G4double plab = std::sqrt(pxRem[remnantIndex]*pxRem[remnantIndex]
109 +pyRem[remnantIndex]*pyRem[remnantIndex]
110 +pzRem[remnantIndex]*pzRem[remnantIndex]);
111 G4double pznorm = pzRem[remnantIndex]/plab;
112 if(pznorm>1.)
113 pznorm = 1.;
114 else if(pznorm<-1.)
115 pznorm = -1.;
117 phi[nParticles] = Math::toDegrees(std::atan2(pyRem[remnantIndex],pxRem[remnantIndex]));
118
119 EKin[nParticles] = EKinRem[remnantIndex];
120 origin[nParticles] = -1; // Origin: cascade
121 parentResonancePDGCode[nParticles] = 0; // No parent resonance
122 parentResonanceID[nParticles] = 0; // No parent resonance
123 history.push_back(""); // history
124 nParticles++;
125// assert(history.size()==(unsigned int)nParticles);
126 }
127}
128
Simple container for output of event results.
#define INCL_WARN(x)
#define INCL_DEBUG(x)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4int getPDGCode() const
Set a PDG Code (MONTE CARLO PARTICLE NUMBERING)
static G4double getTotalBias()
General bias vector function.
G4double arcCos(const G4double x)
Calculates arccos with some tolerance on illegal arguments.
G4double toDegrees(G4double radians)
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
G4int Int_t
G4double Double_t
Short_t S[maxSizeParticles]
Particle strangeness number.
Short_t origin[maxSizeParticles]
Origin of the particle.
Int_t parentResonancePDGCode[maxSizeParticles]
Particle's parent resonance PDG code.
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Short_t A[maxSizeParticles]
Particle mass number.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
Float_t pzPrime[maxSizeParticles]
Particle momentum, z component, in inverse kinematics [MeV/c].
Int_t parentResonanceID[maxSizeParticles]
Particle's parent resonance unique ID identifier.
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Float_t stoppingTime
Cascade stopping time [fm/c].
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
Short_t Z[maxSizeParticles]
Particle charge number.
std::vector< std::string > history
History of the particle.
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
Short_t nParticles
Number of particles in the final state.
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Float_t ParticleBias[maxSizeParticles]
Particle weight due to the bias.
Short_t SRem[maxSizeRemnants]
Remnant strangeness number.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Float_t thetaPrime[maxSizeParticles]
Particle momentum polar angle, in inverse kinematics [radians].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Float_t EKinPrime[maxSizeParticles]
Particle kinetic energy, in inverse kinematics [MeV].
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
static G4ThreadLocal Int_t eventNumber
Number of the event.
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Int_t PDGCode[maxSizeParticles]
PDG numbering of the particles.
void fillInverseKinematics(const Double_t gamma)
Fill the variables describing the reaction in inverse kinematics.
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
void remnantToParticle(const G4int remnantIndex)
Move a remnant to the particle array.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
#define G4ThreadLocal
Definition: tls.hh:77