source: trunk/source/processes/hadronic/models/abrasion/src/G4WilsonAbrasionModel.cc @ 1348

Last change on this file since 1348 was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 32.5 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// *                                                                  *
21// * Parts of this code which have been  developed by QinetiQ Ltd     *
22// * under contract to the European Space Agency (ESA) are the        *
23// * intellectual property of ESA. Rights to use, copy, modify and    *
24// * redistribute this software for general public use are granted    *
25// * in compliance with any licensing, distribution and development   *
26// * policy adopted by the Geant4 Collaboration. This code has been   *
27// * written by QinetiQ Ltd for the European Space Agency, under ESA  *
28// * contract 17191/03/NL/LvH (Aurora Programme).                     *
29// *                                                                  *
30// * By using,  copying,  modifying or  distributing the software (or *
31// * any work based  on the software)  you  agree  to acknowledge its *
32// * use  in  resulting  scientific  publications,  and indicate your *
33// * acceptance of all terms of the Geant4 Software license.          *
34// ********************************************************************
35//
36// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
37//
38// MODULE:              G4WilsonAbrasionModel.cc
39//
40// Version:             1.0
41// Date:                08/12/2009
42// Author:              P R Truscott
43// Organisation:        QinetiQ Ltd, UK
44// Customer:            ESA/ESTEC, NOORDWIJK
45// Contract:            17191/03/NL/LvH
46//
47// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
48//
49// CHANGE HISTORY
50// --------------
51//
52// 6 October 2003, P R Truscott, QinetiQ Ltd, UK
53// Created.
54//
55// 15 March 2004, P R Truscott, QinetiQ Ltd, UK
56// Beta release
57//
58// 18 January 2005, M H Mendenhall, Vanderbilt University, US
59// Pointers to theAbrasionGeometry and products generated by the deexcitation
60// handler deleted to prevent memory leaks.  Also particle change of projectile
61// fragment previously not properly defined.
62//
63// 08 December 2009, P R Truscott, QinetiQ Ltd, Ltd
64// ver 1.0
65// There was originally a possibility of the minimum impact parameter AFTER
66// considering Couloumb repulsion, rm, being too large.  Now if:
67//     rm >= fradius * (rP + rT)
68// where fradius is currently 0.99, then it is assumed the primary track is
69// unchanged.  Additional conditions to escape from while-loop over impact
70// parameter: if the loop counter evtcnt reaches 1000, the primary track is
71// continued.
72// Additional clauses have been included in
73//    G4WilsonAbrasionModel::GetNucleonInducedExcitation
74// Previously it was possible to get sqrt of negative number as Wilson
75// algorithm not properly defined if either:
76//    rT > rP && rsq < rTsq - rPsq) or (rP > rT && rsq < rPsq - rTsq)
77//
78// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
79///////////////////////////////////////////////////////////////////////////////
80
81#include "G4WilsonAbrasionModel.hh"
82#include "G4WilsonRadius.hh"
83#include "G4NuclearAbrasionGeometry.hh"
84#include "G4WilsonAblationModel.hh"
85
86#include "G4ExcitationHandler.hh"
87#include "G4Evaporation.hh"
88#include "G4FermiBreakUp.hh"
89#include "G4StatMF.hh"
90#include "G4ParticleDefinition.hh"
91#include "G4DynamicParticle.hh"
92#include "Randomize.hh"
93#include "G4Fragment.hh"
94#include "G4ReactionProductVector.hh"
95#include "G4LorentzVector.hh"
96#include "G4ParticleMomentum.hh"
97#include "G4Poisson.hh"
98#include "G4ParticleTable.hh"
99#include "G4IonTable.hh"
100#include "globals.hh"
101
102
103G4WilsonAbrasionModel::G4WilsonAbrasionModel (G4bool useAblation1)
104  :G4HadronicInteraction("G4WilsonAbrasion")
105{
106  // Send message to stdout to advise that the G4Abrasion model is being used.
107  PrintWelcomeMessage();
108
109  // Set the default verbose level to 0 - no output.
110  verboseLevel = 0;
111  useAblation  = useAblation1;
112
113  // No de-excitation handler has been supplied - define the default handler.
114
115  theExcitationHandler  = new G4ExcitationHandler;
116  theExcitationHandlerx = new G4ExcitationHandler;
117  if (useAblation)
118  {
119    theAblation = new G4WilsonAblationModel;
120    theAblation->SetVerboseLevel(verboseLevel);
121    theExcitationHandler->SetEvaporation(theAblation);
122    theExcitationHandlerx->SetEvaporation(theAblation);
123  }
124  else
125  {
126    theAblation                      = NULL;
127    G4Evaporation * theEvaporation   = new G4Evaporation;
128    G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp;
129    G4StatMF * theMF                 = new G4StatMF;
130    theExcitationHandler->SetEvaporation(theEvaporation);
131    theExcitationHandler->SetFermiModel(theFermiBreakUp);
132    theExcitationHandler->SetMultiFragmentation(theMF);
133    theExcitationHandler->SetMaxAandZForFermiBreakUp(12, 6);
134    theExcitationHandler->SetMinEForMultiFrag(5.0*MeV);
135
136    theEvaporation  = new G4Evaporation;
137    theFermiBreakUp = new G4FermiBreakUp;
138    theExcitationHandlerx->SetEvaporation(theEvaporation);
139    theExcitationHandlerx->SetFermiModel(theFermiBreakUp);
140    theExcitationHandlerx->SetMaxAandZForFermiBreakUp(12, 6);
141  }
142
143  // Set the minimum and maximum range for the model (despite nomanclature,
144  // this is in energy per nucleon number). 
145
146  SetMinEnergy(70.0*MeV);
147  SetMaxEnergy(10.1*GeV);
148  isBlocked = false;
149
150  // npK, when mutiplied by the nuclear Fermi momentum, determines the range of
151  // momentum over which the secondary nucleon momentum is sampled.
152
153  r0sq = 0.0;
154  npK = 5.0;
155  B = 10.0 * MeV;
156  third = 1.0 / 3.0;
157  fradius = 0.99;
158  conserveEnergy = false;
159  conserveMomentum = true;
160}
161
162
163G4WilsonAbrasionModel::G4WilsonAbrasionModel (G4ExcitationHandler *aExcitationHandler)
164{
165//
166//
167// Send message to stdout to advise that the G4Abrasion model is being used.
168//
169  PrintWelcomeMessage();
170//
171//
172// Set the default verbose level to 0 - no output.
173//
174  verboseLevel = 0;
175//                     
176//
177// The user is able to provide the excitation handler as well as an argument
178// which is provided in this instantiation is used to determine
179// whether the spectators of the interaction are free following the abrasion.
180//
181  theExcitationHandler             = aExcitationHandler;
182  theExcitationHandlerx            = new G4ExcitationHandler;
183  G4Evaporation * theEvaporation   = new G4Evaporation;
184  G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp;
185  theExcitationHandlerx->SetEvaporation(theEvaporation);
186  theExcitationHandlerx->SetFermiModel(theFermiBreakUp);
187  theExcitationHandlerx->SetMaxAandZForFermiBreakUp(12, 6);
188//
189//
190// Set the minimum and maximum range for the model (despite nomanclature, this
191// is in energy per nucleon number). 
192//
193  SetMinEnergy(70.0*MeV);
194  SetMaxEnergy(10.1*GeV);
195  isBlocked = false;
196//
197//
198// npK, when mutiplied by the nuclear Fermi momentum, determines the range of
199// momentum over which the secondary nucleon momentum is sampled.
200//
201  npK              = 5.0;
202  B                = 10.0 * MeV;
203  third            = 1.0 / 3.0;
204  fradius          = 0.99;
205  conserveEnergy   = false;
206  conserveMomentum = true;
207}
208////////////////////////////////////////////////////////////////////////////////
209//
210G4WilsonAbrasionModel::~G4WilsonAbrasionModel ()
211{
212//
213//
214// The destructor doesn't have to do a great deal!
215//
216  if (theExcitationHandler)  delete theExcitationHandler;
217  if (theExcitationHandlerx) delete theExcitationHandlerx;
218  if (useAblation)           delete theAblation;
219//  delete theExcitationHandler;
220//  delete theExcitationHandlerx;
221}
222////////////////////////////////////////////////////////////////////////////////
223//
224G4HadFinalState *G4WilsonAbrasionModel::ApplyYourself (
225  const G4HadProjectile &theTrack, G4Nucleus &theTarget)
226{
227//
228//
229// The secondaries will be returned in G4HadFinalState &theParticleChange -
230// initialise this.  The original track will always be discontinued and
231// secondaries followed.
232//
233  theParticleChange.Clear();
234  theParticleChange.SetStatusChange(stopAndKill);
235//
236//
237// Get relevant information about the projectile and target (A, Z, energy/nuc,
238// momentum, etc).
239//
240  const G4ParticleDefinition *definitionP = theTrack.GetDefinition();
241  const G4double AP   = definitionP->GetBaryonNumber();
242  const G4double ZP   = definitionP->GetPDGCharge();
243  G4LorentzVector pP  = theTrack.Get4Momentum();
244  G4double E          = theTrack.GetKineticEnergy()/AP;
245  G4double AT         = theTarget.GetN();
246  G4double ZT         = theTarget.GetZ();
247  G4double TotalEPre  = theTrack.GetTotalEnergy() +
248    theTarget.AtomicMass(AT, ZT) + theTarget.GetEnergyDeposit();
249  G4double TotalEPost = 0.0;
250//
251//
252// Determine the radii of the projectile and target nuclei.
253//
254  G4WilsonRadius aR;
255  G4double rP   = aR.GetWilsonRadius(AP);
256  G4double rT   = aR.GetWilsonRadius(AT);
257  G4double rPsq = rP * rP;
258  G4double rTsq = rT * rT;
259  if (verboseLevel >= 2)
260  {
261    G4cout <<"########################################"
262           <<"########################################"
263           <<G4endl;
264    G4cout.precision(6);
265    G4cout <<"IN G4WilsonAbrasionModel" <<G4endl;
266    G4cout <<"Initial projectile A=" <<AP
267           <<", Z=" <<ZP
268           <<", radius = " <<rP/fermi <<" fm"
269           <<G4endl; 
270    G4cout <<"Initial target     A=" <<AT
271           <<", Z=" <<ZT
272           <<", radius = " <<rT/fermi <<" fm"
273           <<G4endl;
274    G4cout <<"Projectile momentum and Energy/nuc = " <<pP <<" ," <<E <<G4endl;
275  }
276//
277//
278// The following variables are used to determine the impact parameter in the
279// near-field (i.e. taking into consideration the electrostatic repulsion).
280//
281  G4double rm   = ZP * ZT * elm_coupling / (E * AP);
282  G4double r    = 0.0;
283  G4double rsq  = 0.0;
284//
285//
286// Initialise some of the variables which wll be used to calculate the chord-
287// length for nucleons in the projectile and target, and hence calculate the
288// number of abraded nucleons and the excitation energy.
289//
290  G4NuclearAbrasionGeometry *theAbrasionGeometry = NULL;
291  G4double CT   = 0.0;
292  G4double F    = 0.0;
293  G4int Dabr    = 0;
294//
295//
296// The following loop is performed until the number of nucleons which are
297// abraded by the process is >1, i.e. an interaction MUST occur.
298//
299  while (Dabr == 0)
300  {
301// Added by MHM 20050119 to fix leaking memory on second pass through this loop
302    if (theAbrasionGeometry)
303    {
304      delete theAbrasionGeometry;
305      theAbrasionGeometry = NULL;
306    }
307//
308//
309// Sample the impact parameter.  For the moment, this class takes account of
310// electrostatic effects on the impact parameter, but (like HZETRN AND NUCFRG2)
311// does not make any correction for the effects of nuclear-nuclear repulsion.
312//
313    G4double rPT   = rP + rT;
314    G4double rPTsq = rPT * rPT;
315//
316//
317// This is a "catch" to make sure we don't go into an infinite loop because the
318// energy is too low to overcome nuclear repulsion.  PRT 20091023.  If the
319// value of rm < fradius * rPT then we're unlikely to sample a small enough
320// impact parameter (energy of incident particle is too low).  Return primary
321// and don't complete nuclear interaction analysis.
322//
323    if (rm >= fradius * rPT) {
324      theParticleChange.SetStatusChange(isAlive);
325      theParticleChange.SetEnergyChange(theTrack.GetKineticEnergy());
326      theParticleChange.SetMomentumChange(theTrack.Get4Momentum().vect().unit());
327      if (verboseLevel >= 2) {
328        G4cout <<"Particle energy too low to overcome repulsion." <<G4endl;
329        G4cout <<"Event rejected and original track maintained" <<G4endl;
330        G4cout <<"########################################"
331               <<"########################################"
332               <<G4endl;
333      }
334      return &theParticleChange;
335    }
336//
337//
338// Now sample impact parameter until the criterion is met that projectile
339// and target overlap, but repulsion is taken into consideration.
340//
341    G4int evtcnt   = 0;
342    r              = 1.1 * rPT;
343    while (r > rPT && ++evtcnt < 1000)
344    {
345      G4double bsq = rPTsq * G4UniformRand();
346      r            = (rm + std::sqrt(rm*rm + 4.0*bsq)) / 2.0;
347    }
348//
349//
350// We've tried to sample this 1000 times, but failed.  Assume nuclei do not
351// collide.
352//
353    if (evtcnt >= 1000) {
354      theParticleChange.SetStatusChange(isAlive);
355      theParticleChange.SetEnergyChange(theTrack.GetKineticEnergy());
356      theParticleChange.SetMomentumChange(theTrack.Get4Momentum().vect().unit());
357      if (verboseLevel >= 2) {
358        G4cout <<"Particle energy too low to overcome repulsion." <<G4endl;
359        G4cout <<"Event rejected and original track maintained" <<G4endl;
360        G4cout <<"########################################"
361               <<"########################################"
362               <<G4endl;
363      }
364      return &theParticleChange;
365    }
366
367
368    rsq = r * r;
369//
370//
371// Now determine the chord-length through the target nucleus.
372//
373    if (rT > rP)
374    {
375      G4double x = (rPsq + rsq - rTsq) / 2.0 / r;
376      if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x);
377      else         CT = 2.0 * std::sqrt(rTsq - rsq);
378    }
379    else
380    {
381      G4double x = (rTsq + rsq - rPsq) / 2.0 / r;
382      if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x);
383      else         CT = 2.0 * rT;
384    }
385//
386//
387// Determine the number of abraded nucleons.  Note that the mean number of
388// abraded nucleons is used to sample the Poisson distribution.  The Poisson
389// distribution is sampled only ten times with the current impact parameter,
390// and if it fails after this to find a case for which the number of abraded
391// nucleons >1, the impact parameter is re-sampled.
392//
393    theAbrasionGeometry = new G4NuclearAbrasionGeometry(AP,AT,r);
394    F                   = theAbrasionGeometry->F();
395    G4double lambda     = 16.6*fermi / std::pow(E/MeV,0.26);
396    G4double Mabr       = F * AP * (1.0 - std::exp(-CT/lambda));
397    G4long n            = 0;
398    for (G4int i = 0; i<10; i++)
399    {
400      n = G4Poisson(Mabr);
401      if (n > 0)
402      {
403        if (n>AP) Dabr = (G4int) AP;
404        else      Dabr = (G4int) n;
405        break;
406      }
407    }
408  }
409  if (verboseLevel >= 2)
410  {
411    G4cout <<G4endl;
412    G4cout <<"Impact parameter    = " <<r/fermi <<" fm" <<G4endl;
413    G4cout <<"# Abraded nucleons  = " <<Dabr <<G4endl;
414  }
415//
416//
417// The number of abraded nucleons must be no greater than the number of
418// nucleons in either the projectile or the target.  If AP - Dabr < 2 or
419// AT - Dabr < 2 then either we have only a nucleon left behind in the
420// projectile/target or we've tried to abrade too many nucleons - and Dabr
421// should be limited.
422//
423  if (AP - (G4double) Dabr < 2.0) Dabr = (G4int) AP;
424  if (AT - (G4double) Dabr < 2.0) Dabr = (G4int) AT;
425//
426//
427// Determine the abraded secondary nucleons from the projectile.  *fragmentP
428// is a pointer to the prefragment from the projectile and nSecP is the number
429// of nucleons in theParticleChange which have been abraded.  The total energy
430// from these is determined.
431//
432  G4ThreeVector boost   = pP.findBoostToCM();
433  G4Fragment *fragmentP = GetAbradedNucleons (Dabr, AP, ZP, rP); 
434  G4int nSecP           = theParticleChange.GetNumberOfSecondaries();
435  G4int i               = 0;
436  for (i=0; i<nSecP; i++)
437  {
438    TotalEPost += theParticleChange.GetSecondary(i)->
439      GetParticle()->GetTotalEnergy();
440  }
441//
442//
443// Determine the number of spectators in the interaction region for the
444// projectile.
445//
446  G4int DspcP = (G4int) (AP*F) - Dabr;
447  if (DspcP <= 0)           DspcP = 0;
448  else if (DspcP > AP-Dabr) DspcP = ((G4int) AP) - Dabr;
449//
450//
451// Determine excitation energy associated with excess surface area of the
452// projectile (EsP) and the excitation due to scattering of nucleons which are
453// retained within the projectile (ExP).  Add the total energy from the excited
454// nucleus to the total energy of the secondaries.
455//
456  G4bool excitationAbsorbedByProjectile = false;
457  if (fragmentP != NULL)
458  {
459    G4double EsP = theAbrasionGeometry->GetExcitationEnergyOfProjectile();
460    G4double ExP = 0.0;
461    if (Dabr < AT)
462      excitationAbsorbedByProjectile = G4UniformRand() < 0.5;
463    if (excitationAbsorbedByProjectile)
464      ExP = GetNucleonInducedExcitation(rP, rT, r);
465    G4double xP = EsP + ExP;
466    if (xP > B*(AP-Dabr)) xP = B*(AP-Dabr);
467    G4LorentzVector lorentzVector = fragmentP->GetMomentum();
468    lorentzVector.setE(lorentzVector.e()+xP);
469    fragmentP->SetMomentum(lorentzVector);
470    TotalEPost += lorentzVector.e();
471  }
472  G4double EMassP = TotalEPost;
473//
474//
475// Determine the abraded secondary nucleons from the target.  Note that it's
476// assumed that the same number of nucleons are abraded from the target as for
477// the projectile, and obviously no boost is applied to the products. *fragmentT
478// is a pointer to the prefragment from the target and nSec is the total number
479// of nucleons in theParticleChange which have been abraded.  The total energy
480// from these is determined.
481//
482  G4Fragment *fragmentT = GetAbradedNucleons (Dabr, AT, ZT, rT); 
483  G4int nSec = theParticleChange.GetNumberOfSecondaries();
484  for (i=nSecP; i<nSec; i++)
485  {
486    TotalEPost += theParticleChange.GetSecondary(i)->
487      GetParticle()->GetTotalEnergy();
488  }
489//
490//
491// Determine the number of spectators in the interaction region for the
492// target.
493//
494  G4int DspcT = (G4int) (AT*F) - Dabr;
495  if (DspcT <= 0)           DspcT = 0;
496  else if (DspcT > AP-Dabr) DspcT = ((G4int) AT) - Dabr;
497//
498//
499// Determine excitation energy associated with excess surface area of the
500// target (EsT) and the excitation due to scattering of nucleons which are
501// retained within the target (ExT).  Add the total energy from the excited
502// nucleus to the total energy of the secondaries.
503//
504  if (fragmentT != NULL)
505  {
506    G4double EsT = theAbrasionGeometry->GetExcitationEnergyOfTarget();
507    G4double ExT = 0.0;
508    if (!excitationAbsorbedByProjectile)
509      ExT = GetNucleonInducedExcitation(rT, rP, r);
510    G4double xT = EsT + ExT;
511    if (xT > B*(AT-Dabr)) xT = B*(AT-Dabr);
512    G4LorentzVector lorentzVector = fragmentT->GetMomentum();
513    lorentzVector.setE(lorentzVector.e()+xT);
514    fragmentT->SetMomentum(lorentzVector);
515    TotalEPost += lorentzVector.e();
516  }
517//
518//
519// Now determine the difference between the pre and post interaction
520// energy - this will be used to determine the Lorentz boost if conservation
521// of energy is to be imposed/attempted.
522//
523  G4double deltaE = TotalEPre - TotalEPost;
524  if (deltaE > 0.0 && conserveEnergy)
525  {
526    G4double beta = std::sqrt(1.0 - EMassP*EMassP/std::pow(deltaE+EMassP,2.0));
527    boost = boost / boost.mag() * beta;
528  }
529//
530//
531// Now boost the secondaries from the projectile.
532//
533  G4ThreeVector pBalance = pP.vect();
534  for (i=0; i<nSecP; i++)
535  {
536    G4DynamicParticle *dynamicP = theParticleChange.GetSecondary(i)->
537      GetParticle();
538    G4LorentzVector lorentzVector = dynamicP->Get4Momentum();
539    lorentzVector.boost(-boost);
540    dynamicP->Set4Momentum(lorentzVector);
541    pBalance -= lorentzVector.vect();
542  }
543//
544//
545// Set the boost for the projectile prefragment.  This is now based on the
546// conservation of momentum.  However, if the user selected momentum of the
547// prefragment is not to be conserved this simply boosted to the velocity of the
548// original projectile times the ratio of the unexcited to the excited mass
549// of the prefragment (the excitation increases the effective mass of the
550// prefragment, and therefore modifying the boost is an attempt to prevent
551// the momentum of the prefragment being excessive).
552//
553  if (fragmentP != NULL)
554  {
555    G4LorentzVector lorentzVector = fragmentP->GetMomentum();
556    G4double m                    = lorentzVector.m();
557    if (conserveMomentum)
558      fragmentP->SetMomentum
559        (G4LorentzVector(pBalance,std::sqrt(pBalance.mag2()+m*m+1.0*eV*eV)));
560    else
561    {
562      G4double mg = fragmentP->GetGroundStateMass();
563      fragmentP->SetMomentum(lorentzVector.boost(-boost * mg/m));
564    }
565  }
566//
567//
568// Output information to user if verbose information requested.
569//
570  if (verboseLevel >= 2)
571  {
572    G4cout <<G4endl;
573    G4cout <<"-----------------------------------" <<G4endl;
574    G4cout <<"Secondary nucleons from projectile:" <<G4endl;
575    G4cout <<"-----------------------------------" <<G4endl;
576    G4cout.precision(7);
577    for (i=0; i<nSecP; i++)
578    {
579      G4cout <<"Particle # " <<i <<G4endl;
580      theParticleChange.GetSecondary(i)->GetParticle()->DumpInfo();
581      G4DynamicParticle *dyn = theParticleChange.GetSecondary(i)->GetParticle();
582      G4cout <<"New nucleon (P) " <<dyn->GetDefinition()->GetParticleName()
583             <<" : "              <<dyn->Get4Momentum()
584             <<G4endl;
585    }
586    G4cout <<"---------------------------" <<G4endl;
587    G4cout <<"The projectile prefragment:" <<G4endl;
588    G4cout <<"---------------------------" <<G4endl;
589    if (fragmentP != NULL)
590      G4cout <<*fragmentP <<G4endl;
591    else
592      G4cout <<"(No residual prefragment)" <<G4endl;
593    G4cout <<G4endl;
594    G4cout <<"-------------------------------" <<G4endl;
595    G4cout <<"Secondary nucleons from target:" <<G4endl;
596    G4cout <<"-------------------------------" <<G4endl;
597    G4cout.precision(7);
598    for (i=nSecP; i<nSec; i++)
599    {
600      G4cout <<"Particle # " <<i <<G4endl;
601      theParticleChange.GetSecondary(i)->GetParticle()->DumpInfo();
602      G4DynamicParticle *dyn = theParticleChange.GetSecondary(i)->GetParticle();
603      G4cout <<"New nucleon (T) " <<dyn->GetDefinition()->GetParticleName()
604             <<" : "              <<dyn->Get4Momentum()
605             <<G4endl;
606    }
607    G4cout <<"-----------------------" <<G4endl;
608    G4cout <<"The target prefragment:" <<G4endl;
609    G4cout <<"-----------------------" <<G4endl;
610    if (fragmentT != NULL)
611      G4cout <<*fragmentT <<G4endl;
612    else
613      G4cout <<"(No residual prefragment)" <<G4endl;
614  }
615//
616//
617// Now we can decay the nuclear fragments if present.  The secondaries are
618// collected and boosted as well.  This is performed first for the projectile...
619//
620  if (fragmentP !=NULL)
621  {
622    G4ReactionProductVector *products = NULL;
623    if (fragmentP->GetZ() != fragmentP->GetA())
624      products = theExcitationHandler->BreakItUp(*fragmentP);
625    else
626      products = theExcitationHandlerx->BreakItUp(*fragmentP);     
627    delete fragmentP;
628    fragmentP = NULL;
629 
630    G4ReactionProductVector::iterator iter;
631    for (iter = products->begin(); iter != products->end(); ++iter)
632    {
633      G4DynamicParticle *secondary =
634        new G4DynamicParticle((*iter)->GetDefinition(),
635        (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
636      theParticleChange.AddSecondary (secondary);  // Added MHM 20050118
637      G4String particleName = (*iter)->GetDefinition()->GetParticleName();
638      delete (*iter); // get rid of leftover particle def!  // Added MHM 20050118
639      if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size())
640      {
641        G4cout <<"------------------------" <<G4endl;
642        G4cout <<"The projectile fragment:" <<G4endl;
643        G4cout <<"------------------------" <<G4endl;
644        G4cout <<" fragmentP = " <<particleName
645               <<" Energy    = " <<secondary->GetKineticEnergy()
646               <<G4endl;
647      }
648    }
649    delete products;  // Added MHM 20050118
650  }
651//
652//
653// Now decay the target nucleus - no boost is applied since in this
654// approximation it is assumed that there is negligible momentum transfer from
655// the projectile.
656//
657  if (fragmentT != NULL)
658  {
659    G4ReactionProductVector *products = NULL;
660    if (fragmentT->GetZ() != fragmentT->GetA())
661      products = theExcitationHandler->BreakItUp(*fragmentT);
662    else
663      products = theExcitationHandlerx->BreakItUp(*fragmentT);     
664    delete fragmentT;
665    fragmentT = NULL;
666 
667    G4ReactionProductVector::iterator iter;
668    for (iter = products->begin(); iter != products->end(); ++iter)
669    {
670      G4DynamicParticle *secondary =
671        new G4DynamicParticle((*iter)->GetDefinition(),
672        (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
673      theParticleChange.AddSecondary (secondary);
674      G4String particleName = (*iter)->GetDefinition()->GetParticleName();
675      delete (*iter); // get rid of leftover particle def!  // Added MHM 20050118
676      if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size())
677      {
678        G4cout <<"--------------------" <<G4endl;
679        G4cout <<"The target fragment:" <<G4endl;
680        G4cout <<"--------------------" <<G4endl;
681        G4cout <<" fragmentT = " <<particleName
682               <<" Energy    = " <<secondary->GetKineticEnergy()
683               <<G4endl;
684      }
685    }
686    delete products;  // Added MHM 20050118
687  }
688
689  if (verboseLevel >= 2)
690     G4cout <<"########################################"
691            <<"########################################"
692            <<G4endl;
693 
694  delete theAbrasionGeometry;
695 
696  return &theParticleChange;
697}
698////////////////////////////////////////////////////////////////////////////////
699//
700G4Fragment *G4WilsonAbrasionModel::GetAbradedNucleons (G4int Dabr, G4double A,
701  G4double Z, G4double r)
702{
703//
704//
705// Initialise variables.  tau is the Fermi radius of the nucleus.  The variables
706// p..., C... and g(amma) are used to help sample the secondary nucleon
707// spectrum.
708//
709 
710  G4double pK   = hbarc * std::pow(9.0 * pi / 4.0 * A, third) / (1.29 * r);
711  if (A <= 24.0) pK *= -0.229*std::pow(A,third) + 1.62; 
712  G4double pKsq = pK * pK;
713  G4double p1sq = 2.0/5.0 * pKsq;
714  G4double p2sq = 6.0/5.0 * pKsq;
715  G4double p3sq = 500.0 * 500.0;
716  G4double C1   = 1.0;
717  G4double C2   = 0.03;
718  G4double C3   = 0.0002;
719  G4double g    = 90.0 * MeV;
720  G4double maxn = C1 + C2 + C3;
721//
722//
723// initialise the number of secondary nucleons abraded to zero, and initially set
724// the type of nucleon abraded to proton ... just for now.
725// 
726  G4double Aabr                     = 0.0;
727  G4double Zabr                     = 0.0; 
728  G4ParticleDefinition *typeNucleon = G4Proton::ProtonDefinition();
729  G4DynamicParticle *dynamicNucleon = NULL;
730  G4ParticleMomentum pabr(0.0, 0.0, 0.0);
731//
732//
733// Now go through each abraded nucleon and sample type, spectrum and angle.
734//
735  for (G4int i=0; i<Dabr; i++)
736  {
737//
738//
739// Sample the nucleon momentum distribution by simple rejection techniques.  We
740// reject values of p == 0.0 since this causes bad behaviour in the sinh term.
741//
742    G4double p   = 0.0;
743    G4bool found = false;
744    while (!found)
745    {
746      while (p <= 0.0) p = npK * pK * G4UniformRand();
747      G4double psq = p * p;
748      found = maxn * G4UniformRand() < C1*std::exp(-psq/p1sq/2.0) +
749        C2*std::exp(-psq/p2sq/2.0) + C3*std::exp(-psq/p3sq/2.0) + p/g/std::sinh(p/g);
750    }
751//
752//
753// Determine the type of particle abraded.  Can only be proton or neutron,
754// and the probability is determine to be proportional to the ratio as found
755// in the nucleus at each stage.
756//
757    G4double prob = (Z-Zabr)/(A-Aabr);
758    if (G4UniformRand()<prob)
759    {
760      Zabr++;
761      typeNucleon = G4Proton::ProtonDefinition();
762    }
763    else
764      typeNucleon = G4Neutron::NeutronDefinition();
765    Aabr++;
766//
767//
768// The angular distribution of the secondary nucleons is approximated to an
769// isotropic distribution in the rest frame of the nucleus (this will be Lorentz
770// boosted later.
771//
772    G4double costheta = 2.*G4UniformRand()-1.0;
773    G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
774    G4double phi      = 2.0*pi*G4UniformRand()*rad;
775    G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
776    G4double nucleonMass = typeNucleon->GetPDGMass();
777    G4double E           = std::sqrt(p*p + nucleonMass*nucleonMass)-nucleonMass;
778    dynamicNucleon = new G4DynamicParticle(typeNucleon,direction,E);
779    theParticleChange.AddSecondary (dynamicNucleon);
780    pabr += p*direction;
781  }
782//
783//
784// Next determine the details of the nuclear prefragment .. that is if there
785// is one or more protons in the residue.  (Note that the 1 eV in the total
786// energy is a safety factor to avoid any possibility of negative rest mass
787// energy.)
788//
789  G4Fragment *fragment = NULL;
790  if (Z-Zabr>=1.0)
791  {
792    G4double ionMass = G4ParticleTable::GetParticleTable()->GetIonTable()->
793                       GetIonMass(G4lrint(Z-Zabr),G4lrint(A-Aabr));
794    G4double E       = std::sqrt(pabr.mag2() + ionMass*ionMass);
795    G4LorentzVector lorentzVector = G4LorentzVector(-pabr, E + 1.0*eV);
796    fragment =
797      new G4Fragment((G4int) (A-Aabr), (G4int) (Z-Zabr), lorentzVector);
798  }
799
800  return fragment;
801}
802////////////////////////////////////////////////////////////////////////////////
803//
804G4double G4WilsonAbrasionModel::GetNucleonInducedExcitation
805  (G4double rP, G4double rT, G4double r)
806{
807//
808//
809// Initialise variables.
810//
811  G4double Cl   = 0.0;
812  G4double rPsq = rP * rP;
813  G4double rTsq = rT * rT;
814  G4double rsq  = r * r;
815//
816//
817// Depending upon the impact parameter, a different form of the chord length is
818// is used.
819// 
820  if (r > rT) Cl = 2.0*std::sqrt(rPsq + 2.0*r*rT - rsq - rTsq);
821  else        Cl = 2.0*rP;
822//
823//
824// The next lines have been changed to include a "catch" to make sure if the
825// projectile and target are too close, Ct is set to twice rP or twice rT.
826// Otherwise the standard Wilson algorithm should work fine.
827// PRT 20091023.
828//
829  G4double Ct = 0.0;
830  if      (rT > rP && rsq < rTsq - rPsq) Ct = 2.0 * rP;
831  else if (rP > rT && rsq < rPsq - rTsq) Ct = 2.0 * rT;
832  else {
833    G4double bP = (rPsq+rsq-rTsq)/2.0/r;
834    G4double x  = rPsq - bP*bP;
835    if (x < 0.0) {
836      G4cerr <<"########################################"
837             <<"########################################"
838             <<G4endl;
839      G4cerr <<"ERROR IN G4WilsonAbrasionModel::GetNucleonInducedExcitation"
840             <<G4endl;
841      G4cerr <<"rPsq - bP*bP < 0.0 and cannot be square-rooted" <<G4endl;
842      G4cerr <<"Set to zero instead" <<G4endl;
843      G4cerr <<"########################################"
844             <<"########################################"
845             <<G4endl;
846    }
847    Ct = 2.0*std::sqrt(x);
848  }
849 
850  G4double Ex = 13.0 * Cl / fermi;
851  if (Ct > 1.5*fermi)
852    Ex += 13.0 * Cl / fermi /3.0 * (Ct/fermi - 1.5);
853
854  return Ex;
855}
856////////////////////////////////////////////////////////////////////////////////
857//
858void G4WilsonAbrasionModel::SetUseAblation (G4bool useAblation1)
859{
860  if (useAblation != useAblation1)
861  {
862    useAblation = useAblation1;
863    delete theExcitationHandler;
864    delete theExcitationHandlerx;
865    theExcitationHandler  = new G4ExcitationHandler;
866    theExcitationHandlerx = new G4ExcitationHandler;
867    if (useAblation)
868    {
869      theAblation = new G4WilsonAblationModel;
870      theAblation->SetVerboseLevel(verboseLevel);
871      theExcitationHandler->SetEvaporation(theAblation);
872      theExcitationHandlerx->SetEvaporation(theAblation);
873    }
874    else
875    {
876      theAblation                      = NULL;
877      G4Evaporation * theEvaporation   = new G4Evaporation;
878      G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp;
879      G4StatMF * theMF                 = new G4StatMF;
880      theExcitationHandler->SetEvaporation(theEvaporation);
881      theExcitationHandler->SetFermiModel(theFermiBreakUp);
882      theExcitationHandler->SetMultiFragmentation(theMF);
883      theExcitationHandler->SetMaxAandZForFermiBreakUp(12, 6);
884      theExcitationHandler->SetMinEForMultiFrag(5.0*MeV);
885
886      theEvaporation  = new G4Evaporation;
887      theFermiBreakUp = new G4FermiBreakUp;
888      theExcitationHandlerx->SetEvaporation(theEvaporation);
889      theExcitationHandlerx->SetFermiModel(theFermiBreakUp);
890      theExcitationHandlerx->SetMaxAandZForFermiBreakUp(12, 6);
891    }
892  }
893  return; 
894}
895////////////////////////////////////////////////////////////////////////////////
896//
897void G4WilsonAbrasionModel::PrintWelcomeMessage ()
898{
899  G4cout <<G4endl;
900  G4cout <<" *****************************************************************"
901         <<G4endl;
902  G4cout <<" Nuclear abrasion model for nuclear-nuclear interactions activated"
903         <<G4endl;
904  G4cout <<" (Written by QinetiQ Ltd for the European Space Agency)"
905         <<G4endl;
906  G4cout <<" *****************************************************************"
907         <<G4endl;
908  G4cout << G4endl;
909
910  return;
911}
912////////////////////////////////////////////////////////////////////////////////
913//
Note: See TracBrowser for help on using the repository browser.