Geant4_10
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 // Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28 // Davide Mancusi, CEA
29 // Alain Boudard, CEA
30 // Sylvie Leray, CEA
31 // Joseph Cugnon, University of Liege
32 //
33 #define INCLXX_IN_GEANT4_MODE 1
34 
35 #include "globals.hh"
36 
46 #include "G4INCLEventInfo.hh"
47 #include "G4INCLGlobals.hh"
48 #include "G4INCLParticleTable.hh"
49 #include <cmath>
50 
51 namespace G4INCL {
52 
54 
55 #ifdef INCL_INVERSE_KINEMATICS
56  void EventInfo::fillInverseKinematics(const Double_t gamma) {
57  const Double_t beta = std::sqrt(1.-1./(gamma*gamma));
58  for(Int_t i=0; i<nParticles; ++i) {
59  // determine the particle mass from the kinetic energy and the momentum;
60  // this ensures consistency with the masses uses by the models
61  const Double_t mass = std::max(
62  0.5 * (px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i]-EKin[i]*EKin[i]) / EKin[i],
63  0.0);
64 
65  const Double_t ETot = EKin[i] + mass;
66  const Double_t ETotPrime = gamma*(ETot - beta*pz[i]);
67  EKinPrime[i] = ETotPrime - mass;
68  pzPrime[i] = -gamma*(pz[i] - beta*ETot);
69  const Double_t pPrime = std::sqrt(px[i]*px[i] + py[i]*py[i] + pzPrime[i]*pzPrime[i]);
70  const Double_t cosThetaPrime = pzPrime[i]/pPrime;
71  if(cosThetaPrime>=1.)
72  thetaPrime[i] = 0.;
73  else if(cosThetaPrime<=-1.)
74  thetaPrime[i] = 180.;
75  else
76  thetaPrime[i] = 180.*std::acos(cosThetaPrime)/Math::pi;
77  }
78  }
79 #endif // INCL_INVERSE_KINEMATICS
80 
81  void EventInfo::remnantToParticle(const G4int remnantIndex) {
82  A[nParticles] = ARem[remnantIndex];
83  Z[nParticles] = ZRem[remnantIndex];
85 
86  px[nParticles] = pxRem[remnantIndex];
87  py[nParticles] = pyRem[remnantIndex];
88  pz[nParticles] = pzRem[remnantIndex];
89 
90  const G4double plab = std::sqrt(pxRem[remnantIndex]*pxRem[remnantIndex]
91  +pyRem[remnantIndex]*pyRem[remnantIndex]
92  +pzRem[remnantIndex]*pzRem[remnantIndex]);
93  G4double pznorm = pzRem[remnantIndex]/plab;
94  if(pznorm>1.)
95  pznorm = 1.;
96  else if(pznorm<-1.)
97  pznorm = -1.;
98  theta[nParticles] = 180.*std::acos(pznorm)/G4INCL::Math::pi;
99  phi[nParticles] = 180.*std::atan2(pyRem[remnantIndex],pxRem[remnantIndex])/G4INCL::Math::pi;
100 
101  EKin[nParticles] = EKinRem[remnantIndex];
102  origin[nParticles] = -1; // Origin: cascade
103  history.push_back(""); // history
104  nParticles++;
105 // assert(history.size()==(unsigned int)nParticles);
106  }
107 }
108 
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
const G4double pi
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
#define G4ThreadLocal
Definition: tls.hh:52
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
int G4int
Definition: G4Types.hh:78
Float_t stoppingTime
Cascade stopping time [fm/c].
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
static G4ThreadLocal Int_t eventNumber
Number of the event.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
Simple container for output of event results.
Short_t Z[maxSizeParticles]
Particle charge number.
Short_t nParticles
Number of particles in the final state.
G4int Int_t
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
void remnantToParticle(const G4int remnantIndex)
Move a remnant to the particle array.
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
Short_t A[maxSizeParticles]
Particle mass number.
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
G4double Double_t
double G4double
Definition: G4Types.hh:76
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
Short_t origin[maxSizeParticles]
Origin of the particle.
std::vector< std::string > history
History of the particle.