Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/abrasion/src/G4WilsonAbrasionModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // 12 June 2012, A. Ribon, CERN, Switzerland
 79 // Fixing trivial warning errors of shadowed variables. 
 80 //
 81 // 4 August 2015, A. Ribon, CERN, Switzerland
 82 // Replacing std::exp and std::pow with the faster versions G4Exp and G4Pow.
 83 //
 84 // 7 August 2015, A. Ribon, CERN, Switzerland
 85 // Checking of 'while' loops.
 86 //
 87 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
 88 ///////////////////////////////////////////////////////////////////////////////
 89 
 90 #include "G4WilsonAbrasionModel.hh"
 91 #include "G4WilsonRadius.hh"
 92 #include "G4NuclearAbrasionGeometry.hh"
 93 #include "G4WilsonAblationModel.hh"
 94 
 95 #include "G4PhysicalConstants.hh"
 96 #include "G4SystemOfUnits.hh"
 97 #include "G4ExcitationHandler.hh"
 98 #include "G4Evaporation.hh"
 99 #include "G4ParticleDefinition.hh"
100 #include "G4DynamicParticle.hh"
101 #include "Randomize.hh"
102 #include "G4Fragment.hh"
103 #include "G4ReactionProductVector.hh"
104 #include "G4LorentzVector.hh"
105 #include "G4ParticleMomentum.hh"
106 #include "G4Poisson.hh"
107 #include "G4ParticleTable.hh"
108 #include "G4IonTable.hh"
109 #include "globals.hh"
110 
111 #include "G4Exp.hh"
112 #include "G4Pow.hh"
113 
114 #include "G4PhysicsModelCatalog.hh"
115 
116 
117 G4WilsonAbrasionModel::G4WilsonAbrasionModel(G4bool useAblation1)
118   : G4HadronicInteraction("G4WilsonAbrasion"), secID(-1)
119 {
120   // Send message to stdout to advise that the G4Abrasion model is being used.
121   PrintWelcomeMessage();
122 
123   // Set the default verbose level to 0 - no output.
124   verboseLevel = 0;
125   useAblation  = useAblation1;
126   theAblation  = nullptr;
127 
128   // No de-excitation handler has been supplied - define the default handler.
129 
130   theExcitationHandler = new G4ExcitationHandler();
131   if (useAblation)
132   {
133     theAblation = new G4WilsonAblationModel;
134     theAblation->SetVerboseLevel(verboseLevel);
135     theExcitationHandler->SetEvaporation(theAblation, true);
136   }
137 
138   // Set the minimum and maximum range for the model (despite nomanclature,
139   // this is in energy per nucleon number).  
140 
141   SetMinEnergy(70.0*MeV);
142   SetMaxEnergy(10.1*GeV);
143   isBlocked = false;
144 
145   // npK, when mutiplied by the nuclear Fermi momentum, determines the range of
146   // momentum over which the secondary nucleon momentum is sampled.
147 
148   r0sq = 0.0;
149   npK = 5.0;
150   B = 10.0 * MeV;
151   third = 1.0 / 3.0;
152   fradius = 0.99;
153   conserveEnergy = false;
154   conserveMomentum = true;
155 
156   // Creator model ID for the secondaries created by this model
157   secID = G4PhysicsModelCatalog::GetModelID( "model_" + GetModelName() );  
158 }
159 
160 void G4WilsonAbrasionModel::ModelDescription(std::ostream& outFile) const
161 {
162   outFile << "G4WilsonAbrasionModel is a macroscopic treatment of\n"
163           << "nucleus-nucleus collisions using simple geometric arguments.\n"
164           << "The smaller projectile nucleus gouges out a part of the larger\n"
165           << "target nucleus, leaving a residual nucleus and a fireball\n"
166           << "region where the projectile and target intersect.  The fireball"
167           << "is then treated as a highly excited nuclear fragment.  This\n"
168           << "model is based on the NUCFRG2 model and is valid for all\n"
169           << "projectile energies between 70 MeV/n and 10.1 GeV/n. \n";
170 }
171 
172 G4WilsonAbrasionModel::G4WilsonAbrasionModel(G4ExcitationHandler* aExcitationHandler) :
173  G4HadronicInteraction("G4WilsonAbrasion"), secID(-1)  
174 {
175 // Send message to stdout to advise that the G4Abrasion model is being used.
176 
177   PrintWelcomeMessage();
178 
179 // Set the default verbose level to 0 - no output.
180 
181   verboseLevel = 0;
182 
183   theAblation = nullptr;   //A.R. 26-Jul-2012 Coverity fix.
184   useAblation = false;     //A.R. 14-Aug-2012 Coverity fix.
185                       
186 //
187 // The user is able to provide the excitation handler as well as an argument
188 // which is provided in this instantiation is used to determine
189 // whether the spectators of the interaction are free following the abrasion.
190 //
191   theExcitationHandler = aExcitationHandler;
192 //
193 //
194 // Set the minimum and maximum range for the model (despite nomanclature, this
195 // is in energy per nucleon number).  
196 //
197   SetMinEnergy(70.0*MeV);
198   SetMaxEnergy(10.1*GeV);
199   isBlocked = false;
200 //
201 //
202 // npK, when mutiplied by the nuclear Fermi momentum, determines the range of
203 // momentum over which the secondary nucleon momentum is sampled.
204 //
205   r0sq             = 0.0;
206   npK              = 5.0;
207   B                = 10.0 * MeV;
208   third            = 1.0 / 3.0;
209   fradius          = 0.99;
210   conserveEnergy   = false;
211   conserveMomentum = true;
212 
213   // Creator model ID for the secondaries created by this model
214   secID = G4PhysicsModelCatalog::GetModelID( "model_" + GetModelName() );  
215 }
216 ////////////////////////////////////////////////////////////////////////////////
217 //
218 G4WilsonAbrasionModel::~G4WilsonAbrasionModel()
219 {
220   delete theExcitationHandler;
221 }
222 ////////////////////////////////////////////////////////////////////////////////
223 //
224 G4HadFinalState *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.GetA_asInt();
246   G4double ZT        = theTarget.GetZ_asInt();
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 = nullptr;
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   G4bool skipInteraction = false;  // It will be set true if the two nuclei fail to collide 
300   const G4int maxNumberOfLoops = 1000;
301   G4int loopCounter = -1;
302   while (Dabr == 0 && ++loopCounter < maxNumberOfLoops)  /* Loop checking, 07.08.2015, A.Ribon */
303   {
304 //
305 //
306 // Sample the impact parameter.  For the moment, this class takes account of
307 // electrostatic effects on the impact parameter, but (like HZETRN AND NUCFRG2)
308 // does not make any correction for the effects of nuclear-nuclear repulsion.
309 //
310     G4double rPT   = rP + rT;
311     G4double rPTsq = rPT * rPT;
312 //
313 //
314 // This is a "catch" to make sure we don't go into an infinite loop because the
315 // energy is too low to overcome nuclear repulsion.  PRT 20091023.  If the
316 // value of rm < fradius * rPT then we're unlikely to sample a small enough
317 // impact parameter (energy of incident particle is too low).
318 //
319     if (rm >= fradius * rPT) {
320       skipInteraction = true;
321     }
322 //
323 //
324 // Now sample impact parameter until the criterion is met that projectile
325 // and target overlap, but repulsion is taken into consideration.
326 //
327     G4int evtcnt   = 0;
328     r              = 1.1 * rPT;
329     while (r > rPT && ++evtcnt < 1000)  /* Loop checking, 07.08.2015, A.Ribon */
330     {
331       G4double bsq = rPTsq * G4UniformRand();
332       r            = (rm + std::sqrt(rm*rm + 4.0*bsq)) / 2.0;
333     }
334 //
335 //
336 // We've tried to sample this 1000 times, but failed.  
337 //
338     if (evtcnt >= 1000) {
339       skipInteraction = true;
340     }
341 
342     rsq = r * r;
343 //
344 //
345 // Now determine the chord-length through the target nucleus.
346 //
347     if (rT > rP)
348     {
349       G4double x = (rPsq + rsq - rTsq) / 2.0 / r;
350       if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x);
351       else         CT = 2.0 * std::sqrt(rTsq - rsq);
352     }
353     else
354     {
355       G4double x = (rTsq + rsq - rPsq) / 2.0 / r;
356       if (x > 0.0) CT = 2.0 * std::sqrt(rTsq - x*x);
357       else         CT = 2.0 * rT;
358     }
359 //
360 //
361 // Determine the number of abraded nucleons.  Note that the mean number of
362 // abraded nucleons is used to sample the Poisson distribution.  The Poisson
363 // distribution is sampled only ten times with the current impact parameter,
364 // and if it fails after this to find a case for which the number of abraded
365 // nucleons >1, the impact parameter is re-sampled.
366 //
367     delete theAbrasionGeometry;
368     theAbrasionGeometry = new G4NuclearAbrasionGeometry(AP,AT,r);
369     F                   = theAbrasionGeometry->F();
370     G4double lambda     = 16.6*fermi / G4Pow::GetInstance()->powA(E/MeV,0.26);
371     G4double Mabr       = F * AP * (1.0 - G4Exp(-CT/lambda));
372     G4long n            = 0;
373     for (G4int i = 0; i<10; ++i)
374     {
375       n = G4Poisson(Mabr);
376       if (n > 0)
377       {
378         if (n>AP) Dabr = (G4int) AP;
379         else      Dabr = (G4int) n;
380         break;
381       }
382     }
383   }  // End of while loop
384 
385   if ( loopCounter >= maxNumberOfLoops || skipInteraction ) {
386     // Assume nuclei do not collide and return unchanged primary. 
387     theParticleChange.SetStatusChange(isAlive);
388     theParticleChange.SetEnergyChange(theTrack.GetKineticEnergy());
389     theParticleChange.SetMomentumChange(theTrack.Get4Momentum().vect().unit());
390     if (verboseLevel >= 2) {
391       G4cout <<"Particle energy too low to overcome repulsion." <<G4endl;
392       G4cout <<"Event rejected and original track maintained" <<G4endl;
393       G4cout <<"########################################"
394              <<"########################################"
395              <<G4endl;
396     }
397     delete theAbrasionGeometry;
398     return &theParticleChange;
399   }
400 
401   if (verboseLevel >= 2)
402   {
403     G4cout <<G4endl;
404     G4cout <<"Impact parameter    = " <<r/fermi <<" fm" <<G4endl;
405     G4cout <<"# Abraded nucleons  = " <<Dabr <<G4endl;
406   }
407 //
408 //
409 // The number of abraded nucleons must be no greater than the number of
410 // nucleons in either the projectile or the target.  If AP - Dabr < 2 or 
411 // AT - Dabr < 2 then either we have only a nucleon left behind in the
412 // projectile/target or we've tried to abrade too many nucleons - and Dabr
413 // should be limited.
414 //
415   if (AP - (G4double) Dabr < 2.0) Dabr = (G4int) AP;
416   if (AT - (G4double) Dabr < 2.0) Dabr = (G4int) AT;
417 //
418 //
419 // Determine the abraded secondary nucleons from the projectile.  *fragmentP
420 // is a pointer to the prefragment from the projectile and nSecP is the number
421 // of nucleons in theParticleChange which have been abraded.  The total energy
422 // from these is determined.
423 //
424   G4ThreeVector boost   = pP.findBoostToCM();
425   G4Fragment *fragmentP = GetAbradedNucleons (Dabr, AP, ZP, rP); 
426   G4int nSecP           = (G4int)theParticleChange.GetNumberOfSecondaries();
427   G4int i               = 0;
428   for (i=0; i<nSecP; ++i)
429   {
430     TotalEPost += theParticleChange.GetSecondary(i)->
431       GetParticle()->GetTotalEnergy();
432   }
433 //
434 //
435 // Determine the number of spectators in the interaction region for the
436 // projectile.
437 //
438   G4int DspcP = (G4int) (AP*F) - Dabr;
439   if (DspcP <= 0)           DspcP = 0;
440   else if (DspcP > AP-Dabr) DspcP = ((G4int) AP) - Dabr;
441 //
442 //
443 // Determine excitation energy associated with excess surface area of the
444 // projectile (EsP) and the excitation due to scattering of nucleons which are
445 // retained within the projectile (ExP).  Add the total energy from the excited
446 // nucleus to the total energy of the secondaries.
447 //
448   G4bool excitationAbsorbedByProjectile = false;
449   if (fragmentP != nullptr)
450   {
451     G4double EsP = theAbrasionGeometry->GetExcitationEnergyOfProjectile();
452     G4double ExP = 0.0;
453     if (Dabr < AT)
454       excitationAbsorbedByProjectile = G4UniformRand() < 0.5;
455     if (excitationAbsorbedByProjectile)
456       ExP = GetNucleonInducedExcitation(rP, rT, r);
457     G4double xP = EsP + ExP;
458     if (xP > B*(AP-Dabr)) xP = B*(AP-Dabr);
459     G4LorentzVector lorentzVector = fragmentP->GetMomentum();
460     lorentzVector.setE(lorentzVector.e()+xP);
461     fragmentP->SetMomentum(lorentzVector);
462     TotalEPost += lorentzVector.e();
463   }
464   G4double EMassP = TotalEPost;
465 //
466 //
467 // Determine the abraded secondary nucleons from the target.  Note that it's
468 // assumed that the same number of nucleons are abraded from the target as for
469 // the projectile, and obviously no boost is applied to the products. *fragmentT
470 // is a pointer to the prefragment from the target and nSec is the total number
471 // of nucleons in theParticleChange which have been abraded.  The total energy
472 // from these is determined.
473 //
474   G4Fragment *fragmentT = GetAbradedNucleons (Dabr, AT, ZT, rT); 
475   G4int nSec = (G4int)theParticleChange.GetNumberOfSecondaries();
476   for (i=nSecP; i<nSec; ++i)
477   {
478     TotalEPost += theParticleChange.GetSecondary(i)->
479       GetParticle()->GetTotalEnergy();
480   }
481 //
482 //
483 // Determine the number of spectators in the interaction region for the
484 // target.
485 //
486   G4int DspcT = (G4int) (AT*F) - Dabr;
487   if (DspcT <= 0)           DspcT = 0;
488   else if (DspcT > AP-Dabr) DspcT = ((G4int) AT) - Dabr;
489 //
490 //
491 // Determine excitation energy associated with excess surface area of the
492 // target (EsT) and the excitation due to scattering of nucleons which are
493 // retained within the target (ExT).  Add the total energy from the excited
494 // nucleus to the total energy of the secondaries.
495 //
496   if (fragmentT != nullptr)
497   {
498     G4double EsT = theAbrasionGeometry->GetExcitationEnergyOfTarget();
499     G4double ExT = 0.0;
500     if (!excitationAbsorbedByProjectile)
501       ExT = GetNucleonInducedExcitation(rT, rP, r);
502     G4double xT = EsT + ExT;
503     if (xT > B*(AT-Dabr)) xT = B*(AT-Dabr);
504     G4LorentzVector lorentzVector = fragmentT->GetMomentum();
505     lorentzVector.setE(lorentzVector.e()+xT);
506     fragmentT->SetMomentum(lorentzVector);
507     TotalEPost += lorentzVector.e();
508   }
509 //
510 //
511 // Now determine the difference between the pre and post interaction
512 // energy - this will be used to determine the Lorentz boost if conservation
513 // of energy is to be imposed/attempted.
514 //
515   G4double deltaE = TotalEPre - TotalEPost;
516   if (deltaE > 0.0 && conserveEnergy)
517   {
518     G4double beta = std::sqrt(1.0 - EMassP*EMassP/G4Pow::GetInstance()->powN(deltaE+EMassP,2));
519     boost = boost / boost.mag() * beta;
520   }
521 //
522 //
523 // Now boost the secondaries from the projectile.
524 //
525   G4ThreeVector pBalance = pP.vect();
526   for (i=0; i<nSecP; ++i)
527   {
528     G4DynamicParticle *dynamicP = theParticleChange.GetSecondary(i)->
529       GetParticle();
530     G4LorentzVector lorentzVector = dynamicP->Get4Momentum();
531     lorentzVector.boost(-boost);
532     dynamicP->Set4Momentum(lorentzVector);
533     pBalance -= lorentzVector.vect();
534   }
535 //
536 //
537 // Set the boost for the projectile prefragment.  This is now based on the
538 // conservation of momentum.  However, if the user selected momentum of the
539 // prefragment is not to be conserved this simply boosted to the velocity of the
540 // original projectile times the ratio of the unexcited to the excited mass
541 // of the prefragment (the excitation increases the effective mass of the
542 // prefragment, and therefore modifying the boost is an attempt to prevent
543 // the momentum of the prefragment being excessive).
544 //
545   if (fragmentP != nullptr)
546   {
547     G4LorentzVector lorentzVector = fragmentP->GetMomentum();
548     G4double fragmentM            = lorentzVector.m();
549     if (conserveMomentum)
550       fragmentP->SetMomentum
551         (G4LorentzVector(pBalance,std::sqrt(pBalance.mag2()+fragmentM*fragmentM+1.0*eV*eV)));
552     else
553     {
554       G4double fragmentGroundStateM = fragmentP->GetGroundStateMass();
555       fragmentP->SetMomentum(lorentzVector.boost(-boost * fragmentGroundStateM/fragmentM));
556     }
557   }
558 //
559 //
560 // Output information to user if verbose information requested.
561 //
562   if (verboseLevel >= 2)
563   {
564     G4cout <<G4endl;
565     G4cout <<"-----------------------------------" <<G4endl;
566     G4cout <<"Secondary nucleons from projectile:" <<G4endl;
567     G4cout <<"-----------------------------------" <<G4endl;
568     G4cout.precision(7);
569     for (i=0; i<nSecP; ++i)
570     {
571       G4cout <<"Particle # " <<i <<G4endl;
572       theParticleChange.GetSecondary(i)->GetParticle()->DumpInfo();
573       G4DynamicParticle *dyn = theParticleChange.GetSecondary(i)->GetParticle();
574       G4cout <<"New nucleon (P) " <<dyn->GetDefinition()->GetParticleName()
575              <<" : "              <<dyn->Get4Momentum()
576              <<G4endl;
577     }
578     G4cout <<"---------------------------" <<G4endl;
579     G4cout <<"The projectile prefragment:" <<G4endl;
580     G4cout <<"---------------------------" <<G4endl;
581     if (fragmentP != nullptr)
582       G4cout <<*fragmentP <<G4endl;
583     else
584       G4cout <<"(No residual prefragment)" <<G4endl;
585     G4cout <<G4endl;
586     G4cout <<"-------------------------------" <<G4endl;
587     G4cout <<"Secondary nucleons from target:" <<G4endl;
588     G4cout <<"-------------------------------" <<G4endl;
589     G4cout.precision(7);
590     for (i=nSecP; i<nSec; ++i)
591     {
592       G4cout <<"Particle # " <<i <<G4endl;
593       theParticleChange.GetSecondary(i)->GetParticle()->DumpInfo();
594       G4DynamicParticle *dyn = theParticleChange.GetSecondary(i)->GetParticle();
595       G4cout <<"New nucleon (T) " <<dyn->GetDefinition()->GetParticleName()
596              <<" : "              <<dyn->Get4Momentum()
597              <<G4endl;
598     }
599     G4cout <<"-----------------------" <<G4endl;
600     G4cout <<"The target prefragment:" <<G4endl;
601     G4cout <<"-----------------------" <<G4endl;
602     if (fragmentT != nullptr)
603       G4cout <<*fragmentT <<G4endl;
604     else
605       G4cout <<"(No residual prefragment)" <<G4endl;
606   }
607 //
608 //
609 // Now we can decay the nuclear fragments if present.  The secondaries are
610 // collected and boosted as well.  This is performed first for the projectile...
611 //
612   if (fragmentP !=nullptr)
613   {
614     G4ReactionProductVector *products = nullptr;
615     //    if (fragmentP->GetZ_asInt() != fragmentP->GetA_asInt())
616     products = theExcitationHandler->BreakItUp(*fragmentP);
617     // else
618     //   products = theExcitationHandlerx->BreakItUp(*fragmentP);      
619     delete fragmentP;
620     fragmentP = nullptr;
621   
622     G4ReactionProductVector::iterator iter;
623     for (iter = products->begin(); iter != products->end(); ++iter)
624     {
625       G4DynamicParticle *secondary =
626         new G4DynamicParticle((*iter)->GetDefinition(),
627         (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
628       theParticleChange.AddSecondary (secondary, secID);
629       G4String particleName = (*iter)->GetDefinition()->GetParticleName();
630       delete (*iter); // get rid of leftover particle def!
631       if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size())
632       {
633         G4cout <<"------------------------" <<G4endl;
634         G4cout <<"The projectile fragment:" <<G4endl;
635         G4cout <<"------------------------" <<G4endl;
636         G4cout <<" fragmentP = " <<particleName
637                <<" Energy    = " <<secondary->GetKineticEnergy()
638                <<G4endl;
639       }
640     }
641     delete products;
642   }
643 //
644 //
645 // Now decay the target nucleus - no boost is applied since in this
646 // approximation it is assumed that there is negligible momentum transfer from
647 // the projectile.
648 //
649   if (fragmentT != nullptr)
650   {
651     G4ReactionProductVector *products = nullptr;
652     //    if (fragmentT->GetZ_asInt() != fragmentT->GetA_asInt())
653       products = theExcitationHandler->BreakItUp(*fragmentT);
654     // else
655     //   products = theExcitationHandlerx->BreakItUp(*fragmentT);      
656     // delete fragmentT;
657     fragmentT = nullptr;
658   
659     G4ReactionProductVector::iterator iter;
660     for (iter = products->begin(); iter != products->end(); ++iter)
661     {
662       G4DynamicParticle *secondary =
663         new G4DynamicParticle((*iter)->GetDefinition(),
664         (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
665       theParticleChange.AddSecondary (secondary, secID);
666       G4String particleName = (*iter)->GetDefinition()->GetParticleName();
667       delete (*iter); // get rid of leftover particle def!
668       if (verboseLevel >= 2 && particleName.find("[",0) < particleName.size())
669       {
670         G4cout <<"--------------------" <<G4endl;
671         G4cout <<"The target fragment:" <<G4endl;
672         G4cout <<"--------------------" <<G4endl;
673         G4cout <<" fragmentT = " <<particleName
674                <<" Energy    = " <<secondary->GetKineticEnergy()
675                <<G4endl;
676       }
677     }
678     delete products;
679   }
680 
681   if (verboseLevel >= 2)
682      G4cout <<"########################################"
683             <<"########################################"
684             <<G4endl;
685   
686   delete theAbrasionGeometry;
687   return &theParticleChange;
688 }
689 ////////////////////////////////////////////////////////////////////////////////
690 //
691 G4Fragment *G4WilsonAbrasionModel::GetAbradedNucleons (G4int Dabr, G4double A,
692   G4double Z, G4double r)
693 {
694 //
695 //
696 // Initialise variables.  tau is the Fermi radius of the nucleus.  The variables
697 // p..., C... and gamma are used to help sample the secondary nucleon
698 // spectrum.
699 //
700   
701   G4double pK   = hbarc * G4Pow::GetInstance()->A13(9.0 * pi / 4.0 * A) / (1.29 * r);
702   if (A <= 24.0) pK *= -0.229*G4Pow::GetInstance()->A13(A) + 1.62;
703   G4double pKsq = pK * pK;
704   G4double p1sq = 2.0/5.0 * pKsq;
705   G4double p2sq = 6.0/5.0 * pKsq;
706   G4double p3sq = 500.0 * 500.0;
707   G4double C1   = 1.0;
708   G4double C2   = 0.03;
709   G4double C3   = 0.0002;
710   G4double gamma = 90.0 * MeV;
711   G4double maxn = C1 + C2 + C3;
712 //
713 //
714 // initialise the number of secondary nucleons abraded to zero, and initially set
715 // the type of nucleon abraded to proton ... just for now.
716 //  
717   G4double Aabr                     = 0.0;
718   G4double Zabr                     = 0.0; 
719   G4ParticleDefinition *typeNucleon = G4Proton::ProtonDefinition();
720   G4DynamicParticle *dynamicNucleon = nullptr;
721   G4ParticleMomentum pabr(0.0, 0.0, 0.0);
722 //
723 //
724 // Now go through each abraded nucleon and sample type, spectrum and angle.
725 //
726   G4bool isForLoopExitAnticipated = false;
727   for (G4int i=0; i<Dabr; ++i)
728   {
729 //
730 //
731 // Sample the nucleon momentum distribution by simple rejection techniques.  We
732 // reject values of p == 0.0 since this causes bad behaviour in the sinh term.
733 //
734     G4double p   = 0.0;
735     G4bool found = false;
736     const G4int maxNumberOfLoops = 100000;
737     G4int loopCounter = -1;
738     while (!found && ++loopCounter < maxNumberOfLoops)  /* Loop checking, 07.08.2015, A.Ribon */
739     {
740       while (p <= 0.0) p = npK * pK * G4UniformRand();  /* Loop checking, 07.08.2015, A.Ribon */
741       G4double psq = p * p;
742       found = maxn * G4UniformRand() < C1*G4Exp(-psq/p1sq/2.0) +
743         C2*G4Exp(-psq/p2sq/2.0) + C3*G4Exp(-psq/p3sq/2.0) + p/gamma/(0.5*(G4Exp(p/gamma)-G4Exp(-p/gamma)));
744     }
745     if ( loopCounter >= maxNumberOfLoops )
746     {
747       isForLoopExitAnticipated = true;
748       break;
749     }
750 //
751 //
752 // Determine the type of particle abraded.  Can only be proton or neutron,
753 // and the probability is determine to be proportional to the ratio as found
754 // in the nucleus at each stage.
755 //
756     G4double prob = (Z-Zabr)/(A-Aabr);
757     if (G4UniformRand()<prob)
758     {
759       Zabr++;
760       typeNucleon = G4Proton::ProtonDefinition();
761     }
762     else
763       typeNucleon = G4Neutron::NeutronDefinition();
764     Aabr++;
765 //
766 //
767 // The angular distribution of the secondary nucleons is approximated to an
768 // isotropic distribution in the rest frame of the nucleus (this will be Lorentz
769 // boosted later.
770 //
771     G4double costheta = 2.*G4UniformRand()-1.0;
772     G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
773     G4double phi      = 2.0*pi*G4UniformRand()*rad;
774     G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
775     G4double nucleonMass = typeNucleon->GetPDGMass();
776     G4double E           = std::sqrt(p*p + nucleonMass*nucleonMass)-nucleonMass;
777     dynamicNucleon = new G4DynamicParticle(typeNucleon,direction,E);
778     theParticleChange.AddSecondary (dynamicNucleon, secID);
779     pabr += p*direction;
780   }
781 //
782 //
783 // Next determine the details of the nuclear prefragment .. that is if there
784 // is one or more protons in the residue.  (Note that the 1 eV in the total
785 // energy is a safety factor to avoid any possibility of negative rest mass
786 // energy.)
787 //
788   G4Fragment *fragment = nullptr;
789   if ( ! isForLoopExitAnticipated && Z-Zabr>=1.0 )
790   {
791     G4double ionMass = G4ParticleTable::GetParticleTable()->GetIonTable()->
792                        GetIonMass(G4lrint(Z-Zabr),G4lrint(A-Aabr));
793     G4double E       = std::sqrt(pabr.mag2() + ionMass*ionMass);
794     G4LorentzVector lorentzVector = G4LorentzVector(-pabr, E + 1.0*eV);
795     fragment =
796       new G4Fragment((G4int) (A-Aabr), (G4int) (Z-Zabr), lorentzVector);
797   }
798 
799   return fragment;
800 }
801 ////////////////////////////////////////////////////////////////////////////////
802 //
803 G4double G4WilsonAbrasionModel::GetNucleonInducedExcitation
804   (G4double rP, G4double rT, G4double r)
805 {
806 //
807 //
808 // Initialise variables.
809 //
810   G4double Cl   = 0.0;
811   G4double rPsq = rP * rP;
812   G4double rTsq = rT * rT;
813   G4double rsq  = r * r;
814 //
815 //
816 // Depending upon the impact parameter, a different form of the chord length is
817 // is used.
818 //  
819   if (r > rT) Cl = 2.0*std::sqrt(rPsq + 2.0*r*rT - rsq - rTsq);
820   else        Cl = 2.0*rP;
821 //
822 //
823 // The next lines have been changed to include a "catch" to make sure if the
824 // projectile and target are too close, Ct is set to twice rP or twice rT.
825 // Otherwise the standard Wilson algorithm should work fine.
826 // PRT 20091023.
827 //
828   G4double Ct = 0.0;
829   if      (rT > rP && rsq < rTsq - rPsq) Ct = 2.0 * rP;
830   else if (rP > rT && rsq < rPsq - rTsq) Ct = 2.0 * rT;
831   else {
832     G4double bP = (rPsq+rsq-rTsq)/2.0/r;
833     G4double x  = rPsq - bP*bP;
834     if (x < 0.0) {
835       G4cerr <<"########################################"
836              <<"########################################"
837              <<G4endl;
838       G4cerr <<"ERROR IN G4WilsonAbrasionModel::GetNucleonInducedExcitation"
839              <<G4endl;
840       G4cerr <<"rPsq - bP*bP < 0.0 and cannot be square-rooted" <<G4endl;
841       G4cerr <<"Set to zero instead" <<G4endl;
842       G4cerr <<"########################################"
843              <<"########################################"
844              <<G4endl;
845     }
846     Ct = 2.0*std::sqrt(x);
847   }
848   
849   G4double Ex = 13.0 * Cl / fermi;
850   if (Ct > 1.5*fermi)
851     Ex += 13.0 * Cl / fermi /3.0 * (Ct/fermi - 1.5);
852 
853   return Ex;
854 }
855 ////////////////////////////////////////////////////////////////////////////////
856 //
857 void G4WilsonAbrasionModel::SetUseAblation (G4bool useAblation1)
858 {
859   if (useAblation != useAblation1)
860   {
861     useAblation = useAblation1;
862     if (useAblation)
863     {
864       theAblation = new G4WilsonAblationModel;
865       theAblation->SetVerboseLevel(verboseLevel);
866       theExcitationHandler->SetEvaporation(theAblation);
867     }
868     else
869     {
870       delete theExcitationHandler;
871       theAblation                      = nullptr;
872       theExcitationHandler  = new G4ExcitationHandler();
873     }
874   }
875   return; 
876 }
877 ////////////////////////////////////////////////////////////////////////////////
878 //
879 void G4WilsonAbrasionModel::PrintWelcomeMessage ()
880 {
881   G4cout <<G4endl;
882   G4cout <<" *****************************************************************"
883          <<G4endl;
884   G4cout <<" Nuclear abrasion model for nuclear-nuclear interactions activated"
885          <<G4endl;
886   G4cout <<" (Written by QinetiQ Ltd for the European Space Agency)"
887          <<G4endl;
888   G4cout <<" *****************************************************************"
889          <<G4endl;
890   G4cout << G4endl;
891 
892   return;
893 }
894 ////////////////////////////////////////////////////////////////////////////////
895 //
896