Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BinaryCascade.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//
28// -------------------------------------------------------------------
29// GEANT4 Class file
30//
31//
32// File name: G4BinaryCascade.hh
33//
34// Authors: Alessandro Brunengo (Alessandro.Brunengo@ge.infn.it)
35// Hans-Peter Wellisch
36// Gunter Folger
37//
38// Creation date: 8 June 2000
39//
40// -----------------------------------------------------------------------------
41
42#ifndef G4BinaryCascade_hh
43#define G4BinaryCascade_hh
44
48#include "G4ListOfCollisions.hh"
49#include "G4V3DNucleus.hh"
50#include "G4Fancy3DNucleus.hh"
51#include "G4Fragment.hh"
53#include "G4VScatterer.hh"
54#include "G4LorentzVector.hh"
55#include "G4LorentzRotation.hh"
56
57#include "G4BCDecay.hh"
58#include "G4BCLateParticle.hh"
59#include "G4BCAction.hh"
60
62
63#include "G4Threading.hh"
64
66
67class G4Track;
68class G4KineticTrack;
69class G43DNucleus;
70class G4Scatterer;
71
73{
74public:
75
77
78 virtual ~G4BinaryCascade();
79
81 G4Nucleus& theNucleus);
83 G4V3DNucleus *);
84
85 virtual void ModelDescription(std::ostream&) const;
86 virtual void PropagateModelDescription(std::ostream&) const;
87
88private:
89
90 G4BinaryCascade(const G4BinaryCascade & right);
91 const G4BinaryCascade& operator=(G4BinaryCascade & right);
92 G4bool operator==(G4BinaryCascade& right) {return (this == &right);}
93 G4bool operator!=(G4BinaryCascade& right) {return (this != &right);}
94
95// Implementation
96 void PrintWelcomeMessage();
97 void BuildTargetList();
98 G4bool BuildLateParticleCollisions(G4KineticTrackVector * secondaries);
99 void FindCollisions(G4KineticTrackVector *);
100 void FindDecayCollision(G4KineticTrack *);
101 void FindLateParticleCollision(G4KineticTrack *);
102
103// void PropagateInit();
104// void Cascade();
105 G4ReactionProductVector * DeExcite();
106 G4ReactionProductVector * DecayVoidNucleus();
107 G4ReactionProductVector * ProductsAddFinalState (G4ReactionProductVector * products, G4KineticTrackVector & finalState );
108 G4ReactionProductVector * ProductsAddPrecompound(G4ReactionProductVector * products ,G4ReactionProductVector * preco );
109
110 G4bool ApplyCollision(G4CollisionInitialState *);
111 G4bool Capture(G4bool verbose=false);
112 G4bool Absorb();
113 G4bool CheckPauliPrinciple(G4KineticTrackVector *);
114 G4double GetExcitationEnergy();
115 void CorrectFinalPandE();
116
117 G4double CorrectShortlivedPrimaryForFermi(
118 G4KineticTrack* primary,G4KineticTrackVector target_collection);
119 G4bool CorrectShortlivedFinalsForFermi(
120 G4KineticTrackVector * products, G4double initial_Efermi);
121
122 void UpdateTracksAndCollisions(G4KineticTrackVector * oldSecondaries,
123 G4KineticTrackVector * oldTarget,
124 G4KineticTrackVector * newSecondaries);
125 G4bool DoTimeStep(G4double timeStep);
126 G4KineticTrackVector* CorrectBarionsOnBoundary(G4KineticTrackVector *in,
128 G4Fragment * FindFragments();
129 void StepParticlesOut();
130
131 G4LorentzVector GetFinal4Momentum();
132 G4LorentzVector GetFinalNucleusMomentum();
133
135 G4V3DNucleus *);
136 G4double GetIonMass(G4int Z, G4int A);
137
138 G4int GetTotalCharge(std::vector<G4KineticTrack *> & aV)
139 {
140 G4int result = 0;
141 std::vector<G4KineticTrack *>::iterator i;
142 for(i = aV.begin(); i != aV.end(); ++i)
143 {
144 result += G4lrint((*i)->GetDefinition()->GetPDGCharge());
145 }
146 return result;
147 }
148 G4int GetTotalBaryonCharge(std::vector<G4KineticTrack *> & aV)
149 {
150 G4int result = 0;
151 std::vector<G4KineticTrack *>::iterator i;
152 for(i = aV.begin(); i != aV.end(); ++i)
153 {
154 if ( (*i)->GetDefinition()->GetBaryonNumber() != 0 ){
155 result += G4lrint((*i)->GetDefinition()->GetPDGCharge());
156 }
157 }
158 return result;
159 }
160
161 G4ReactionProductVector * HighEnergyModelFSProducts(G4ReactionProductVector *,
162 G4KineticTrackVector * secondaries); // add secondaries of string model to G4RPV
163 G4ReactionProductVector * FillVoidNucleusProducts(G4ReactionProductVector * ); // nucleus destroyed, add secondaries to G4RPV
164// utility methods
165
166 G4ThreeVector GetSpherePoint(G4double r, const G4LorentzVector & momentumdirection);
167
168 void ClearAndDestroy(G4KineticTrackVector * ktv);
169 void ClearAndDestroy(G4ReactionProductVector * rpv);
170
171// for debugging purpose
172
173 G4ReactionProductVector * ProductsAddFakeGamma(G4ReactionProductVector * products );
174
175 void PrintKTVector(G4KineticTrackVector * ktv, std::string comment=std::string(""));
176 void PrintKTVector(G4KineticTrack* kt, std::string comment=std::string(""));
177 void DebugApplyCollisionFail(G4CollisionInitialState * collision,
178 G4KineticTrackVector * products);
179 void DebugApplyCollision(G4CollisionInitialState * collision,
180 G4KineticTrackVector *products);
181 G4bool DebugFinalEpConservation(const G4HadProjectile & aTrack, G4ReactionProductVector* products);
182 G4bool DebugEpConservation(const G4String where);
183 G4bool CheckChargeAndBaryonNumber(G4String where);
184
185private:
186 G4KineticTrackVector theProjectileList; // replaced by theProjectile4Momentum
187 G4KineticTrackVector theTargetList; // list of nucleons in Nucleus
188 G4KineticTrackVector theSecondaryList; // particles being followed
189 G4KineticTrackVector theCapturedList; // captured particles
190 G4KineticTrackVector theFinalState; // particles for final state
191
192
193 G4ExcitationHandler * theExcitationHandler;
194 G4CollisionManager * theCollisionMgr;
195
196 G4Scatterer * theH1Scatterer;
197
198 std::vector<G4BCAction *> theImR;
199 G4BCDecay * theDecay;
200 G4BCLateParticle * theLateParticle;
201 G4VFieldPropagation * thePropagator;
202 G4DecayKineticTracks decayKTV;
203
204 G4double theCurrentTime;
205 G4double theBCminP;
206 G4double theCutOnP;
207 G4double theCutOnPAbsorb;
208 G4LorentzVector theInitial4Mom; // suppress?
209 G4LorentzVector theProjectile4Momentum;
210 G4int currentA, currentZ, lateA, lateZ;
211 G4int initialZ, initialA,projectileA,projectileZ;
212 G4double massInNucleus, initial_nuclear_mass;
213 G4double currentInitialEnergy; // for debugging
214 G4LorentzRotation precompoundLorentzboost;
215 G4double theOuterRadius;
216 G4bool thePrimaryEscape;
217 const G4ParticleDefinition * thePrimaryType;
218 G4ThreeVector theMomentumTransfer;
219 static G4int theBIC_ID;
220};
221
222#endif
223
224
225
226
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
virtual void PropagateModelDescription(std::ostream &) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *, G4V3DNucleus *)
virtual void ModelDescription(std::ostream &) const
virtual ~G4BinaryCascade()
int G4lrint(double ad)
Definition: templates.hh:134