Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCLClusterDecay.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// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// INCL++ intra-nuclear cascade model
27// Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28// Davide Mancusi, CEA
29// Alain Boudard, CEA
30// Sylvie Leray, CEA
31// Joseph Cugnon, University of Liege
32//
33// INCL++ revision: v5.1.8
34//
35#define INCLXX_IN_GEANT4_MODE 1
36
37#include "globals.hh"
38
39/** \file G4INCLClusterDecay.cc
40 * \brief Static class for carrying out cluster decays
41 *
42 * \date 6th July 2011
43 * \author Davide Mancusi
44 */
45
46#include "G4INCLClusterDecay.hh"
49#include "G4INCLRandom.hh"
50// #include <cassert>
51#include <algorithm>
52
53namespace G4INCL {
54
56 ParticleList decayProducts;
57 recursiveDecay(c, &decayProducts);
58
59 // Correctly update the particle type
60 if(c->getA()==1) {
61// assert(c->getZ()==1 || c->getZ()==0);
62 if(c->getZ()==1)
63 c->setType(Proton);
64 else
65 c->setType(Neutron);
66 c->setTableMass();
67 }
68
69 return decayProducts;
70 }
71
72 void ClusterDecay::recursiveDecay(Cluster * const c, ParticleList *decayProducts) {
73 const G4int Z = c->getZ();
74 const G4int A = c->getA();
75// assert(c->getExcitationEnergy()>-1.e-5);
76 if(c->getExcitationEnergy()<0.)
78
81
82 switch(theDecayMode) {
83 default:
84 ERROR("Unrecognized cluster-decay mode: " << theDecayMode << std::endl
85 << c->print());
87 // For stable clusters, just return
88 return;
89 break;
93 // Two-body decays
94 twoBodyDecay(c, theDecayMode, decayProducts);
95 break;
98 // Three-body decays
99 threeBodyDecay(c, theDecayMode, decayProducts);
100 break;
103 // Phase-space decays
104 phaseSpaceDecay(c, theDecayMode, decayProducts);
105 break;
106 }
107
108 // Calls itself recursively in case the produced remnant is still unstable.
109 // Sneaky, isn't it.
110 recursiveDecay(c,decayProducts);
111
112 } else {
113 // The cluster is too large for our decay-mode table. Decompose it only
114 // if Z==0 || Z==A.
115 DEBUG("Cluster is outside the decay-mode table." << c->print() << std::endl);
116 if(Z==A) {
117 DEBUG("Z==A, will decompose it in free protons." << std::endl);
118 phaseSpaceDecay(c, ParticleTable::ProtonUnbound, decayProducts);
119 } else if(Z==0) {
120 DEBUG("Z==0, will decompose it in free neutrons." << std::endl);
121 phaseSpaceDecay(c, ParticleTable::NeutronUnbound, decayProducts);
122 }
123 }
124 }
125
126 void ClusterDecay::twoBodyDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts) {
127 Particle *decayParticle = 0;
128 const ThreeVector mom(0.0, 0.0, 0.0);
129 const ThreeVector pos = c->getPosition();
130
131 // Create the emitted particle
132 switch(theDecayMode) {
134 decayParticle = new Particle(Proton, mom, pos);
135 break;
137 decayParticle = new Particle(Neutron, mom, pos);
138 break;
140 decayParticle = new Cluster(2,4);
141 break;
142 default:
143 ERROR("Unrecognized cluster-decay mode in two-body decay: " << theDecayMode << std::endl
144 << c->print());
145 return;
146 }
147 decayParticle->makeParticipant();
148 decayParticle->setNumberOfDecays(1);
149 decayParticle->setPosition(c->getPosition());
150 decayParticle->setEmissionTime(c->getEmissionTime());
151 decayParticle->setTableMass();
152
153 // Save some variables of the mother cluster
154 G4double motherMass = c->getMass();
155 const ThreeVector velocity = -c->boostVector();
156
157 // Characteristics of the daughter particle
158 const G4int daughterZ = c->getZ() - decayParticle->getZ();
159 const G4int daughterA = c->getA() - decayParticle->getA();
160 const G4double daughterMass = ParticleTable::getTableMass(daughterA,daughterZ);
161
162 // The mother cluster becomes the daughter
163 c->setZ(daughterZ);
164 c->setA(daughterA);
165 c->setMass(daughterMass);
166 c->setExcitationEnergy(0.);
167
168 // Decay kinematics in the mother rest frame
169 const G4double decayMass = decayParticle->getMass();
170// assert(motherMass-daughterMass-decayMass>-1.e-5); // Q-value should be >0
171 G4double pCM = 0.;
172 if(motherMass-daughterMass-decayMass>0.)
173 pCM = KinematicsUtils::momentumInCM(motherMass, daughterMass, decayMass);
174 const ThreeVector momentum = Random::normVector(pCM);
175 c->setMomentum(momentum);
176 c->adjustEnergyFromMomentum();
177 decayParticle->setMomentum(-momentum);
178 decayParticle->adjustEnergyFromMomentum();
179
180 // Boost to the lab frame
181 decayParticle->boost(velocity);
182 c->boost(velocity);
183
184 // Add the decay particle to the list of decay products
185 decayProducts->push_back(decayParticle);
186 }
187
188 void ClusterDecay::threeBodyDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts) {
189 Particle *decayParticle1 = 0;
190 Particle *decayParticle2 = 0;
191 const ThreeVector mom(0.0, 0.0, 0.0);
192 const ThreeVector pos = c->getPosition();
193
194 // Create the emitted particles
195 switch(theDecayMode) {
197 decayParticle1 = new Particle(Proton, mom, pos);
198 decayParticle2 = new Particle(Proton, mom, pos);
199 break;
201 decayParticle1 = new Particle(Neutron, mom, pos);
202 decayParticle2 = new Particle(Neutron, mom, pos);
203 break;
204 default:
205 ERROR("Unrecognized cluster-decay mode in three-body decay: " << theDecayMode << std::endl
206 << c->print());
207 return;
208 }
209 decayParticle1->makeParticipant();
210 decayParticle2->makeParticipant();
211 decayParticle1->setNumberOfDecays(1);
212 decayParticle2->setNumberOfDecays(1);
213 decayParticle1->setTableMass();
214 decayParticle2->setTableMass();
215
216 // Save some variables of the mother cluster
217 const G4double motherMass = c->getMass();
218 const ThreeVector velocity = -c->boostVector();
219
220 // Masses and charges of the daughter particle and of the decay products
221 const G4int decayZ1 = decayParticle1->getZ();
222 const G4int decayA1 = decayParticle1->getA();
223 const G4int decayZ2 = decayParticle2->getZ();
224 const G4int decayA2 = decayParticle2->getA();
225 const G4int decayZ = decayZ1 + decayZ2;
226 const G4int decayA = decayA1 + decayA2;
227 const G4int daughterZ = c->getZ() - decayZ;
228 const G4int daughterA = c->getA() - decayA;
229 const G4double decayMass1 = decayParticle1->getMass();
230 const G4double decayMass2 = decayParticle2->getMass();
231 const G4double daughterMass = ParticleTable::getTableMass(daughterA,daughterZ);
232
233 // Q-values
234 G4double qValue = motherMass - daughterMass - decayMass1 - decayMass2;
235// assert(qValue > -1e-5); // Q-value should be >0
236 if(qValue<0.)
237 qValue=0.;
238 const G4double qValueB = qValue * Random::shoot();
239
240 // The decay particles behave as if they had more mass until the second
241 // decay
242 const G4double decayMass = decayMass1 + decayMass2 + qValueB;
243
244 /* Stage A: mother --> daughter + (decay1+decay2) */
245
246 // The mother cluster becomes the daughter
247 c->setZ(daughterZ);
248 c->setA(daughterA);
249 c->setMass(daughterMass);
250 c->setExcitationEnergy(0.);
251
252 // Decay kinematics in the mother rest frame
253 const G4double pCMA = KinematicsUtils::momentumInCM(motherMass, daughterMass, decayMass);
254 const ThreeVector momentumA = Random::normVector(pCMA);
255 c->setMomentum(momentumA);
256 c->adjustEnergyFromMomentum();
257 const ThreeVector decayBoostVector = momentumA/std::sqrt(decayMass*decayMass + momentumA.mag2());
258
259 /* Stage B: (decay1+decay2) --> decay1 + decay2 */
260
261 // Decay kinematics in the (decay1+decay2) rest frame
262 const G4double pCMB = KinematicsUtils::momentumInCM(decayMass, decayMass1, decayMass2);
263 const ThreeVector momentumB = Random::normVector(pCMB);
264 decayParticle1->setMomentum(momentumB);
265 decayParticle2->setMomentum(-momentumB);
266 decayParticle1->adjustEnergyFromMomentum();
267 decayParticle2->adjustEnergyFromMomentum();
268
269 // Boost decay1 and decay2 to the Stage-A decay frame
270 decayParticle1->boost(decayBoostVector);
271 decayParticle2->boost(decayBoostVector);
272
273 // Boost all particles to the lab frame
274 decayParticle1->boost(velocity);
275 decayParticle2->boost(velocity);
276 c->boost(velocity);
277
278 // Add the decay particles to the list of decay products
279 decayProducts->push_back(decayParticle1);
280 decayProducts->push_back(decayParticle2);
281 }
282
283 void ClusterDecay::phaseSpaceDecay(Cluster * const c, ParticleTable::ClusterDecayType theDecayMode, ParticleList *decayProducts) {
284 const G4int theA = c->getA();
285 const G4int theZ = c->getZ();
286 const ThreeVector mom(0.0, 0.0, 0.0);
287 const ThreeVector pos = c->getPosition();
288
289 G4int theZStep;
290 ParticleType theEjectileType;
291 switch(theDecayMode) {
293 theZStep = 1;
294 theEjectileType = Proton;
295 break;
297 theZStep = 0;
298 theEjectileType = Neutron;
299 break;
300 default:
301 ERROR("Unrecognized cluster-decay mode in phase-space decay: " << theDecayMode << std::endl
302 << c->print());
303 return;
304 }
305
306 // Find the daughter cluster (first cluster which is not
307 // proton/neutron-unbound, in the sense of the table)
308 G4int finalDaughterZ, finalDaughterA;
310 finalDaughterZ=theZ;
311 finalDaughterA=theA;
312 while(ParticleTable::clusterDecayMode[finalDaughterZ][finalDaughterA]==theDecayMode) {
313 finalDaughterA--;
314 finalDaughterZ -= theZStep;
315 }
316 } else {
317 finalDaughterA = 1;
318 if(theDecayMode==ParticleTable::ProtonUnbound)
319 finalDaughterZ = 1;
320 else
321 finalDaughterZ = 0;
322 }
323// assert(finalDaughterZ<=theZ && finalDaughterA<theA && finalDaughterA>0 && finalDaughterZ>=0);
324 const G4double finalDaughterMass = ParticleTable::getTableMass(finalDaughterA, finalDaughterZ);
325
326 // Compute the available decay energy
327 const G4int nSplits = theA-finalDaughterA;
328 const G4double ejectileMass = ParticleTable::getTableMass(1, theZStep);
329 // c->getMass() can possibly contain some excitation energy, too
330 G4double availableEnergy = c->getMass() - finalDaughterMass - nSplits*ejectileMass;
331// assert(availableEnergy>-1.e-5);
332 if(availableEnergy<0.)
333 availableEnergy=0.;
334
335 // Compute an estimate of the maximum event weight
336 G4double maximumWeight = 1.;
337 G4double eMax = finalDaughterMass + availableEnergy;
338 G4double eMin = finalDaughterMass - ejectileMass;
339 for(G4int iSplit=0; iSplit<nSplits; ++iSplit) {
340 eMax += ejectileMass;
341 eMin += ejectileMass;
342 maximumWeight *= KinematicsUtils::momentumInCM(eMax, eMin, ejectileMass);
343 }
344
345 // Sample decays until the weight cutoff is satisfied
346 G4double weight;
347 std::vector<G4double> theCMMomenta;
348 std::vector<G4double> invariantMasses;
349 G4int nTries=0;
350 /* Maximum number of trials dependent on nSplits. 50 trials seems to be
351 * sufficient for small nSplits. For nSplits>=5, maximumWeight is a gross
352 * overestimate of the actual maximum weight, leading to unreasonably high
353 * rejection rates. For these cases, we set nSplits=1000, although the sane
354 * thing to do would be to improve the importance sampling (maybe by
355 * parametrising maximumWeight?).
356 */
357 G4int maxTries;
358 if(nSplits<5)
359 maxTries=50;
360 else
361 maxTries=1000;
362 do {
363 if(nTries++>maxTries) {
364 WARN("Phase-space decay exceeded the maximum number of rejections (" << maxTries
365 << "). Z=" << theZ << ", A=" << theA << ", E*=" << c->getExcitationEnergy()
366 << ", availableEnergy=" << availableEnergy
367 << ", nSplits=" << nSplits
368 << std::endl);
369 break;
370 }
371
372 // Construct a sorted vector of random numbers
373 std::vector<G4double> randomNumbers;
374 for(G4int iSplit=0; iSplit<nSplits-1; ++iSplit)
375 randomNumbers.push_back(Random::shoot0());
376 std::sort(randomNumbers.begin(), randomNumbers.end());
377
378 // Divide the available decay energy in the right number of steps
379 invariantMasses.clear();
380 invariantMasses.push_back(finalDaughterMass);
381 for(G4int iSplit=0; iSplit<nSplits-1; ++iSplit)
382 invariantMasses.push_back(finalDaughterMass + (iSplit+1)*ejectileMass + randomNumbers.at(iSplit)*availableEnergy);
383 invariantMasses.push_back(c->getMass());
384
385 weight = 1.;
386 theCMMomenta.clear();
387 for(G4int iSplit=0; iSplit<nSplits; ++iSplit) {
388 G4double motherMass = invariantMasses.at(nSplits-iSplit);
389 const G4double daughterMass = invariantMasses.at(nSplits-iSplit-1);
390// assert(motherMass-daughterMass-ejectileMass>-1.e-5);
391 G4double pCM = 0.;
392 if(motherMass-daughterMass-ejectileMass>0.)
393 pCM = KinematicsUtils::momentumInCM(motherMass, daughterMass, ejectileMass);
394 theCMMomenta.push_back(pCM);
395 weight *= pCM;
396 }
397 } while(maximumWeight*Random::shoot()>weight);
398
399 for(G4int iSplit=0; iSplit<nSplits; ++iSplit) {
400 ThreeVector const velocity = -c->boostVector();
401
402#if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
403 const G4double motherMass = c->getMass();
404#endif
405 c->setA(c->getA() - 1);
406 c->setZ(c->getZ() - theZStep);
407 c->setMass(invariantMasses.at(nSplits-iSplit-1));
408
409 Particle *ejectile = new Particle(theEjectileType, mom, pos);
410 ejectile->setTableMass();
411
412// assert(motherMass-c->getMass()-ejectileMass>-1.e-5);
413 ThreeVector momentum;
414 momentum = Random::normVector(theCMMomenta.at(iSplit));
415 c->setMomentum(momentum);
416 c->adjustEnergyFromMomentum();
417 ejectile->setMomentum(-momentum);
418 ejectile->adjustEnergyFromMomentum();
419
420 // Boost to the lab frame
421 c->boost(velocity);
422 ejectile->boost(velocity);
423
424 // Add the decay particle to the list of decay products
425 decayProducts->push_back(ejectile);
426 }
427// assert(std::abs(c->getTableMass()-c->getMass())<1.e-3);
428 c->setExcitationEnergy(0.);
429 }
430}
431
Static class for carrying out cluster decays.
#define WARN(x)
#define DEBUG(x)
#define ERROR(x)
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static ParticleList decay(Cluster *const c)
Carries out a cluster decay.
void setExcitationEnergy(const G4double e)
Set the excitation energy of the cluster.
G4double getExcitationEnergy() const
Get the excitation energy of the cluster.
std::string print() const
static G4double momentumInCM(Particle const *const p1, Particle const *const p2)
gives the momentum in the CM frame of two particles.
static NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
static const G4int clusterTableZSize
static const G4int clusterTableASize
static const ClusterDecayType clusterDecayMode[clusterTableZSize][clusterTableASize]
G4int getZ() const
Returns the charge number.
void setType(ParticleType t)
void setTableMass()
Set the mass of the Particle to its table mass.
G4int getA() const
Returns the baryon number.
static ThreeVector normVector(G4double norm=1.)
Definition: G4INCLRandom.cc:73
static G4double shoot()
Definition: G4INCLRandom.hh:99
static G4double shoot0()
std::list< G4INCL::Particle * > ParticleList