Geant4  10.02.p02
G4ElementaryParticleCollider.hh
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 // $Id: G4ElementaryParticleCollider.hh 71940 2013-06-28 19:04:44Z mkelsey $
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 // 20100714 M. Kelsey -- Switch to new G4CascadeColliderBase class
43 // 20100726 M. Kelsey -- Move remaining std::vector<> buffers here
44 // 20100804 M. Kelsey -- Add printFinalStateTables() function.
45 // 20110923 M. Kelsey -- Add optional stream& to printFinalStateTables().
46 // 20130129 M. Kelsey -- Add static arrays and interpolators for two-body
47 // angular distributions (addresses MT thread-local issue)
48 // 20130131 D. Wright -- Use new *AngDst classes for gamma-N two-body
49 // 20130220 M. Kelsey -- Replace two-body angular code with new *AngDst classes
50 // 20130221 M. Kelsey -- Move two-body angular dist classes to factory
51 // 20130306 M. Kelsey -- Move printFinalStateTables() to table factory
52 // 20130307 M. Kelsey -- Reverse order of dimensions for rmn array
53 // 20130422 M. Kelsey -- Move kinematics to G4CascadeFinalStateAlgorithm
54 // 20130508 D. Wright -- Add muon capture, with absorption on quasideuterons
55 // 20130620 Address Coverity complaint about missing copy actions
56 // 20130628 Add function to split dibaryon into particle_kinds list
57 // 20141009 M. Kelsey -- Add pion absorption by single nucleons, with
58 // nuclear recoil. Improves pi- capture performance.
59 
60 #ifndef G4ELEMENTARY_PARTICLE_COLLIDER_HH
61 #define G4ELEMENTARY_PARTICLE_COLLIDER_HH
62 
63 #include "G4CascadeColliderBase.hh"
65 #include "G4CascadeInterpolator.hh"
67 #include "G4LorentzVector.hh"
68 #include <iosfwd>
69 #include <vector>
70 
71 class G4CollisionOutput;
72 
73 
75 public:
78 
79  void collide(G4InuclParticle* bullet, G4InuclParticle* target,
80  G4CollisionOutput& output);
81 
82  void setNucleusState(G4int a, G4int z) { // Nucleus to use for recoil
83  nucleusA = a; nucleusZ = z;
84  }
85 
86 private:
87  G4int generateMultiplicity(G4int is, G4double ekin) const;
88 
89  void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin);
90 
91  void generateSCMfinalState(G4double ekin, G4double etot_scm,
92  G4InuclElementaryParticle* particle1,
93  G4InuclElementaryParticle* particle2);
94 
95  // Pion (or photon) absorption on a dibaryon
96  void generateSCMpionAbsorption(G4double etot_scm,
97  G4InuclElementaryParticle* particle1,
98  G4InuclElementaryParticle* particle2);
99 
100  // Muon absorption on a dibaryon (with outgoing neutrino)
101  void generateSCMmuonAbsorption(G4double etot_scm,
102  G4InuclElementaryParticle* particle1,
103  G4InuclElementaryParticle* particle2);
104 
105  // Pion absorption on a single nucleon (charge exchange)
106  void generateSCMpionNAbsorption(G4double etot_scm,
107  G4InuclElementaryParticle* particle1,
108  G4InuclElementaryParticle* particle2);
109 
111 
112  G4bool splitQuasiDeuteron(G4int qdtype); // Fill kinds with NN components
113 
114  void fillOutgoingMasses(); // Fill mass arrays from particle types
115 
116  // Utility class to generate final-state kinematics
118 
119  // Internal buffers for lists of secondaries
120  std::vector<G4InuclElementaryParticle> particles;
121  std::vector<G4LorentzVector> scm_momentums;
122  std::vector<G4double> modules;
123  std::vector<G4double> masses;
124  std::vector<G4double> masses2;
125  std::vector<G4int> particle_kinds;
126 
127  // Nuclear environment (to do pion-nucleon absorption)
129 
130 private:
131  // Copying of modules is forbidden
134 };
135 
136 #endif /* G4ELEMENTARY_PARTICLE_COLLIDER_HH */
137 
138 
void generateSCMpionNAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
void generateSCMpionAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
G4CascadeFinalStateGenerator fsGenerator
G4double z
Definition: TRTMaterials.hh:39
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
void generateSCMmuonAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
std::vector< G4LorentzVector > scm_momentums
void generateSCMfinalState(G4double ekin, G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
std::vector< G4InuclElementaryParticle > particles
bool G4bool
Definition: G4Types.hh:79
void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin)
G4int generateMultiplicity(G4int is, G4double ekin) const
G4ElementaryParticleCollider & operator=(const G4ElementaryParticleCollider &)
double G4double
Definition: G4Types.hh:76
G4bool pionNucleonAbsorption(G4double ekin) const