Geant4  10.01.p02
G4WilsonAblationModel.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 // * *
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: G4WilsonAblationModel.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 // 08 December 2009, P R Truscott, QinetiQ Ltd, UK
59 // Ver 1.0
60 // Updated as a result of changes in the G4Evaporation classes. These changes
61 // affect mostly SelectSecondariesByEvaporation, and now you have variables
62 // associated with the evaporation model which can be changed:
63 // OPTxs to select the inverse cross-section
64 // OPTxs = 0 => Dostrovski's parameterization
65 // OPTxs = 1 or 2 => Chatterjee's paramaterization
66 // OPTxs = 3 or 4 => Kalbach's parameterization
67 // useSICB => use superimposed Coulomb Barrier for inverse cross
68 // sections
69 // Other problem found with G4Fragment definition using Lorentz vector and
70 // **G4ParticleDefinition**. This does not allow A and Z to be defined for the
71 // fragment for some reason. Now the fragment is defined more explicitly:
72 // G4Fragment *fragment = new G4Fragment(A, Z, lorentzVector);
73 // to avoid this quirk. Bug found in SelectSecondariesByDefault: *type is now
74 // equated to evapType[i] whereas previously it was equated to fragType[i].
75 //
76 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
78 //
79 #include <iomanip>
80 #include <numeric>
81 
82 #include "G4WilsonAblationModel.hh"
83 #include "G4PhysicalConstants.hh"
84 #include "G4SystemOfUnits.hh"
85 #include "Randomize.hh"
86 #include "G4ParticleTable.hh"
87 #include "G4IonTable.hh"
88 #include "G4Alpha.hh"
89 #include "G4He3.hh"
90 #include "G4Triton.hh"
91 #include "G4Deuteron.hh"
92 #include "G4Proton.hh"
93 #include "G4Neutron.hh"
100 #include "G4LorentzVector.hh"
101 #include "G4VEvaporationChannel.hh"
102 
104 //
106 {
107 //
108 //
109 // Send message to stdout to advise that the G4Abrasion model is being used.
110 //
112 //
113 //
114 // Set the default verbose level to 0 - no output.
115 //
116  verboseLevel = 0;
117 //
118 //
119 // Set the binding energy per nucleon .... did I mention that this is a crude
120 // model for nuclear de-excitation?
121 //
122  B = 10.0 * MeV;
123 //
124 //
125 // It is possuble to switch off secondary particle production (other than the
126 // final nuclear fragment). The default is on.
127 //
128  produceSecondaries = true;
129 //
130 //
131 // Now we need to define the decay modes. We're using the G4Evaporation model
132 // to help determine the kinematics of the decay.
133 //
134  nFragTypes = 6;
135  fragType[0] = G4Alpha::Alpha();
136  fragType[1] = G4He3::He3();
137  fragType[2] = G4Triton::Triton();
139  fragType[4] = G4Proton::Proton();
141 //
142 //
143 // Set verboseLevel default to no output.
144 //
145  verboseLevel = 0;
148 //
149 //
150 // Set defaults for evaporation classes. These can be overridden by user
151 // "set" methods.
152 //
153  OPTxs = 3;
154  useSICB = false;
155  fragmentVector = 0;
156 }
158 //
160 {
161  if (theChannels != 0) theChannels = 0;
162  if (theChannelFactory != 0) delete theChannelFactory;
163 }
165 //
167  (const G4Fragment &theNucleus)
168 {
169 //
170 //
171 // Initilise the pointer to the G4FragmentVector used to return the information
172 // about the breakup.
173 //
174  fragmentVector = new G4FragmentVector;
175  fragmentVector->clear();
176 //
177 //
178 // Get the A, Z and excitation of the nucleus.
179 //
180  G4int A = theNucleus.GetA_asInt();
181  G4int Z = theNucleus.GetZ_asInt();
182  G4double ex = theNucleus.GetExcitationEnergy();
183  if (verboseLevel >= 2)
184  {
185  G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
186  <<"oooooooooooooooooooooooooooooooooooooooo"
187  <<G4endl;
188  G4cout.precision(6);
189  G4cout <<"IN G4WilsonAblationModel" <<G4endl;
190  G4cout <<"Initial prefragment A=" <<A
191  <<", Z=" <<Z
192  <<", excitation energy = " <<ex/MeV <<" MeV"
193  <<G4endl;
194  }
195 //
196 //
197 // Check that there is a nucleus to speak of. It's possible there isn't one
198 // or its just a proton or neutron. In either case, the excitation energy
199 // (from the Lorentz vector) is not used.
200 //
201  if (A == 0)
202  {
203  if (verboseLevel >= 2)
204  {
205  G4cout <<"No nucleus to decay" <<G4endl;
206  G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
207  <<"oooooooooooooooooooooooooooooooooooooooo"
208  <<G4endl;
209  }
210  return fragmentVector;
211  }
212  else if (A == 1)
213  {
214  G4LorentzVector lorentzVector = theNucleus.GetMomentum();
215  lorentzVector.setE(lorentzVector.e()-ex+10.0*eV);
216  if (Z == 0)
217  {
218  G4Fragment *fragment = new G4Fragment(lorentzVector,G4Neutron::Neutron());
219  fragmentVector->push_back(fragment);
220  }
221  else
222  {
223  G4Fragment *fragment = new G4Fragment(lorentzVector,G4Proton::Proton());
224  fragmentVector->push_back(fragment);
225  }
226  if (verboseLevel >= 2)
227  {
228  G4cout <<"Final fragment is in fact only a nucleon) :" <<G4endl;
229  G4cout <<(*fragmentVector)[0] <<G4endl;
230  G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
231  <<"oooooooooooooooooooooooooooooooooooooooo"
232  <<G4endl;
233  }
234  return fragmentVector;
235  }
236 //
237 //
238 // Then the number of nucleons ablated (either as nucleons or light nuclear
239 // fragments) is based on a simple argument for the binding energy per nucleon.
240 //
241  G4int DAabl = (G4int) (ex / B);
242  if (DAabl > A) DAabl = A;
243 // The following lines are no longer accurate given we now treat the final fragment
244 // if (verboseLevel >= 2)
245 // G4cout <<"Number of nucleons ejected = " <<DAabl <<G4endl;
246 
247 //
248 //
249 // Determine the nuclear fragment from the ablation process by sampling the
250 // Rudstam equation.
251 //
252  G4int AF = A - DAabl;
253  G4int ZF = 0;
254  if (AF > 0)
255  {
256  G4double AFd = (G4double) AF;
257  G4double R = 11.8 / std::pow(AFd, 0.45);
258  G4int minZ = Z - DAabl;
259  if (minZ <= 0) minZ = 1;
260 //
261 //
262 // Here we define an integral probability distribution based on the Rudstam
263 // equation assuming a constant AF.
264 //
265  G4double sig[100];
266  G4double sum = 0.0;
267  for (G4int ii=minZ; ii<= Z; ii++)
268  {
269  sum += std::exp(-R*std::pow(std::abs(ii - 0.486*AFd + 3.8E-04*AFd*AFd),1.5));
270  sig[ii] = sum;
271  }
272 //
273 //
274 // Now sample that distribution to determine a value for ZF.
275 //
276  G4double xi = G4UniformRand();
277  G4int iz = minZ;
278  G4bool found = false;
279  while (iz <= Z && !found)
280  {
281  found = (xi <= sig[iz]/sum);
282  if (!found) iz++;
283  }
284  if (iz > Z)
285  ZF = Z;
286  else
287  ZF = iz;
288  }
289  G4int DZabl = Z - ZF;
290 //
291 //
292 // Now determine the nucleons or nuclei which have bee ablated. The preference
293 // is for the production of alphas, then other nuclei in order of decreasing
294 // binding energy. The energies assigned to the products of the decay are
295 // provisional for the moment (the 10eV is just to avoid errors with negative
296 // excitation energies due to rounding).
297 //
298  G4double totalEpost = 0.0;
299  evapType.clear();
300  for (G4int ift=0; ift<nFragTypes; ift++)
301  {
302  G4ParticleDefinition *type = fragType[ift];
303  G4double n = std::floor((G4double) DAabl / type->GetBaryonNumber() + 1.0E-10);
304  G4double n1 = 1.0E+10;
305  if (fragType[ift]->GetPDGCharge() > 0.0)
306  n1 = std::floor((G4double) DZabl / type->GetPDGCharge() + 1.0E-10);
307  if (n > n1) n = n1;
308  if (n > 0.0)
309  {
310  G4double mass = type->GetPDGMass();
311  for (G4int j=0; j<(G4int) n; j++)
312  {
313  totalEpost += mass;
314  evapType.push_back(type);
315  }
316  DAabl -= (G4int) (n * type->GetBaryonNumber() + 1.0E-10);
317  DZabl -= (G4int) (n * type->GetPDGCharge() + 1.0E-10);
318  }
319  }
320 //
321 //
322 // Determine the properties of the final nuclear fragment. Note that if
323 // the final fragment is predicted to have a nucleon number of zero, then
324 // really it's the particle last in the vector evapType which becomes the
325 // final fragment. Therefore delete this from the vector if this is the
326 // case.
327 //
328  G4double massFinalFrag = 0.0;
329  if (AF > 0)
330  massFinalFrag = G4ParticleTable::GetParticleTable()->GetIonTable()->
331  GetIonMass(ZF,AF);
332  else
333  {
334  G4ParticleDefinition *type = evapType[evapType.size()-1];
335  AF = type->GetBaryonNumber();
336  ZF = (G4int) (type->GetPDGCharge() + 1.0E-10);
337  evapType.erase(evapType.end()-1);
338  }
339  totalEpost += massFinalFrag;
340 //
341 //
342 // Provide verbose output on the nuclear fragment if requested.
343 //
344  if (verboseLevel >= 2)
345  {
346  G4cout <<"Final fragment A=" <<AF
347  <<", Z=" <<ZF
348  <<G4endl;
349  for (G4int ift=0; ift<nFragTypes; ift++)
350  {
351  G4ParticleDefinition *type = fragType[ift];
352  G4int n = std::count(evapType.begin(),evapType.end(),type);
353  if (n > 0)
354  G4cout <<"Particle type: " <<std::setw(10) <<type->GetParticleName()
355  <<", number of particles emitted = " <<n <<G4endl;
356  }
357  }
358 //
359 // Add the total energy from the fragment. Note that the fragment is assumed
360 // to be de-excited and does not undergo photo-evaporation .... I did mention
361 // this is a bit of a crude model?
362 //
363  G4double massPreFrag = theNucleus.GetGroundStateMass();
364  G4double totalEpre = massPreFrag + ex;
365  G4double excess = totalEpre - totalEpost;
366 // G4Fragment *resultNucleus(theNucleus);
367  G4Fragment *resultNucleus = new G4Fragment(A, Z, theNucleus.GetMomentum());
368  G4ThreeVector boost(0.0,0.0,0.0);
369  G4int nEvap = 0;
370  if (produceSecondaries && evapType.size()>0)
371  {
372  if (excess > 0.0)
373  {
374  SelectSecondariesByEvaporation (resultNucleus);
375  nEvap = fragmentVector->size();
376  boost = resultNucleus->GetMomentum().findBoostToCM();
377  if (evapType.size() > 0)
378  SelectSecondariesByDefault (boost);
379  }
380  else
381  SelectSecondariesByDefault(G4ThreeVector(0.0,0.0,0.0));
382  }
383 
384  if (AF > 0)
385  {
387  GetIonMass(ZF,AF);
388  G4double e = mass + 10.0*eV;
389  G4double p = std::sqrt(e*e-mass*mass);
390  G4ThreeVector direction(0.0,0.0,1.0);
391  G4LorentzVector lorentzVector = G4LorentzVector(direction*p, e);
392  lorentzVector.boost(-boost);
393  G4Fragment* frag = new G4Fragment(AF, ZF, lorentzVector);
394  fragmentVector->push_back(frag);
395  }
396  delete resultNucleus;
397 //
398 //
399 // Provide verbose output on the ablation products if requested.
400 //
401  if (verboseLevel >= 2)
402  {
403  if (nEvap > 0)
404  {
405  G4cout <<"----------------------" <<G4endl;
406  G4cout <<"Evaporated particles :" <<G4endl;
407  G4cout <<"----------------------" <<G4endl;
408  }
409  G4int ie = 0;
410  G4FragmentVector::iterator iter;
411  for (iter = fragmentVector->begin(); iter != fragmentVector->end(); iter++)
412  {
413  if (ie == nEvap)
414  {
415 // G4cout <<*iter <<G4endl;
416  G4cout <<"---------------------------------" <<G4endl;
417  G4cout <<"Particles from default emission :" <<G4endl;
418  G4cout <<"---------------------------------" <<G4endl;
419  }
420  G4cout <<*iter <<G4endl;
421  }
422  G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
423  <<"oooooooooooooooooooooooooooooooooooooooo"
424  <<G4endl;
425  }
426 
427  return fragmentVector;
428 }
430 //
432  (G4Fragment *intermediateNucleus)
433 {
434  G4Fragment theResidualNucleus = *intermediateNucleus;
435  G4bool evaporate = true;
436  while (evaporate && evapType.size() != 0)
437  {
438 //
439 //
440 // Here's the cheaky bit. We're hijacking the G4Evaporation model, in order to
441 // more accurately sample to kinematics, but the species of the nuclear
442 // fragments will be the ones of our choosing as above.
443 //
444  std::vector <G4VEvaporationChannel*> theChannels1;
445  theChannels1.clear();
446  std::vector <G4VEvaporationChannel*>::iterator i;
447  VectorOfFragmentTypes::iterator iter;
448  std::vector <VectorOfFragmentTypes::iterator> iters;
449  iters.clear();
450  iter = std::find(evapType.begin(), evapType.end(), G4Alpha::Alpha());
451  if (iter != evapType.end())
452  {
453  theChannels1.push_back(new G4AlphaEvaporationChannel);
454  i = theChannels1.end() - 1;
455  (*i)->SetOPTxs(OPTxs);
456  (*i)->UseSICB(useSICB);
457 // (*i)->Initialize(theResidualNucleus);
458  iters.push_back(iter);
459  }
460  iter = std::find(evapType.begin(), evapType.end(), G4He3::He3());
461  if (iter != evapType.end())
462  {
463  theChannels1.push_back(new G4He3EvaporationChannel);
464  i = theChannels1.end() - 1;
465  (*i)->SetOPTxs(OPTxs);
466  (*i)->UseSICB(useSICB);
467 // (*i)->Initialize(theResidualNucleus);
468  iters.push_back(iter);
469  }
470  iter = std::find(evapType.begin(), evapType.end(), G4Triton::Triton());
471  if (iter != evapType.end())
472  {
473  theChannels1.push_back(new G4TritonEvaporationChannel);
474  i = theChannels1.end() - 1;
475  (*i)->SetOPTxs(OPTxs);
476  (*i)->UseSICB(useSICB);
477 // (*i)->Initialize(theResidualNucleus);
478  iters.push_back(iter);
479  }
480  iter = std::find(evapType.begin(), evapType.end(), G4Deuteron::Deuteron());
481  if (iter != evapType.end())
482  {
483  theChannels1.push_back(new G4DeuteronEvaporationChannel);
484  i = theChannels1.end() - 1;
485  (*i)->SetOPTxs(OPTxs);
486  (*i)->UseSICB(useSICB);
487 // (*i)->Initialize(theResidualNucleus);
488  iters.push_back(iter);
489  }
490  iter = std::find(evapType.begin(), evapType.end(), G4Proton::Proton());
491  if (iter != evapType.end())
492  {
493  theChannels1.push_back(new G4ProtonEvaporationChannel);
494  i = theChannels1.end() - 1;
495  (*i)->SetOPTxs(OPTxs);
496  (*i)->UseSICB(useSICB);
497 // (*i)->Initialize(theResidualNucleus);
498  iters.push_back(iter);
499  }
500  iter = std::find(evapType.begin(), evapType.end(), G4Neutron::Neutron());
501  if (iter != evapType.end())
502  {
503  theChannels1.push_back(new G4NeutronEvaporationChannel);
504  i = theChannels1.end() - 1;
505  (*i)->SetOPTxs(OPTxs);
506  (*i)->UseSICB(useSICB);
507 // (*i)->Initialize(theResidualNucleus);
508  iters.push_back(iter);
509  }
510  G4int nChannels = theChannels1.size();
511 
512  G4double totalProb = 0.0;
513  G4int ich = 0;
514  G4double probEvapType[6] = {0.0};
515  std::vector<G4VEvaporationChannel*>::iterator iterEv;
516  for (iterEv=theChannels1.begin(); iterEv!=theChannels1.end(); iterEv++) {
517  totalProb += (*iterEv)->GetEmissionProbability(intermediateNucleus);
518  probEvapType[ich] = totalProb;
519  ++ich;
520  }
521  if (totalProb > 0.0) {
522 //
523 //
524 // The emission probability for at least one of the evaporation channels is
525 // positive, therefore work out which one should be selected and decay
526 // the nucleus.
527 //
528  G4double xi = totalProb*G4UniformRand();
529  G4int ii = 0;
530  for (ii=0; ii<nChannels; ii++) {
531  if (xi < probEvapType[ii]) { break; }
532  }
533  if (ii >= nChannels) { ii = nChannels - 1; }
534  G4FragmentVector *evaporationResult = theChannels1[ii]->
535  BreakUp(*intermediateNucleus);
536  fragmentVector->push_back((*evaporationResult)[0]);
537  *intermediateNucleus = *(*evaporationResult)[1];
538  //delete evaporationResult->back();
539  delete evaporationResult;
540  //evapType.erase(iters[ii]);
541  }
542  else
543  {
544 //
545 //
546 // Probability for further evaporation is nil so have to escape from this
547 // routine and set the energies of the secondaries to 10eV.
548 //
549  evaporate = false;
550  }
551  }
552 
553  return;
554 }
556 //
558 {
559  for (unsigned i=0; i<evapType.size(); i++)
560  {
561  G4ParticleDefinition *type = evapType[i];
562  G4double mass = type->GetPDGMass();
563  G4double e = mass + 10.0*eV;
564  G4double p = std::sqrt(e*e-mass*mass);
565  G4double costheta = 2.0*G4UniformRand() - 1.0;
566  G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
567  G4double phi = twopi * G4UniformRand() * rad;
568  G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
569  G4LorentzVector lorentzVector = G4LorentzVector(direction*p, e);
570  lorentzVector.boost(-boost);
571 // Possibility that the following line is not correctly carrying over A and Z
572 // from particle definition. Force values. PRT 03/12/2009.
573 // G4Fragment *fragment =
574 // new G4Fragment(lorentzVector, type);
575  G4int A = type->GetBaryonNumber();
576  G4int Z = (G4int) (type->GetPDGCharge() + 1.0E-10);
577  G4Fragment *fragment =
578  new G4Fragment(A, Z, lorentzVector);
579 
580  fragmentVector->push_back(fragment);
581  }
582 }
584 //
586 {
587  G4cout <<G4endl;
588  G4cout <<" *****************************************************************"
589  <<G4endl;
590  G4cout <<" Nuclear ablation model for nuclear-nuclear interactions activated"
591  <<G4endl;
592  G4cout <<" (Written by QinetiQ Ltd for the European Space Agency)"
593  <<G4endl;
594  G4cout <<" !!! WARNING: This model is not well validation and should not be used for accurate simulation !!!"
595  <<G4endl;
596  G4cout <<" *****************************************************************"
597  <<G4endl;
598  G4cout << G4endl;
599 
600  return;
601 }
603 //
virtual std::vector< G4VEvaporationChannel * > * GetChannel()=0
static const double MeV
Definition: G4SIunits.hh:193
void SelectSecondariesByDefault(G4ThreeVector)
VectorOfFragmentTypes evapType
CLHEP::Hep3Vector G4ThreeVector
G4ParticleDefinition * fragType[6]
void SelectSecondariesByEvaporation(G4Fragment *)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
std::vector< G4VEvaporationChannel * > * theChannels
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
G4int GetA_asInt() const
Definition: G4Fragment.hh:243
bool G4bool
Definition: G4Types.hh:79
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:276
G4double iz
Definition: TRTMaterials.hh:39
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:265
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4int n
static const G4double A[nN]
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static const double rad
Definition: G4SIunits.hh:130
static const double eV
Definition: G4SIunits.hh:194
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4int GetZ_asInt() const
Definition: G4Fragment.hh:248
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4VEvaporationFactory * theChannelFactory
static G4He3 * He3()
Definition: G4He3.cc:94
G4FragmentVector * fragmentVector
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:260
CLHEP::HepLorentzVector G4LorentzVector
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)