Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4IntraNucleiCascader.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// 20100315 M. Kelsey -- Remove "using" directory and unnecessary #includes.
28// 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
29// 20100517 M. Kelsey -- Inherit from common base class, make other colliders
30// simple data members
31// 20100617 M. Kelsey -- Make G4NucleiModel a data member, instead of
32// creating and deleting on every cycle.
33// 20100623 M. Kelsey -- Undo change from 0617. G4NucleiModel not reusable.
34// 20100714 M. Kelsey -- Switch to new G4CascadeColliderBase class
35// 20100716 M. Kelsey -- Eliminate inter_case; use base-class functionality,
36// add function to compute recoil nuclear mass on the fly
37// 20100720 M. Kelsey -- Make EPCollider pointer member
38// 20100722 M. Kelsey -- Move cascade output buffers to .hh file
39// 20100728 M. Kelsey -- Move G4NucleiModel here, as pointer member
40// 20100907 M. Kelsey -- Add new "makeResidualFragment" to create
41// G4InuclNuclei at current stage of cascade
42// 20100909 M. Kelsey -- Drop makeResidualFragment(), getResidualMass() and
43// local G4InuclNuclei object, replace with new RecoilMaker.
44// Move goodCase() to RecoilMaker.
45// 20100916 M. Kelsey -- Add functions to handle trapped particles, and to
46// decay hyperons.
47// 20110224 M. Kelsey -- Add ::rescatter() function which takes a list of
48// pre-existing secondaries as input. Split ::collide() into
49// separate utility functions. Move cascade parameters to static
50// data members. Add setVerboseLevel().
51// 20110303 M. Kelsey -- Add more cascade functions to support rescattering
52// 20110304 M. Kelsey -- Modify rescatter to use original Propagate() input
53// 20110316 M. Kelsey -- Add function to do G4KineticTrack conversion, decay
54// rescattering resonances in situ.
55// 20110324 M. Kelsey -- Add list of nucleon hit locations for rescatter().
56// 20110721 M. Kelsey -- Drop decayTrappedParticle(G4KineticTrack*).
57// 20110722 M. Kelsey -- Deprecate "output_particles" list in favor of using
58// output directly (will help with pre-cascade issues).
59// 20110729 M. Kelsey -- Replace convertKineticToCascade() to reduce churn.
60// 20110801 M. Kelsey -- Add local target buffers for rescattering, to avoid
61// memory leak.
62// 20110919 M. Kelsey -- Add optional final-state clustering
63// 20130304 M. Kelsey -- Add new G4CascadeHistory for cacasde structure reporting
64// 20130620 Address Coverity complaint about missing copy actions
65// 20141204 M. Kelsey -- Add function to test for non-interacting particles
66
67#ifndef G4INTRA_NUCLEI_CASCADER_HH
68#define G4INTRA_NUCLEI_CASCADER_HH
69
71#include "G4CollisionOutput.hh"
72#include "G4ThreeVector.hh"
73#include <vector>
74
81class G4InuclParticle;
82class G4KineticTrack;
84class G4NucleiModel;
85class G4V3DNucleus;
86
87
89public:
91 virtual ~G4IntraNucleiCascader();
92
93 void collide(G4InuclParticle* bullet, G4InuclParticle* target,
94 G4CollisionOutput& globalOutput);
95
96 // For use with Propagate to preload a set of secondaries
97 void rescatter(G4InuclParticle* bullet, G4KineticTrackVector* theSecondaries,
98 G4V3DNucleus* theNucleus, G4CollisionOutput& globalOutput);
99
100 void setVerboseLevel(G4int verbose=0);
101
102private:
103 static const G4int itry_max; // Maximum number of attempts
104 static const G4int reflection_cut; // Maximum number of reflections
105 static const G4double small_ekin; // Tolerance for round-off zero
106 static const G4double quasielast_cut; // To recover elastic scatters
107
108protected:
110
111 void newCascade(G4int itry); // Clear buffers for next attempt
112 void setupCascade(); // Fill cascade using nuclear model
113 void generateCascade(); // Track secondaries through nucleus
114 G4bool finishCascade(); // Clean up output, check consistency
115
116 void finalize(G4int itry, // Transfer final state for return
117 G4InuclParticle* bullet, G4InuclParticle* target,
118 G4CollisionOutput& globalOutput);
119
121
122 // Functions to transfer input high-energy cascade for propagation
123 void preloadCascade(G4V3DNucleus* theNucleus,
124 G4KineticTrackVector* theSecondaries);
125 void copyWoundedNucleus(G4V3DNucleus* theNucleus);
126 void copySecondaries(G4KineticTrackVector* theSecondaries);
127 void processSecondary(const G4KineticTrack* aSecondary);
128 void releaseSecondary(const G4KineticTrack* aSecondary);
129
130 // Functions to handle, e.g., low-energy hyperons stuck inside potential
131 void processTrappedParticle(const G4CascadParticle& trapped);
132 void decayTrappedParticle(const G4CascadParticle& trapped);
133
134 // Test if particle is able to interact in nucleus
135 G4bool particleCanInteract(const G4CascadParticle& cpart) const;
136
137private:
138 G4NucleiModel* model;
139 G4ElementaryParticleCollider* theElementaryParticleCollider;
140 G4CascadeRecoilMaker* theRecoilMaker;
141 G4CascadeCoalescence* theClusterMaker;
142 G4CascadeHistory* theCascadeHistory;
143
144 // Buffers and parameters for cascade attempts
145 G4InuclNuclei* tnuclei; // Target nucleus (must be non-zero)
146 G4InuclNuclei* bnuclei; // Non-zero if ion-ion collision
147 G4InuclElementaryParticle* bparticle; // Non-zero if hadron-ion collision
148
149 G4double minimum_recoil_A; // Require fragment with this mass
150 G4double coulombBarrier;
151
152 // Buffers for creation (and reuse) of rescattering targets
153 G4InuclNuclei* nucleusTarget;
154 G4InuclElementaryParticle* protonTarget;
155
156 // Buffers for collecting result of cascade (reset on each iteration)
157 G4CollisionOutput output;
158 std::vector<G4CascadParticle> cascad_particles;
159 std::vector<G4CascadParticle> new_cascad_particles;
160 G4ExitonConfiguration theExitonConfiguration;
161
162 std::vector<G4ThreeVector> hitNucleons; // Nucleons hit before rescatter
163
164private:
165 // Copying of modules is forbidden
168};
169
170#endif /* G4INTRA_NUCLEI_CASCADER_HH */
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void processTrappedParticle(const G4CascadParticle &trapped)
void releaseSecondary(const G4KineticTrack *aSecondary)
void copyWoundedNucleus(G4V3DNucleus *theNucleus)
G4bool particleCanInteract(const G4CascadParticle &cpart) const
void rescatter(G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
void setVerboseLevel(G4int verbose=0)
void decayTrappedParticle(const G4CascadParticle &trapped)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
G4bool initialize(G4InuclParticle *bullet, G4InuclParticle *target)
G4InuclParticle * createTarget(G4V3DNucleus *theNucleus)
void preloadCascade(G4V3DNucleus *theNucleus, G4KineticTrackVector *theSecondaries)
void copySecondaries(G4KineticTrackVector *theSecondaries)
void processSecondary(const G4KineticTrack *aSecondary)
void finalize(G4int itry, G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)