source: trunk/source/processes/hadronic/models/pre_equilibrium/exciton_model/src/G4LowEIonFragmentation.cc @ 1228

Last change on this file since 1228 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 7.7 KB
Line 
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#include "G4LowEIonFragmentation.hh"
27#include <algorithm>
28
29G4int G4LowEIonFragmentation::hits = 0;
30G4int G4LowEIonFragmentation::totalTries = 0;
31G4double G4LowEIonFragmentation::area = 0;
32
33G4HadFinalState * G4LowEIonFragmentation::
34ApplyYourself(const G4HadProjectile & thePrimary, G4Nucleus & theNucleus)
35{
36  area = 0;
37  // initialize the particle change
38  theResult.Clear();
39  theResult.SetStatusChange( stopAndKill );
40  theResult.SetEnergyChange( 0.0 );
41
42  // Get Target A, Z
43  G4double aTargetA = theNucleus.GetN();
44  G4double aTargetZ = theNucleus.GetZ();
45
46  // Get Projectile A, Z
47  G4double aProjectileA = thePrimary.GetDefinition()->GetBaryonNumber();
48  G4double aProjectileZ = thePrimary.GetDefinition()->GetPDGCharge();
49
50  // Get Maximum radius of both
51 
52  G4Fancy3DNucleus aPrim;
53  aPrim.Init(aProjectileA, aProjectileZ);
54  G4double projectileOuterRadius = aPrim.GetOuterRadius();
55 
56  G4Fancy3DNucleus aTarg;
57  aTarg.Init(aTargetA, aTargetZ);
58  G4double targetOuterRadius = aTarg.GetOuterRadius();
59
60  // Get the Impact parameter
61  G4int particlesFromProjectile = 0;
62  G4int chargedFromProjectile = 0;
63  G4double impactParameter = 0;
64  G4double x,y;
65  G4Nucleon * pNucleon;
66  // need at lease one particle from the projectile model beyond the
67  // projectileHorizon.
68  while(0==particlesFromProjectile)
69  {
70    do
71    {
72      x = 2*G4UniformRand() - 1;
73      y = 2*G4UniformRand() - 1;
74    }
75    while(x*x + y*y > 1);
76    impactParameter = std::sqrt(x*x+y*y)*(targetOuterRadius+projectileOuterRadius);
77    totalTries++;
78    area = pi*(targetOuterRadius+projectileOuterRadius)*
79              (targetOuterRadius+projectileOuterRadius);
80    G4double projectileHorizon = impactParameter-targetOuterRadius; 
81   
82    // Empirical boundary transparency.
83    G4double empirical = G4UniformRand();
84    if(projectileHorizon/projectileOuterRadius>empirical) continue;
85   
86    // Calculate the number of nucleons involved in collision
87    // From projectile
88    aPrim.StartLoop();
89    while((pNucleon = aPrim.GetNextNucleon()))
90    {
91      if(pNucleon->GetPosition().y()>projectileHorizon)
92      {
93        // We have one
94        particlesFromProjectile++;
95        if(pNucleon->GetParticleType()==G4Proton::ProtonDefinition()) 
96        {
97          chargedFromProjectile++;
98        } 
99      }
100    }
101  }
102  hits ++;
103
104  // From target:
105  G4double targetHorizon = impactParameter-projectileOuterRadius;
106  G4int chargedFromTarget = 0;
107  G4int particlesFromTarget = 0;
108  aTarg.StartLoop(); 
109  while((pNucleon = aTarg.GetNextNucleon()))
110  {
111    if(pNucleon->GetPosition().y()>targetHorizon)
112    {
113      // We have one
114      particlesFromTarget++;
115      if(pNucleon->GetParticleType()==G4Proton::ProtonDefinition()) 
116      {
117        chargedFromTarget++;
118      }
119    }
120  }
121 
122  // Energy sharing between projectile and target. Note that this is a quite simplistic kinetically.
123  G4ThreeVector exciton3Momentum = thePrimary.Get4Momentum().vect();
124  exciton3Momentum *= particlesFromProjectile/aProjectileA;
125 
126  G4double compoundEnergy = thePrimary.GetTotalEnergy()*particlesFromProjectile/aProjectileA; 
127  G4double targetMass = G4ParticleTable::GetParticleTable()
128                        ->GetIonTable()->GetIonMass(static_cast<G4int>(aTargetZ) ,static_cast<G4int>(aTargetA));
129  compoundEnergy += targetMass;
130  G4LorentzVector fragment4Momentum(exciton3Momentum, compoundEnergy);
131 
132  // take the nucleons and fill the Fragments
133  G4Fragment anInitialState;
134  anInitialState.SetA(aTargetA+particlesFromProjectile);
135  anInitialState.SetZ(aTargetZ+chargedFromProjectile);
136  anInitialState.SetNumberOfParticles(particlesFromProjectile);
137  anInitialState.SetNumberOfHoles(particlesFromTarget);
138  anInitialState.SetNumberOfCharged(chargedFromProjectile + chargedFromTarget);
139  anInitialState.SetMomentum(fragment4Momentum);
140
141  // Fragment the Fragment using Pre-compound
142  G4ReactionProductVector* thePreCompoundResult;
143  thePreCompoundResult = theModel->DeExcite(anInitialState);
144 
145  // De-excite the projectile using ExcitationHandler
146 
147  G4ReactionProductVector * theExcitationResult = 0; 
148  if(particlesFromProjectile != aProjectileA)
149  {
150    G4ThreeVector residual3Momentum = thePrimary.Get4Momentum().vect();
151    residual3Momentum -= exciton3Momentum;
152    G4double residualEnergy = thePrimary.GetTotalEnergy()*(1.-particlesFromProjectile/aProjectileA);
153    G4LorentzVector residual4Momentum(residual3Momentum, residualEnergy); 
154 
155    G4Fragment initialState2;
156    initialState2.SetA(aProjectileA-particlesFromProjectile);
157    initialState2.SetZ(aProjectileZ-chargedFromProjectile);
158    initialState2.SetNumberOfHoles(static_cast<G4int>((aProjectileA-particlesFromProjectile)/2.0));
159    initialState2.SetNumberOfParticles(static_cast<G4int>((aProjectileZ-chargedFromProjectile)/2.0));
160    initialState2.SetNumberOfCharged(static_cast<G4int>((aProjectileZ-chargedFromProjectile)/2.0));
161
162
163    initialState2.SetMomentum(residual4Momentum);
164    theExcitationResult = theHandler->BreakItUp(initialState2);
165  }
166
167  // Fill the particle change
168  G4int nSecondaries = 0;
169  if(theExcitationResult) nSecondaries+=theExcitationResult->size();
170  if(thePreCompoundResult) nSecondaries+=thePreCompoundResult->size();
171 
172  unsigned int k;
173  if(theExcitationResult!=0)
174  {
175    for(k=0; k<theExcitationResult->size(); k++)
176    {
177      G4DynamicParticle* p0 = new G4DynamicParticle;
178      p0->SetDefinition( theExcitationResult->operator[](k)->GetDefinition() );
179      p0->SetMomentum( theExcitationResult->operator[](k)->GetMomentum() );
180      theResult.AddSecondary(p0);
181    }
182  }
183 
184  for(k=0; k<thePreCompoundResult->size(); k++)
185  {
186    G4DynamicParticle* p0 = new G4DynamicParticle;
187    p0->SetDefinition(thePreCompoundResult->operator[](k)->GetDefinition());
188    p0->SetMomentum(thePreCompoundResult->operator[](k)->GetMomentum());
189    theResult.AddSecondary(p0);
190  }
191 
192  // clean up
193  std::for_each(thePreCompoundResult->begin(), thePreCompoundResult->end(), DeleteReactionProduct());
194  if(theExcitationResult) 
195  {
196    std::for_each(theExcitationResult->begin(), theExcitationResult->end(), DeleteReactionProduct());
197  }
198  delete thePreCompoundResult;
199  if(theExcitationResult) delete theExcitationResult;
200
201  // return the particle change
202  return &theResult;
203 
204}
Note: See TracBrowser for help on using the repository browser.