Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
47 #include "G4INCLEventInfo.hh"
48 #include "G4INCLGlobals.hh"
49 #include "G4INCLParticleTable.hh"
50 #include <cmath>
51 
52 namespace G4INCL {
53 
55 
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  Double_t mass;
62  if(EKin[i]>0.) {
63  mass = std::max(
64  0.5 * (px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i]-EKin[i]*EKin[i]) / EKin[i],
65  0.0);
66  } else {
67  INCL_WARN("Particle with null kinetic energy in fillInverseKinematics, cannot determine its mass:\n"
68  << " A=" << A[i] << ", Z=" << Z[i] << '\n'
69  << " EKin=" << EKin[i] << ", px=" << px[i] << ", py=" << py[i] << ", pz=" << pz[i] << '\n'
70  << " Falling back to the mass from the INCL ParticleTable" << '\n');
71  mass = ParticleTable::getRealMass(A[i], Z[i]);
72  }
73 
74  const Double_t ETot = EKin[i] + mass;
75  const Double_t ETotPrime = gamma*(ETot - beta*pz[i]);
76  EKinPrime[i] = ETotPrime - mass;
77  pzPrime[i] = -gamma*(pz[i] - beta*ETot);
78  const Double_t pPrime = std::sqrt(px[i]*px[i] + py[i]*py[i] + pzPrime[i]*pzPrime[i]);
79  const Double_t cosThetaPrime = (pPrime>0.) ? (pzPrime[i]/pPrime) : 1.;
80  if(cosThetaPrime>=1.)
81  thetaPrime[i] = 0.;
82  else if(cosThetaPrime<=-1.)
83  thetaPrime[i] = 180.;
84  else
85  thetaPrime[i] = Math::toDegrees(Math::arcCos(cosThetaPrime));
86  }
87  }
88 
89  void EventInfo::remnantToParticle(const G4int remnantIndex) {
90  A[nParticles] = ARem[remnantIndex];
91  Z[nParticles] = ZRem[remnantIndex];
93 
94  px[nParticles] = pxRem[remnantIndex];
95  py[nParticles] = pyRem[remnantIndex];
96  pz[nParticles] = pzRem[remnantIndex];
97 
98  const G4double plab = std::sqrt(pxRem[remnantIndex]*pxRem[remnantIndex]
99  +pyRem[remnantIndex]*pyRem[remnantIndex]
100  +pzRem[remnantIndex]*pzRem[remnantIndex]);
101  G4double pznorm = pzRem[remnantIndex]/plab;
102  if(pznorm>1.)
103  pznorm = 1.;
104  else if(pznorm<-1.)
105  pznorm = -1.;
107  phi[nParticles] = Math::toDegrees(std::atan2(pyRem[remnantIndex],pxRem[remnantIndex]));
108 
109  EKin[nParticles] = EKinRem[remnantIndex];
110  origin[nParticles] = -1; // Origin: cascade
111  history.push_back(""); // history
112  nParticles++;
113 // assert(history.size()==(unsigned int)nParticles);
114  }
115 }
116 
G4double arcCos(const G4double x)
Calculates arccos with some tolerance on illegal arguments.
Float_t thetaPrime[maxSizeParticles]
Particle momentum polar angle, in inverse kinematics [radians].
void fillInverseKinematics(const Double_t gamma)
Fill the variables describing the reaction in inverse kinematics.
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
#define INCL_WARN(x)
G4double toDegrees(G4double radians)
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
#define G4ThreadLocal
Definition: tls.hh:89
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].
Float_t EKinPrime[maxSizeParticles]
Particle kinetic energy, in inverse kinematics [MeV].
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
static G4ThreadLocal Int_t eventNumber
Number of the event.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
Float_t pzPrime[maxSizeParticles]
Particle momentum, z component, in inverse kinematics [MeV/c].
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.