Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LowEIonFragmentation.cc
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//---------------------------------------------------------------------------
29//
30// ClassName: G4LowEIonFragmentation
31//
32// Author: H.P. Wellisch
33//
34// Modified:
35// 02 Jun 2010 M. A. Cortes Giraldo fix: particlesFromTarget must be
36// accounted for as particles of initial compound nucleus
37// 28 Oct 2010 V.Ivanchenko complete migration to integer Z and A;
38// use updated G4Fragment methods
39
40#include <algorithm>
41
44#include "G4SystemOfUnits.hh"
45#include "G4Fancy3DNucleus.hh"
46#include "G4Proton.hh"
47#include "G4NucleiProperties.hh"
48
49G4int G4LowEIonFragmentation::hits = 0;
50G4int G4LowEIonFragmentation::totalTries = 0;
51G4double G4LowEIonFragmentation::area = 0;
52
54{
55 theHandler = value;
56 theModel = new G4PreCompoundModel(theHandler);
57 proton = G4Proton::Proton();
58}
59
61{
62 theHandler = new G4ExcitationHandler;
63 theModel = new G4PreCompoundModel(theHandler);
64 proton = G4Proton::Proton();
65}
66
68{
69 delete theModel;
70}
71
73ApplyYourself(const G4HadProjectile & thePrimary, G4Nucleus & theNucleus)
74{
75 area = 0;
76 // initialize the particle change
77 theResult.Clear();
78 theResult.SetStatusChange( stopAndKill );
79 theResult.SetEnergyChange( 0.0 );
80
81 // Get Target A, Z
82 G4int aTargetA = theNucleus.GetA_asInt();
83 G4int aTargetZ = theNucleus.GetZ_asInt();
84
85 // Get Projectile A, Z
86 G4int aProjectileA = thePrimary.GetDefinition()->GetBaryonNumber();
87 G4int aProjectileZ = G4lrint(thePrimary.GetDefinition()->GetPDGCharge()/eplus);
88
89 // Get Maximum radius of both
90
91 G4Fancy3DNucleus aPrim;
92 aPrim.Init(aProjectileA, aProjectileZ);
93 G4double projectileOuterRadius = aPrim.GetOuterRadius();
94
95 G4Fancy3DNucleus aTarg;
96 aTarg.Init(aTargetA, aTargetZ);
97 G4double targetOuterRadius = aTarg.GetOuterRadius();
98
99 // Get the Impact parameter
100 G4int particlesFromProjectile = 0;
101 G4int chargedFromProjectile = 0;
102 G4double impactParameter = 0;
103 G4double x,y;
104 G4Nucleon * pNucleon;
105 // need at lease one particle from the projectile model beyond the
106 // projectileHorizon.
107 while(0==particlesFromProjectile)
108 {
109 do
110 {
111 x = 2*G4UniformRand() - 1;
112 y = 2*G4UniformRand() - 1;
113 }
114 while(x*x + y*y > 1);
115 impactParameter = std::sqrt(x*x+y*y)*(targetOuterRadius+projectileOuterRadius);
116 ++totalTries;
117 area = pi*(targetOuterRadius+projectileOuterRadius)*
118 (targetOuterRadius+projectileOuterRadius);
119 G4double projectileHorizon = impactParameter-targetOuterRadius;
120
121 // Empirical boundary transparency.
122 G4double empirical = G4UniformRand();
123 if(projectileHorizon > empirical*projectileOuterRadius) { continue; }
124
125 // Calculate the number of nucleons involved in collision
126 // From projectile
127 aPrim.StartLoop();
128 while((pNucleon = aPrim.GetNextNucleon()))
129 {
130 if(pNucleon->GetPosition().y()>projectileHorizon)
131 {
132 // We have one
133 ++particlesFromProjectile;
134 if(pNucleon->GetParticleType() == proton)
135 {
136 ++chargedFromProjectile;
137 }
138 }
139 }
140 }
141 ++hits;
142
143 // From target:
144 G4double targetHorizon = impactParameter-projectileOuterRadius;
145 G4int chargedFromTarget = 0;
146 G4int particlesFromTarget = 0;
147 aTarg.StartLoop();
148 while((pNucleon = aTarg.GetNextNucleon()))
149 {
150 if(pNucleon->GetPosition().y()>targetHorizon)
151 {
152 // We have one
153 ++particlesFromTarget;
154 if(pNucleon->GetParticleType() == proton)
155 {
156 ++chargedFromTarget;
157 }
158 }
159 }
160
161 // Energy sharing between projectile and target.
162 // Note that this is a quite simplistic kinetically.
163 G4ThreeVector momentum = thePrimary.Get4Momentum().vect();
164 G4double w = (G4double)particlesFromProjectile/(G4double)aProjectileA;
165
166 G4double projTotEnergy = thePrimary.GetTotalEnergy();
167 G4double targetMass = G4NucleiProperties::GetNuclearMass(aTargetA, aTargetZ);
168 G4LorentzVector fragment4Momentum(momentum*w, projTotEnergy*w + targetMass);
169
170 // take the nucleons and fill the Fragments
171 G4Fragment anInitialState(aTargetA+particlesFromProjectile,
172 aTargetZ+chargedFromProjectile,
173 fragment4Momentum);
174 // M.A. Cortes fix
175 //anInitialState.SetNumberOfParticles(particlesFromProjectile);
176 anInitialState.SetNumberOfExcitedParticle(particlesFromProjectile+particlesFromTarget,
177 chargedFromProjectile + chargedFromTarget);
178 anInitialState.SetNumberOfHoles(particlesFromProjectile+particlesFromTarget,
179 chargedFromProjectile + chargedFromTarget);
180 G4double time = thePrimary.GetGlobalTime();
181 anInitialState.SetCreationTime(time);
182
183 // Fragment the Fragment using Pre-compound
184 G4ReactionProductVector* thePreCompoundResult = theModel->DeExcite(anInitialState);
185
186 // De-excite the projectile using ExcitationHandler
187 G4ReactionProductVector * theExcitationResult = 0;
188 if(particlesFromProjectile < aProjectileA)
189 {
190 G4LorentzVector residual4Momentum(momentum*(1.0-w), projTotEnergy*(1.0-w));
191
192 G4Fragment initialState2(aProjectileA-particlesFromProjectile,
193 aProjectileZ-chargedFromProjectile,
194 residual4Momentum );
195
196 // half of particles are excited (?!)
197 G4int pinit = (aProjectileA-particlesFromProjectile)/2;
198 G4int cinit = (aProjectileZ-chargedFromProjectile)/2;
199
200 initialState2.SetNumberOfExcitedParticle(pinit,cinit);
201 initialState2.SetNumberOfHoles(pinit,cinit);
202 initialState2.SetCreationTime(time);
203
204 theExcitationResult = theHandler->BreakItUp(initialState2);
205 }
206
207 // Fill the particle change and clear intermediate vectors
208 G4int nexc = 0;
209 G4int npre = 0;
210 if(theExcitationResult) { nexc = theExcitationResult->size(); }
211 if(thePreCompoundResult) { npre = thePreCompoundResult->size();}
212
213 if(nexc > 0) {
214 for(G4int k=0; k<nexc; ++k) {
215 G4ReactionProduct* p = (*theExcitationResult)[k];
217 delete p;
218 }
219 }
220
221 if(npre > 0) {
222 for(G4int k=0; k<npre; ++k) {
223 G4ReactionProduct* p = (*thePreCompoundResult)[k];
225 delete p;
226 }
227 }
228
229 delete thePreCompoundResult;
230 delete theExcitationResult;
231
232 // return the particle change
233 return &theResult;
234
235}
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
double y() const
Hep3Vector vect() const
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState) const
G4Nucleon * GetNextNucleon()
G4double GetOuterRadius()
void Init(G4int theA, G4int theZ)
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:383
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:335
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:316
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
G4double GetTotalEnergy() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ParticleDefinition * GetParticleType() const
Definition: G4Nucleon.hh:84
virtual const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:68
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetPDGCharge() const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
int G4lrint(double ad)
Definition: templates.hh:163