Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
28// 20100315 M. Kelsey -- Remove "using" directive and unnecessary #includes.
29// 20100407 M. Kelsey -- Eliminate return-by-value std::vector<> by creating
30// three data buffers for particles, momenta, and particle types.
31// The following functions now return void and are non-const:
32// ::generateSCMfinalState()
33// ::generateMomModules() (also remove input vector)
34// ::generateStrangeChannelPartTypes()
35// ::generateSCMpionAbsorption()
36// 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide(); merge
37// public vs. private ::collide() functions.
38// 20100511 M. Kelsey -- Remove G4PionSampler and G4NucleonSampler. Expand
39// particle-types selector to all modes, not just strangeness.
40// 20100517 M. Kelsey -- Inherit from common base class, make arrays static
41// 20100714 M. Kelsey -- Switch to new G4CascadeColliderBase class
42// 20100726 M. Kelsey -- Move remaining std::vector<> buffers here
43// 20100804 M. Kelsey -- Add printFinalStateTables() function.
44// 20110923 M. Kelsey -- Add optional stream& to printFinalStateTables().
45// 20130129 M. Kelsey -- Add static arrays and interpolators for two-body
46// angular distributions (addresses MT thread-local issue)
47// 20130131 D. Wright -- Use new *AngDst classes for gamma-N two-body
48// 20130220 M. Kelsey -- Replace two-body angular code with new *AngDst classes
49// 20130221 M. Kelsey -- Move two-body angular dist classes to factory
50// 20130306 M. Kelsey -- Move printFinalStateTables() to table factory
51// 20130307 M. Kelsey -- Reverse order of dimensions for rmn array
52// 20130422 M. Kelsey -- Move kinematics to G4CascadeFinalStateAlgorithm
53// 20130508 D. Wright -- Add muon capture, with absorption on quasideuterons
54// 20130620 Address Coverity complaint about missing copy actions
55// 20130628 Add function to split dibaryon into particle_kinds list
56// 20141009 M. Kelsey -- Add pion absorption by single nucleons, with
57// nuclear recoil. Improves pi- capture performance.
58
59#ifndef G4ELEMENTARY_PARTICLE_COLLIDER_HH
60#define G4ELEMENTARY_PARTICLE_COLLIDER_HH
61
66#include "G4LorentzVector.hh"
67#include <iosfwd>
68#include <vector>
69
71
72
74public:
77
78 void collide(G4InuclParticle* bullet, G4InuclParticle* target,
79 G4CollisionOutput& output);
80
81 void setNucleusState(G4int a, G4int z) { // Nucleus to use for recoil
82 nucleusA = a; nucleusZ = z;
83 }
84
85private:
86 G4int generateMultiplicity(G4int is, G4double ekin) const;
87
88 void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin);
89
90 void generateSCMfinalState(G4double ekin, G4double etot_scm,
92 G4InuclElementaryParticle* particle2);
93
94 // Pion (or photon) absorption on a dibaryon
95 void generateSCMpionAbsorption(G4double etot_scm,
97 G4InuclElementaryParticle* particle2);
98
99 // Muon absorption on a dibaryon (with outgoing neutrino)
100 void generateSCMmuonAbsorption(G4double etot_scm,
101 G4InuclElementaryParticle* particle1,
102 G4InuclElementaryParticle* particle2);
103
104 // Pion absorption on a single nucleon (charge exchange)
105 void generateSCMpionNAbsorption(G4double etot_scm,
106 G4InuclElementaryParticle* particle1,
107 G4InuclElementaryParticle* particle2);
108
109 G4bool pionNucleonAbsorption(G4double ekin) const;
110
111 void fillOutgoingMasses(); // Fill mass arrays from particle types
112
113 // Utility class to generate final-state kinematics
115
116 // Internal buffers for lists of secondaries
117 std::vector<G4InuclElementaryParticle> particles;
118 std::vector<G4LorentzVector> scm_momentums;
119 std::vector<G4double> modules;
120 std::vector<G4double> masses;
121 std::vector<G4double> masses2;
122 std::vector<G4int> particle_kinds;
123
124 // Nuclear environment (to do pion-nucleon absorption)
125 G4int nucleusA, nucleusZ;
126
127private:
128 // Copying of modules is forbidden
131};
132
133#endif /* G4ELEMENTARY_PARTICLE_COLLIDER_HH */
134
135
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)