|
||||
File indexing completed on 2025-01-18 09:58:11
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly 0028 // 20100315 M. Kelsey -- Remove "using" directive and unnecessary #includes. 0029 // 20100407 M. Kelsey -- Eliminate return-by-value std::vector<> by creating 0030 // three data buffers for particles, momenta, and particle types. 0031 // The following functions now return void and are non-const: 0032 // ::generateSCMfinalState() 0033 // ::generateMomModules() (also remove input vector) 0034 // ::generateStrangeChannelPartTypes() 0035 // ::generateSCMpionAbsorption() 0036 // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide(); merge 0037 // public vs. private ::collide() functions. 0038 // 20100511 M. Kelsey -- Remove G4PionSampler and G4NucleonSampler. Expand 0039 // particle-types selector to all modes, not just strangeness. 0040 // 20100517 M. Kelsey -- Inherit from common base class, make arrays static 0041 // 20100714 M. Kelsey -- Switch to new G4CascadeColliderBase class 0042 // 20100726 M. Kelsey -- Move remaining std::vector<> buffers here 0043 // 20100804 M. Kelsey -- Add printFinalStateTables() function. 0044 // 20110923 M. Kelsey -- Add optional stream& to printFinalStateTables(). 0045 // 20130129 M. Kelsey -- Add static arrays and interpolators for two-body 0046 // angular distributions (addresses MT thread-local issue) 0047 // 20130131 D. Wright -- Use new *AngDst classes for gamma-N two-body 0048 // 20130220 M. Kelsey -- Replace two-body angular code with new *AngDst classes 0049 // 20130221 M. Kelsey -- Move two-body angular dist classes to factory 0050 // 20130306 M. Kelsey -- Move printFinalStateTables() to table factory 0051 // 20130307 M. Kelsey -- Reverse order of dimensions for rmn array 0052 // 20130422 M. Kelsey -- Move kinematics to G4CascadeFinalStateAlgorithm 0053 // 20130508 D. Wright -- Add muon capture, with absorption on quasideuterons 0054 // 20130620 Address Coverity complaint about missing copy actions 0055 // 20130628 Add function to split dibaryon into particle_kinds list 0056 // 20141009 M. Kelsey -- Add pion absorption by single nucleons, with 0057 // nuclear recoil. Improves pi- capture performance. 0058 0059 #ifndef G4ELEMENTARY_PARTICLE_COLLIDER_HH 0060 #define G4ELEMENTARY_PARTICLE_COLLIDER_HH 0061 0062 #include "G4CascadeColliderBase.hh" 0063 #include "G4CascadeFinalStateGenerator.hh" 0064 #include "G4CascadeInterpolator.hh" 0065 #include "G4InuclElementaryParticle.hh" 0066 #include "G4LorentzVector.hh" 0067 #include <iosfwd> 0068 #include <vector> 0069 0070 class G4CollisionOutput; 0071 0072 0073 class G4ElementaryParticleCollider : public G4CascadeColliderBase { 0074 public: 0075 G4ElementaryParticleCollider(); 0076 virtual ~G4ElementaryParticleCollider() {}; 0077 0078 void collide(G4InuclParticle* bullet, G4InuclParticle* target, 0079 G4CollisionOutput& output); 0080 0081 void setNucleusState(G4int a, G4int z) { // Nucleus to use for recoil 0082 nucleusA = a; nucleusZ = z; 0083 } 0084 0085 private: 0086 G4int generateMultiplicity(G4int is, G4double ekin) const; 0087 0088 void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin); 0089 0090 void generateSCMfinalState(G4double ekin, G4double etot_scm, 0091 G4InuclElementaryParticle* particle1, 0092 G4InuclElementaryParticle* particle2); 0093 0094 // Pion (or photon) absorption on a dibaryon 0095 void generateSCMpionAbsorption(G4double etot_scm, 0096 G4InuclElementaryParticle* particle1, 0097 G4InuclElementaryParticle* particle2); 0098 0099 // Muon absorption on a dibaryon (with outgoing neutrino) 0100 void generateSCMmuonAbsorption(G4double etot_scm, 0101 G4InuclElementaryParticle* particle1, 0102 G4InuclElementaryParticle* particle2); 0103 0104 // Pion absorption on a single nucleon (charge exchange) 0105 void generateSCMpionNAbsorption(G4double etot_scm, 0106 G4InuclElementaryParticle* particle1, 0107 G4InuclElementaryParticle* particle2); 0108 0109 G4bool pionNucleonAbsorption(G4double ekin) const; 0110 0111 void fillOutgoingMasses(); // Fill mass arrays from particle types 0112 0113 // Utility class to generate final-state kinematics 0114 G4CascadeFinalStateGenerator fsGenerator; 0115 0116 // Internal buffers for lists of secondaries 0117 std::vector<G4InuclElementaryParticle> particles; 0118 std::vector<G4LorentzVector> scm_momentums; 0119 std::vector<G4double> modules; 0120 std::vector<G4double> masses; 0121 std::vector<G4double> masses2; 0122 std::vector<G4int> particle_kinds; 0123 0124 // Nuclear environment (to do pion-nucleon absorption) 0125 G4int nucleusA, nucleusZ; 0126 0127 private: 0128 // Copying of modules is forbidden 0129 G4ElementaryParticleCollider(const G4ElementaryParticleCollider&); 0130 G4ElementaryParticleCollider& operator=(const G4ElementaryParticleCollider&); 0131 }; 0132 0133 #endif /* G4ELEMENTARY_PARTICLE_COLLIDER_HH */ 0134 0135
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |