Geant4  10.01.p03
G4EMDissociation.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: G4EMDissociation.cc
39 //
40 // Version: B.1
41 // Date: 15/04/04
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 // 17 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 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
60 //
61 #include "G4EMDissociation.hh"
62 #include "G4PhysicalConstants.hh"
63 #include "G4SystemOfUnits.hh"
64 #include "G4Evaporation.hh"
65 #include "G4FermiBreakUp.hh"
66 #include "G4StatMF.hh"
67 #include "G4ParticleDefinition.hh"
68 #include "G4LorentzVector.hh"
69 #include "G4PhysicsFreeVector.hh"
71 #include "G4Proton.hh"
72 #include "G4Neutron.hh"
73 #include "G4IonTable.hh"
75 #include "G4DecayProducts.hh"
76 #include "G4DynamicParticle.hh"
77 #include "G4Fragment.hh"
79 #include "Randomize.hh"
80 #include "globals.hh"
81 
83 
84  // Send message to stdout to advise that the G4EMDissociation model is being
85  // used.
87 
88  // No de-excitation handler has been supplied - define the default handler.
90  G4Evaporation* theEvaporation = new G4Evaporation;
91  G4FermiBreakUp* theFermiBreakUp = new G4FermiBreakUp;
92  G4StatMF* theMF = new G4StatMF;
93  theExcitationHandler->SetEvaporation(theEvaporation);
94  theExcitationHandler->SetFermiModel(theFermiBreakUp);
99 
100  // This EM dissociation model needs access to the cross-sections held in
101  // G4EMDissociationCrossSection.
104 
105  // Set the minimum and maximum range for the model (despite nomanclature, this
106  // is in energy per nucleon number).
107  SetMinEnergy(100.0*MeV);
108  SetMaxEnergy(500.0*GeV);
109 
110  // Set the default verbose level to 0 - no output.
111  verboseLevel = 0;
112 }
113 
114 /*
115 G4EMDissociation::G4EMDissociation(const G4EMDissociation& emd)
116  : G4HadronicInteraction(emd)
117 {
118  if (emd.theExcitationHandler != 0) {
119  theExcitationHandler = new G4ExcitationHandler;
120  *theExcitationHandler = *emd.theExcitationHandler;
121  }
122 
123  handlerDefinedInternally = emd.handlerDefinedInternally;
124 
125  if (emd.dissociationCrossSection != 0) {
126  dissociationCrossSection = new G4EMDissociationCrossSection;
127  *dissociationCrossSection = *emd.dissociationCrossSection;
128  }
129 
130  if (emd.thePhotonSpectrum !- 0) {
131  thePhotonSpectrum = new G4EMDissociationSpectrum;
132  *thePhotonSpectrum = *emd.thePhotonSpectrum;
133 }
134 */
135 
137 {
138  // Send message to stdout to advise that the G4EMDissociation model is being
139  // used.
141 
142  theExcitationHandler = aExcitationHandler;
143  handlerDefinedInternally = false;
144 
145  // This EM dissociation model needs access to the cross-sections held in
146  // G4EMDissociationCrossSection.
149 
150  // Set the minimum and maximum range for the model (despite nomanclature, this
151  // is in energy per nucleon number)
152  SetMinEnergy(100.0*MeV);
153  SetMaxEnergy(500.0*GeV);
154  verboseLevel = 0;
155 }
156 
157 
160  // delete dissociationCrossSection;
161  // Cross section deleted by G4CrossSectionRegistry; don't do it here
162  // Bug reported by Gong Ding in Bug Report #1339
163  delete thePhotonSpectrum;
164 }
165 
166 
168  (const G4HadProjectile &theTrack, G4Nucleus &theTarget)
169 {
170  // The secondaries will be returned in G4HadFinalState &theParticleChange -
171  // initialise this.
172 
173  theParticleChange.Clear();
174  theParticleChange.SetStatusChange(stopAndKill);
175 
176  // Get relevant information about the projectile and target (A, Z) and
177  // energy/nuc, momentum, velocity, Lorentz factor and rest-mass of the
178  // projectile.
179 
180  const G4ParticleDefinition *definitionP = theTrack.GetDefinition();
181  const G4double AP = definitionP->GetBaryonNumber();
182  const G4double ZP = definitionP->GetPDGCharge();
183  G4LorentzVector pP = theTrack.Get4Momentum();
184  G4double E = theTrack.GetKineticEnergy()/AP;
185  G4double MP = theTrack.GetTotalEnergy() - E*AP;
186  G4double b = pP.beta();
187  G4double AT = theTarget.GetA_asInt();
188  G4double ZT = theTarget.GetZ_asInt();
190 
191  // Depending upon the verbosity level, output the initial information on the
192  // projectile and target
193  if (verboseLevel >= 2) {
194  G4cout.precision(6);
195  G4cout <<"########################################"
196  <<"########################################"
197  <<G4endl;
198  G4cout <<"IN G4EMDissociation" <<G4endl;
199  G4cout <<"Initial projectile A=" <<AP
200  <<", Z=" <<ZP
201  <<G4endl;
202  G4cout <<"Initial target A=" <<AT
203  <<", Z=" <<ZT
204  <<G4endl;
205  G4cout <<"Projectile momentum and Energy/nuc = " <<pP <<" ," <<E <<G4endl;
206  }
207 
208  // Initialise the variables which will be used with the phase-space decay and
209  // to boost the secondaries from the interaction.
210 
211  G4ParticleDefinition *typeNucleon = NULL;
212  G4ParticleDefinition *typeDaughter = NULL;
213  G4double Eg = 0.0;
214  G4double mass = 0.0;
215  G4ThreeVector boost = G4ThreeVector(0.0, 0.0, 0.0);
216 
217  // Determine the cross-sections at the giant dipole and giant quadrupole
218  // resonance energies for the projectile and then target. The information is
219  // initially provided in the G4PhysicsFreeVector individually for the E1
220  // and E2 fields. These are then summed.
221 
222  G4double bmin = thePhotonSpectrum->GetClosestApproach(AP, ZP, AT, ZT, b);
223  G4PhysicsFreeVector *crossSectionP = dissociationCrossSection->
224  GetCrossSectionForProjectile(AP, ZP, AT, ZT, b, bmin);
225  G4PhysicsFreeVector *crossSectionT = dissociationCrossSection->
226  GetCrossSectionForTarget(AP, ZP, AT, ZT, b, bmin);
227 
228  G4double totCrossSectionP = (*crossSectionP)[0]+(*crossSectionP)[1];
229  G4double totCrossSectionT = (*crossSectionT)[0]+(*crossSectionT)[1];
230 
231  // Now sample whether the interaction involved EM dissociation of the projectile
232  // or the target.
233 
234  if (G4UniformRand() <
235  totCrossSectionP / (totCrossSectionP + totCrossSectionT)) {
236 
237  // It was the projectile which underwent EM dissociation. Define the Lorentz
238  // boost to be applied to the secondaries, and sample whether a proton or a
239  // neutron was ejected. Then determine the energy of the virtual gamma ray
240  // which passed from the target nucleus ... this will be used to define the
241  // excitation of the projectile.
242 
243  mass = MP;
244  if (G4UniformRand() < dissociationCrossSection->
245  GetWilsonProbabilityForProtonDissociation (AP, ZP))
246  {
247  if (verboseLevel >= 2)
248  G4cout <<"Projectile underwent EM dissociation producing a proton"
249  <<G4endl;
250  typeNucleon = G4Proton::ProtonDefinition();
251  typeDaughter = G4IonTable::GetIonTable()->
252  GetIon((G4int) ZP-1, (G4int) AP-1, 0.0);
253  }
254  else
255  {
256  if (verboseLevel >= 2)
257  G4cout <<"Projectile underwent EM dissociation producing a neutron"
258  <<G4endl;
259  typeNucleon = G4Neutron::NeutronDefinition();
260  typeDaughter = G4IonTable::GetIonTable()->
261  GetIon((G4int) ZP, (G4int) AP-1, 0.0);
262  }
263  if (G4UniformRand() < (*crossSectionP)[0]/totCrossSectionP)
264  {
265  Eg = crossSectionP->GetLowEdgeEnergy(0);
266  if (verboseLevel >= 2)
267  G4cout <<"Transition type was E1" <<G4endl;
268  }
269  else
270  {
271  Eg = crossSectionP->GetLowEdgeEnergy(1);
272  if (verboseLevel >= 2)
273  G4cout <<"Transition type was E2" <<G4endl;
274  }
275 
276  // We need to define a Lorentz vector with the original momentum, but total
277  // energy includes the projectile and virtual gamma. This is then used
278  // to calculate the boost required for the secondaries.
279 
280  pP.setE(pP.e()+Eg);
281  boost = pP.findBoostToCM();
282  }
283  else
284  {
285  // It was the target which underwent EM dissociation. Sample whether a
286  // proton or a neutron was ejected. Then determine the energy of the virtual
287  // gamma ray which passed from the projectile nucleus ... this will be used to
288  // define the excitation of the target.
289 
290  mass = MT;
291  if (G4UniformRand() < dissociationCrossSection->
292  GetWilsonProbabilityForProtonDissociation (AT, ZT))
293  {
294  if (verboseLevel >= 2)
295  G4cout <<"Target underwent EM dissociation producing a proton"
296  <<G4endl;
297  typeNucleon = G4Proton::ProtonDefinition();
298  typeDaughter = G4IonTable::GetIonTable()->
299  GetIon((G4int) ZT-1, (G4int) AT-1, 0.0);
300  }
301  else
302  {
303  if (verboseLevel >= 2)
304  G4cout <<"Target underwent EM dissociation producing a neutron"
305  <<G4endl;
306  typeNucleon = G4Neutron::NeutronDefinition();
307  typeDaughter = G4IonTable::GetIonTable()->
308  GetIon((G4int) ZT, (G4int) AT-1, 0.0);
309  }
310  if (G4UniformRand() < (*crossSectionT)[0]/totCrossSectionT)
311  {
312  Eg = crossSectionT->GetLowEdgeEnergy(0);
313  if (verboseLevel >= 2)
314  G4cout <<"Transition type was E1" <<G4endl;
315  }
316  else
317  {
318  Eg = crossSectionT->GetLowEdgeEnergy(1);
319  if (verboseLevel >= 2)
320  G4cout <<"Transition type was E2" <<G4endl;
321  }
322 
323  // Add the projectile to theParticleChange, less the energy of the
324  // not-so-virtual gamma-ray. Not that at the moment, no lateral momentum
325  // is transferred between the projectile and target nuclei.
326 
327  G4ThreeVector v = pP.vect();
328  v.setMag(1.0);
329  G4DynamicParticle *changedP = new G4DynamicParticle (definitionP, v, E*AP-Eg);
330  theParticleChange.AddSecondary (changedP);
331  if (verboseLevel >= 2)
332  {
333  G4cout <<"Projectile change:" <<G4endl;
334  changedP->DumpInfo();
335  }
336  }
337 
338  // Perform a two-body decay based on the restmass energy of the parent and
339  // gamma-ray, and the masses of the daughters. In the frame of reference of
340  // the nucles, the angular distribution is sampled isotropically, but the
341  // the nucleon and secondary nucleus are boosted if they've come from the
342  // projectile.
343 
344  G4double e = mass + Eg;
345  G4double mass1 = typeNucleon->GetPDGMass();
346  G4double mass2 = typeDaughter->GetPDGMass();
347  G4double pp = (e+mass1+mass2)*(e+mass1-mass2)*
348  (e-mass1+mass2)*(e-mass1-mass2)/(4.0*e*e);
349  if (pp < 0.0) {
350  pp = 1.0*eV;
351 // if (verboseLevel >`= 1)
352 // {
353 // G4cout <<"IN G4EMDissociation::ApplyYoursef" <<G4endl;
354 // G4cout <<"Error in mass of secondaries compared with primary:" <<G4endl;
355 // G4cout <<"Rest mass of primary = " <<mass <<" MeV" <<G4endl;
356 // G4cout <<"Virtual gamma energy = " <<Eg <<" MeV" <<G4endl;
357 // G4cout <<"Rest mass of secondary #1 = " <<mass1 <<" MeV" <<G4endl;
358 // G4cout <<"Rest mass of secondary #2 = " <<mass2 <<" MeV" <<G4endl;
359 // }
360  }
361  else
362  pp = std::sqrt(pp);
363  G4double costheta = 2.*G4UniformRand()-1.0;
364  G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
365  G4double phi = 2.0*pi*G4UniformRand()*rad;
366  G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
367  G4DynamicParticle *dynamicNucleon =
368  new G4DynamicParticle(typeNucleon, direction*pp);
369  dynamicNucleon->Set4Momentum(dynamicNucleon->Get4Momentum().boost(-boost));
370  G4DynamicParticle *dynamicDaughter =
371  new G4DynamicParticle(typeDaughter, -direction*pp);
372  dynamicDaughter->Set4Momentum(dynamicDaughter->Get4Momentum().boost(-boost));
373 
374  // The "decay" products have to be transferred to the G4HadFinalState object.
375  // Furthermore, the residual nucleus should be de-excited.
376 
377  theParticleChange.AddSecondary (dynamicNucleon);
378  if (verboseLevel >= 2) {
379  G4cout <<"Nucleon from the EMD process:" <<G4endl;
380  dynamicNucleon->DumpInfo();
381  }
382 
383  G4Fragment* theFragment = new
384  G4Fragment((G4int) typeDaughter->GetBaryonNumber(),
385  (G4int) typeDaughter->GetPDGCharge(), dynamicDaughter->Get4Momentum());
386 
387  if (verboseLevel >= 2) {
388  G4cout <<"Dynamic properties of the prefragment:" <<G4endl;
389  G4cout.precision(6);
390  dynamicDaughter->DumpInfo();
391  G4cout <<"Nuclear properties of the prefragment:" <<G4endl;
392  G4cout <<theFragment <<G4endl;
393  }
394 
395  G4ReactionProductVector* products =
396  theExcitationHandler->BreakItUp(*theFragment);
397  delete theFragment;
398  theFragment = NULL;
399 
400  G4DynamicParticle* secondary = 0;
401  G4ReactionProductVector::iterator iter;
402  for (iter = products->begin(); iter != products->end(); ++iter) {
403  secondary = new G4DynamicParticle((*iter)->GetDefinition(),
404  (*iter)->GetTotalEnergy(), (*iter)->GetMomentum());
405  theParticleChange.AddSecondary (secondary);
406  }
407  delete products;
408 
409  if (verboseLevel >= 2)
410  G4cout <<"########################################"
411  <<"########################################"
412  <<G4endl;
413 
414  return &theParticleChange;
415 }
416 
417 
419 {
420  G4cout <<G4endl;
421  G4cout <<" ****************************************************************"
422  <<G4endl;
423  G4cout <<" EM dissociation model for nuclear-nuclear interactions activated"
424  <<G4endl;
425  G4cout <<" (Written by QinetiQ Ltd for the European Space Agency)"
426  <<G4endl;
427  G4cout <<" ****************************************************************"
428  <<G4endl;
429  G4cout << G4endl;
430 
431  return;
432 }
433 
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4bool handlerDefinedInternally
static const double MeV
Definition: G4SIunits.hh:193
static G4double GetNuclearMass(const G4double A, const G4double Z)
CLHEP::Hep3Vector G4ThreeVector
const G4double pi
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
void DumpInfo(G4int mode=0) const
void SetMinEForMultiFrag(G4double anE)
int G4int
Definition: G4Types.hh:78
G4ExcitationHandler * theExcitationHandler
std::vector< G4ReactionProduct * > G4ReactionProductVector
void SetMinEnergy(G4double anEnergy)
G4EMDissociationCrossSection * dissociationCrossSection
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void SetFermiModel(G4VFermiBreakUp *ptr)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &, G4Nucleus &)
G4ErrorTarget * theTarget
Definition: errprop.cc:59
void SetMultiFragmentation(G4VMultiFragmentation *ptr)
static const double GeV
Definition: G4SIunits.hh:196
const G4LorentzVector & Get4Momentum() const
G4LorentzVector Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static const double rad
Definition: G4SIunits.hh:130
void Set4Momentum(const G4LorentzVector &momentum)
static const double eV
Definition: G4SIunits.hh:194
G4double GetPDGMass() const
void SetMaxAandZForFermiBreakUp(G4int anA, G4int aZ)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4EMDissociationSpectrum * thePhotonSpectrum
void SetEvaporation(G4VEvaporation *ptr)
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:99
G4double GetTotalEnergy() const
CLHEP::HepLorentzVector G4LorentzVector