Changeset 1315 for trunk/source/processes/hadronic/models/cascade/cascade/include/G4ElementaryParticleCollider.hh
- Timestamp:
- Jun 18, 2010, 11:42:07 AM (14 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/source/processes/hadronic/models/cascade/cascade/include/G4ElementaryParticleCollider.hh
r1196 r1315 23 23 // * acceptance of all terms of the Geant4 Software license. * 24 24 // ******************************************************************** 25 // $Id: G4ElementaryParticleCollider.hh,v 1.31 2010/05/21 17:56:34 mkelsey Exp $ 26 // Geant4 tag: $Name: geant4-09-04-beta-cand-01 $ 25 27 // 28 // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly 29 // 20100315 M. Kelsey -- Remove "using" directive and unnecessary #includes. 30 // 20100407 M. Kelsey -- Eliminate return-by-value std::vector<> by creating 31 // three data buffers for particles, momenta, and particle types. 32 // The following functions now return void and are non-const: 33 // ::generateSCMfinalState() 34 // ::generateMomModules() (also remove input vector) 35 // ::generateStrangeChannelPartTypes() 36 // ::generateSCMpionAbsorption() 37 // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide(); merge 38 // public vs. private ::collide() functions. 39 // 20100511 M. Kelsey -- Remove G4PionSampler and G4NucleonSampler. Expand 40 // particle-types selector to all modes, not just strangeness. 41 // 20100517 M. Kelsey -- Inherit from common base class, make arrays static 42 26 43 #ifndef G4ELEMENTARY_PARTICLE_COLLIDER_HH 27 44 #define G4ELEMENTARY_PARTICLE_COLLIDER_HH 28 45 29 #include "G4 Collider.hh"46 #include "G4VCascadeCollider.hh" 30 47 #include "G4InuclElementaryParticle.hh" 31 #include "G4InuclSpecialFunctions.hh" 32 #include "G4CascadSpecialFunctions.hh" 33 #include "G4LorentzConvertor.hh" 34 #include "G4NucleonSampler.hh" 35 #include "G4PionSampler.hh" 48 #include "G4LorentzVector.hh" 49 #include <vector> 36 50 37 38 using namespace G4InuclSpecialFunctions; 39 using namespace G4CascadSpecialFunctions; 51 class G4LorentzConvertor; 52 class G4CollisionOutput; 40 53 41 54 42 class G4ElementaryParticleCollider { 43 55 class G4ElementaryParticleCollider : public G4VCascadeCollider { 44 56 public: 45 46 57 G4ElementaryParticleCollider(); 47 48 G4CollisionOutput collide(G4InuclParticle* bullet, 49 G4InuclParticle* target); 58 virtual ~G4ElementaryParticleCollider() {}; 59 60 void collide(G4InuclParticle* bullet, G4InuclParticle* target, 61 G4CollisionOutput& output); 50 62 51 63 private: 52 53 G4int verboseLevel;54 55 64 void initializeArrays(); 56 65 57 66 G4int generateMultiplicity(G4int is, G4double ekin) const; 58 67 59 void collide(G4InuclElementaryParticle* bullet, 60 G4InuclElementaryParticle* target, 61 G4CollisionOutput& output); 68 void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin); 62 69 63 64 std::vector<G4InuclElementaryParticle> 65 generateSCMfinalState(G4double ekin, G4double etot_scm, G4double pscm, 66 G4InuclElementaryParticle* particle1, 67 G4InuclElementaryParticle* particle2, 68 G4LorentzConvertor* toSCM) const; 70 void generateSCMfinalState(G4double ekin, G4double etot_scm, G4double pscm, 71 G4InuclElementaryParticle* particle1, 72 G4InuclElementaryParticle* particle2, 73 G4LorentzConvertor* toSCM); 74 75 void generateSCMpionAbsorption(G4double etot_scm, 76 G4InuclElementaryParticle* particle1, 77 G4InuclElementaryParticle* particle2); 78 79 void generateMomModules(G4int mult, G4int is, G4double ekin, 80 G4double etot_cm); 81 82 // Samples the CM momentum for elastic and charge exchange scattering 83 // 84 G4LorentzVector 85 sampleCMmomentumFor2to2(G4int is, G4int kw, G4double ekin, 86 G4double pscm) const; 69 87 70 88 71 std::vector<G4double> 72 generateMomModules(const std::vector<G4int>& kinds, G4int mult, 73 G4int is, G4double ekin, G4double etot_cm) const; 89 // Samples cos(theta) in the CM for elastic and charge exchange scattering 90 // 91 G4double sampleCMcosFor2to2(G4double pscm, G4double pFrac, 92 G4double pA, G4double pC, G4double pCos) const; 74 93 75 94 76 G4CascadeMomentum 77 particleSCMmomentumFor2to2(G4int is, G4int kw, G4double ekin, 78 G4double pscm) const; 79 80 81 G4int getElasticCase(G4int is, G4int kw, G4double ekin) const; 82 83 84 std::vector<G4int> 85 generateStrangeChannelPartTypes(G4int is, G4int mult, 86 G4double ekin) const; 87 88 89 G4double 90 getMomModuleFor2toMany(G4int is, G4int mult, G4int knd, 91 G4double ekin) const; 95 G4double getMomModuleFor2toMany(G4int is, G4int mult, G4int knd, 96 G4double ekin) const; 92 97 93 98 94 99 G4bool satisfyTriangle(const std::vector<G4double>& modules) const; 95 100 96 G4 CascadeMomentum101 G4LorentzVector 97 102 particleSCMmomentumFor2to3(G4int is, G4int knd, G4double ekin, 98 103 G4double pmod) const; 99 104 100 101 std::pair<G4double, G4double> 102 adjustIntervalForElastic(G4double ekin, G4double ak, G4double ae, 103 G4int k, G4int l, const std::vector<G4double>& ssv, 104 G4double st) const; 105 106 std::vector<G4InuclElementaryParticle> 107 generateSCMpionAbsorption(G4double etot_scm, 108 G4InuclElementaryParticle* particle1, 109 G4InuclElementaryParticle* particle2) const; 110 111 G4NucleonSampler nucSampler; 112 G4PionSampler piSampler; 105 // Internal buffers for lists of secondaries 106 std::vector<G4InuclElementaryParticle> particles; 107 std::vector<G4double> modules; 108 std::vector<G4int> particle_kinds; 113 109 114 110 // Parameter arrays 115 116 G4double rmn[14][10][2]; 117 G4double ang[4][4][13]; 118 G4double abn[4][4][4]; 119 111 static const G4double rmn[14][10][2]; 112 static const G4double abn[4][4][4]; 120 113 }; 121 114 122 #endif // G4ELEMENTARY_PARTICLE_COLLIDER_HH115 #endif /* G4ELEMENTARY_PARTICLE_COLLIDER_HH */ 123 116 124 117
Note: See TracChangeset
for help on using the changeset viewer.