Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // 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  Double_t mass;
60  if(A[i]>0) {
61  mass = ParticleTable::getTableMass(A[i],Z[i]);
62  } else if(origin[i]==-1) { // cascade particles with A=0, must be pions
63  if(Z[i]==1)
65  else if(Z[i]==0)
67  else
69  } else // gamma rays
70  mass = 0.;
71 
72  const Double_t ETot = EKin[i] + mass;
73  const Double_t ETotPrime = gamma*(ETot - beta*pz[i]);
74  /* Using the invariant mass here avoids negative kinetic energies with
75  * particles produced by the de-excitation models, which do not
76  * necessarily use the same mass look-up tables as INCL.
77  */
78  Double_t invariantMass;
79  if(A[i]>0 || origin[i]==-1) { // massive particles
80  invariantMass = std::sqrt(ETot*ETot - px[i]*px[i] - py[i]*py[i] - pz[i]*pz[i]);
81  } else { // gamma rays
82  invariantMass = 0.;
83  }
84  EKinPrime[i] = ETotPrime - invariantMass;
85  pzPrime[i] = -gamma*(pz[i] - beta*ETot);
86  const Double_t pPrime = std::sqrt(px[i]*px[i] + py[i]*py[i] + pzPrime[i]*pzPrime[i]);
87  const Double_t cosThetaPrime = pzPrime[i]/pPrime;
88  if(cosThetaPrime>=1.)
89  thetaPrime[i] = 0.;
90  else if(cosThetaPrime<=-1.)
91  thetaPrime[i] = 180.;
92  else
93  thetaPrime[i] = 180.*std::acos(cosThetaPrime)/Math::pi;
94  }
95  }
96 #endif // INCL_INVERSE_KINEMATICS
97 }
98