ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4INCLEventInfo.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4INCLEventInfo.cc
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 "G4INCLParticle.hh"
51 #include <cmath>
52 
53 namespace 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 
98  ParticleSpecies pt(A[nParticles],Z[nParticles],S[nParticles]);
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  history.push_back(""); // history
122  nParticles++;
123 // assert(history.size()==(unsigned int)nParticles);
124  }
125 }
126