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