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