Geant4 Cross Reference |
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.hh 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 #ifndef G4INCLEVENTINFO_HH_HH 48 #define G4INCLEVENTINFO_HH_HH 1 49 50 #include "G4INCLParticleType.hh" 51 #ifdef INCL_ROOT_USE 52 #include <Rtypes.h> 53 #endif 54 #include <string> 55 #include <vector> 56 #include <algorithm> 57 58 namespace G4INCL { 59 #ifndef INCL_ROOT_USE 60 typedef G4int Int_t; 61 typedef short Short_t; 62 typedef G4float Float_t; 63 typedef G4double Double_t; 64 typedef G4bool Bool_t; 65 #endif 66 67 struct EventInfo { 68 EventInfo() : 69 nParticles(0), 70 event(0), 71 eventBias((Float_t)0.0), 72 nRemnants(0), 73 projectileType(0), 74 At(0), 75 Zt(0), 76 St(0), 77 Ap(0), 78 Zp(0), 79 Sp(0), 80 Ep((Float_t)0.0), 81 impactParameter((Float_t)0.0), 82 nCollisions(0), 83 stoppingTime((Float_t)0.0), 84 EBalance((Float_t)0.0), 85 firstEBalance((Float_t)0.0), 86 pLongBalance((Float_t)0.0), 87 pTransBalance((Float_t)0.0), 88 nCascadeParticles(0), 89 transparent(false), 90 annihilationP(false), 91 annihilationN(false), 92 forcedCompoundNucleus(false), 93 nucleonAbsorption(false), 94 pionAbsorption(false), 95 nDecays(0), 96 nSrcCollisions(0), 97 nSrcPairs(0), 98 nBlockedCollisions(0), 99 nBlockedDecays(0), 100 effectiveImpactParameter((Float_t)0.0), 101 deltasInside(false), 102 sigmasInside(false), 103 kaonsInside(false), 104 antikaonsInside(false), 105 lambdasInside(false), 106 forcedDeltasInside(false), 107 forcedDeltasOutside(false), 108 forcedPionResonancesOutside(false), 109 absorbedStrangeParticle(false), 110 forcedSigmaOutside(false), 111 forcedStrangeInside(false), 112 emitLambda(0), 113 emitKaon(false), 114 clusterDecay(false), 115 firstCollisionTime((Float_t)0.0), 116 firstCollisionXSec((Float_t)0.0), 117 firstCollisionSpectatorPosition((Float_t)0.0), 118 firstCollisionSpectatorMomentum((Float_t)0.0), 119 firstCollisionIsElastic(false), 120 nReflectionAvatars(0), 121 nCollisionAvatars(0), 122 nDecayAvatars(0), 123 nUnmergedSpectators(0), 124 nEnergyViolationInteraction(0) 125 126 { 127 std::fill_n(A, maxSizeParticles, 0); 128 std::fill_n(Z, maxSizeParticles, 0); 129 std::fill_n(S, maxSizeParticles, 0); 130 std::fill_n(J, maxSizeParticles, 0); 131 std::fill_n(PDGCode, maxSizeParticles, 0); 132 std::fill_n(ParticleBias, maxSizeParticles, (Float_t)0.0); 133 std::fill_n(EKin, maxSizeParticles, (Float_t)0.0); 134 std::fill_n(px, maxSizeParticles, (Float_t)0.0); 135 std::fill_n(py, maxSizeParticles, (Float_t)0.0); 136 std::fill_n(pz, maxSizeParticles, (Float_t)0.0); 137 std::fill_n(theta, maxSizeParticles, (Float_t)0.0); 138 std::fill_n(phi, maxSizeParticles, (Float_t)0.0); 139 std::fill_n(origin, maxSizeParticles, 0); 140 std::fill_n(parentResonancePDGCode, maxSizeParticles, 0); 141 std::fill_n(parentResonanceID, maxSizeParticles, 0); 142 std::fill_n(emissionTime, maxSizeParticles, (Float_t)0.0); 143 std::fill_n(ARem, maxSizeRemnants, 0); 144 std::fill_n(ZRem, maxSizeRemnants, 0); 145 std::fill_n(SRem, maxSizeRemnants, 0); 146 std::fill_n(EStarRem, maxSizeRemnants, (Float_t)0.0); 147 std::fill_n(JRem, maxSizeRemnants, (Float_t)0.0); 148 std::fill_n(EKinRem, maxSizeRemnants, (Float_t)0.0); 149 std::fill_n(pxRem, maxSizeRemnants, (Float_t)0.0); 150 std::fill_n(pyRem, maxSizeRemnants, (Float_t)0.0); 151 std::fill_n(pzRem, maxSizeRemnants, (Float_t)0.0); 152 std::fill_n(thetaRem, maxSizeRemnants, (Float_t)0.0); 153 std::fill_n(phiRem, maxSizeRemnants, (Float_t)0.0); 154 std::fill_n(jxRem, maxSizeRemnants, (Float_t)0.0); 155 std::fill_n(jyRem, maxSizeRemnants, (Float_t)0.0); 156 std::fill_n(jzRem, maxSizeRemnants, (Float_t)0.0); 157 std::fill_n(EKinPrime, maxSizeParticles, (Float_t)0.0); 158 std::fill_n(pzPrime, maxSizeParticles, (Float_t)0.0); 159 std::fill_n(thetaPrime, maxSizeParticles, (Float_t)0.0); 160 } 161 162 /** \brief Number of the event */ 163 static G4ThreadLocal Int_t eventNumber; 164 165 /** \brief Maximum array size for remnants */ 166 static const Short_t maxSizeRemnants = 10; 167 168 /** \brief Maximum array size for emitted particles */ 169 static const Short_t maxSizeParticles = 1000; 170 171 /** \brief Number of particles in the final state */ 172 Short_t nParticles; 173 /** \brief Sequential number of the event in the event loop */ 174 Int_t event; 175 /** \brief Particle mass number */ 176 Short_t A[maxSizeParticles]; 177 /** \brief Particle charge number */ 178 Short_t Z[maxSizeParticles]; 179 /** \brief Particle strangeness number */ 180 Short_t S[maxSizeParticles]; 181 /** \brief Particle angular momemtum */ 182 Short_t J[maxSizeParticles]; 183 /** \brief PDG numbering of the particles */ 184 Int_t PDGCode[maxSizeParticles]; 185 /** \brief Event bias */ 186 Float_t eventBias; 187 /** \brief Particle weight due to the bias */ 188 Float_t ParticleBias[maxSizeParticles]; 189 /** \brief Particle kinetic energy [MeV] */ 190 Float_t EKin[maxSizeParticles]; 191 /** \brief Particle momentum, x component [MeV/c] */ 192 Float_t px[maxSizeParticles]; 193 /** \brief Particle momentum, y component [MeV/c] */ 194 Float_t py[maxSizeParticles]; 195 /** \brief Particle momentum, z component [MeV/c] */ 196 Float_t pz[maxSizeParticles]; 197 /** \brief Particle momentum polar angle [radians] */ 198 Float_t theta[maxSizeParticles]; 199 /** \brief Particle momentum azimuthal angle [radians] */ 200 Float_t phi[maxSizeParticles]; 201 /** \brief Origin of the particle 202 * 203 * Should be -1 for cascade particles, or the number of the remnant for 204 * de-excitation particles. */ 205 Short_t origin[maxSizeParticles]; 206 /** \brief Particle's parent resonance PDG code */ 207 Int_t parentResonancePDGCode[maxSizeParticles]; 208 /** \brief Particle's parent resonance unique ID identifier */ 209 Int_t parentResonanceID[maxSizeParticles]; 210 /** \brief History of the particle 211 * 212 * Condensed information about the de-excitation chain of a particle. For 213 * cascade particles, it is just an empty string. For particles arising 214 * from the de-excitation of a cascade remnant, it is a string of 215 * characters. Each character represents one or more identical steps in 216 * the de-excitation process. The currently defined possible character 217 * values and their meanings are the following: 218 * 219 * e: evaporation product 220 * E: evaporation residue 221 * m: multifragmentation 222 * a: light partner in asymmetric fission or IMF emission 223 * A: heavy partner in asymmetric fission or IMF emission 224 * f: light partner in fission 225 * F: heavy partner in fission 226 * s: saddle-to-scission emission 227 * n: non-statistical emission (decay) */ 228 std::vector<std::string> history; 229 /** \brief Number of remnants */ 230 Short_t nRemnants; 231 /** \brief Projectile particle type */ 232 Int_t projectileType; 233 /** \brief Mass number of the target nucleus */ 234 Short_t At; 235 /** \brief Charge number of the target nucleus */ 236 Short_t Zt; 237 /** \brief Strangeness number of the target nucleus */ 238 Short_t St; 239 /** \brief Mass number of the projectile nucleus */ 240 Short_t Ap; 241 /** \brief Charge number of the projectile nucleus */ 242 Short_t Zp; 243 /** \brief Strangeness number of the projectile nucleus */ 244 Short_t Sp; 245 /** \brief Projectile kinetic energy given as input */ 246 Float_t Ep; 247 /** \brief Impact parameter [fm] */ 248 Float_t impactParameter; 249 /** \brief Number of accepted two-body collisions */ 250 Int_t nCollisions; 251 /** \brief Cascade stopping time [fm/c] */ 252 Float_t stoppingTime; 253 /** \brief Energy-conservation balance [MeV] */ 254 Float_t EBalance; 255 /** \brief First value for the energy-conservation balance [MeV] */ 256 Float_t firstEBalance; 257 /** \brief Longitudinal momentum-conservation balance [MeV/c] */ 258 Float_t pLongBalance; 259 /** \brief Transverse momentum-conservation balance [MeV/c] */ 260 Float_t pTransBalance; 261 /** \brief Number of cascade particles */ 262 Short_t nCascadeParticles; 263 /** \brief True if the event is transparent */ 264 Bool_t transparent; 265 /** \brief True if annihilation at rest on a proton */ 266 Bool_t annihilationP; 267 /** \brief True if annihilation at rest on a neutron */ 268 Bool_t annihilationN; 269 /** \brief True if the event is a forced CN */ 270 Bool_t forcedCompoundNucleus; 271 /** \brief True if the event is a nucleon absorption */ 272 Bool_t nucleonAbsorption; 273 /** \brief True if the event is a pion absorption */ 274 Bool_t pionAbsorption; 275 /** \brief Number of accepted Delta decays */ 276 Int_t nDecays; 277 /** \brief Number of accepted SRC collisions */ 278 Int_t nSrcCollisions; 279 /** \brief Number of src pairs */ 280 Int_t nSrcPairs; 281 /** \brief Number of two-body collisions blocked by Pauli or CDPP */ 282 Int_t nBlockedCollisions; 283 /** \brief Number of decays blocked by Pauli or CDPP */ 284 Int_t nBlockedDecays; 285 /** \brief Effective (Coulomb-distorted) impact parameter [fm] */ 286 Float_t effectiveImpactParameter; 287 /** \brief Event involved deltas in the nucleus at the end of the cascade */ 288 Bool_t deltasInside; 289 /** \brief Event involved sigmas in the nucleus at the end of the cascade */ 290 Bool_t sigmasInside; 291 /** \brief Event involved kaons in the nucleus at the end of the cascade */ 292 Bool_t kaonsInside; 293 /** \brief Event involved antikaons in the nucleus at the end of the cascade */ 294 Bool_t antikaonsInside; 295 /** \brief Event involved lambdas in the nucleus at the end of the cascade */ 296 Bool_t lambdasInside; 297 /** \brief Event involved forced delta decays inside the nucleus */ 298 Bool_t forcedDeltasInside; 299 /** \brief Event involved forced delta decays outside the nucleus */ 300 Bool_t forcedDeltasOutside; 301 /** \brief Event involved forced eta/omega decays outside the nucleus */ 302 Bool_t forcedPionResonancesOutside; 303 /** \brief Event involved forced strange absorption inside the nucleus */ 304 Bool_t absorbedStrangeParticle; 305 /** \brief Event involved forced Sigma Zero decays outside the nucleus */ 306 Bool_t forcedSigmaOutside; 307 /** \brief Event involved forced antiKaon/Sigma absorption inside the nucleus */ 308 Bool_t forcedStrangeInside; 309 /** \brief Number of forced Lambda emit out of the nucleus */ 310 Int_t emitLambda; 311 /** \brief Event involved forced Kaon emission */ 312 Bool_t emitKaon; 313 /** \brief Event involved cluster decay */ 314 Bool_t clusterDecay; 315 /** \brief Time of the first collision [fm/c] */ 316 Float_t firstCollisionTime; 317 /** \brief Cross section of the first collision (mb) */ 318 Float_t firstCollisionXSec; 319 /** \brief Position of the spectator on the first collision (fm) */ 320 Float_t firstCollisionSpectatorPosition; 321 /** \brief Momentum of the spectator on the first collision (fm) */ 322 Float_t firstCollisionSpectatorMomentum; 323 /** \brief True if the first collision was elastic */ 324 Bool_t firstCollisionIsElastic; 325 /** \brief Number of reflection avatars */ 326 Int_t nReflectionAvatars; 327 /** \brief Number of collision avatars */ 328 Int_t nCollisionAvatars; 329 /** \brief Number of decay avatars */ 330 Int_t nDecayAvatars; 331 /** \brief Number of dynamical spectators that were merged back into the projectile remnant */ 332 Int_t nUnmergedSpectators; 333 /** \brief Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a solution. */ 334 Int_t nEnergyViolationInteraction; 335 /** \brief Emission time [fm/c] */ 336 Float_t emissionTime[maxSizeParticles]; 337 /** \brief Remnant mass number */ 338 Short_t ARem[maxSizeRemnants]; 339 /** \brief Remnant charge number */ 340 Short_t ZRem[maxSizeRemnants]; 341 /** \brief Remnant strangeness number */ 342 Short_t SRem[maxSizeRemnants]; 343 /** \brief Remnant excitation energy [MeV] */ 344 Float_t EStarRem[maxSizeRemnants]; 345 /** \brief Remnant spin [\f$\hbar\f$] */ 346 Float_t JRem[maxSizeRemnants]; 347 /** \brief Remnant kinetic energy [MeV] */ 348 Float_t EKinRem[maxSizeRemnants]; 349 /** \brief Remnant momentum, x component [MeV/c] */ 350 Float_t pxRem[maxSizeRemnants]; 351 /** \brief Remnant momentum, y component [MeV/c] */ 352 Float_t pyRem[maxSizeRemnants]; 353 /** \brief Remnant momentum, z component [MeV/c] */ 354 Float_t pzRem[maxSizeRemnants]; 355 /** \brief Remnant momentum polar angle [radians] */ 356 Float_t thetaRem[maxSizeRemnants]; 357 /** \brief Remnant momentum azimuthal angle [radians] */ 358 Float_t phiRem[maxSizeRemnants]; 359 /** \brief Remnant angular momentum, x component [\f$\hbar\f$] */ 360 Float_t jxRem[maxSizeRemnants]; 361 /** \brief Remnant angular momentum, y component [\f$\hbar\f$] */ 362 Float_t jyRem[maxSizeRemnants]; 363 /** \brief Remnant angular momentum, z component [\f$\hbar\f$] */ 364 Float_t jzRem[maxSizeRemnants]; 365 /** \brief Particle kinetic energy, in inverse kinematics [MeV] */ 366 Float_t EKinPrime[maxSizeParticles]; 367 /** \brief Particle momentum, z component, in inverse kinematics [MeV/c] */ 368 Float_t pzPrime[maxSizeParticles]; 369 /** \brief Particle momentum polar angle, in inverse kinematics [radians] */ 370 Float_t thetaPrime[maxSizeParticles]; 371 372 /** \brief Reset the EventInfo members */ 373 void reset() { 374 nParticles = 0; 375 event = 0; 376 eventBias = (Float_t)0.0; 377 history.clear(); 378 nRemnants = 0; 379 projectileType = 0; 380 At = 0; 381 Zt = 0; 382 St = 0; 383 Ap = 0; 384 Zp = 0; 385 Sp = 0; 386 Ep = (Float_t)0.0; 387 impactParameter = (Float_t)0.0; 388 nCollisions = 0; 389 stoppingTime = (Float_t)0.0; 390 EBalance = (Float_t)0.0; 391 firstEBalance = (Float_t)0.0; 392 pLongBalance = (Float_t)0.0; 393 pTransBalance = (Float_t)0.0; 394 nCascadeParticles = 0; 395 transparent = false; 396 annihilationP = false; 397 annihilationN = false; 398 forcedCompoundNucleus = false; 399 nucleonAbsorption = false; 400 pionAbsorption = false; 401 nDecays = 0; 402 nSrcCollisions = 0; 403 nSrcPairs = 0; 404 nBlockedCollisions = 0; 405 nBlockedDecays = 0; 406 effectiveImpactParameter = (Float_t)0.0; 407 deltasInside = false; 408 sigmasInside = false; 409 kaonsInside = false; 410 antikaonsInside = false; 411 lambdasInside = false; 412 forcedDeltasInside = false; 413 forcedDeltasOutside = false; 414 forcedPionResonancesOutside = false; 415 absorbedStrangeParticle = false; 416 forcedSigmaOutside = false; 417 forcedStrangeInside = false; 418 emitLambda = 0; 419 emitKaon = false; 420 clusterDecay = false; 421 firstCollisionTime = (Float_t)0.0; 422 firstCollisionXSec = (Float_t)0.0; 423 firstCollisionSpectatorPosition = (Float_t)0.0; 424 firstCollisionSpectatorMomentum = (Float_t)0.0; 425 firstCollisionIsElastic = false; 426 nReflectionAvatars = 0; 427 nCollisionAvatars = 0; 428 nDecayAvatars = 0; 429 nUnmergedSpectators = 0; 430 nEnergyViolationInteraction = 0; 431 432 } 433 434 /// \brief Move a remnant to the particle array 435 void remnantToParticle(const G4int remnantIndex); 436 437 /// \brief Fill the variables describing the reaction in inverse kinematics 438 void fillInverseKinematics(const Double_t gamma); 439 }; 440 } 441 442 #endif /* G4INCLEVENTINFO_HH_HH */ 443