Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4INCLStandardPropagationModel.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 #define INCLXX_IN_GEANT4_MODE 1
34 
35 #include "globals.hh"
36 
37 /*
38  * StandardPropagationModel.cpp
39  *
40  * \date 4 juin 2009
41  * \author Pekka Kaitaniemi
42  */
43 
45 #include "G4INCLSurfaceAvatar.hh"
47 #include "G4INCLDecayAvatar.hh"
48 #include "G4INCLCrossSections.hh"
49 #include "G4INCLRandom.hh"
50 #include <iostream>
51 #include <algorithm>
52 #include "G4INCLLogger.hh"
53 #include "G4INCLGlobals.hh"
54 #include "G4INCLKinematicsUtils.hh"
58 #include "G4INCLIntersection.hh"
59 
60 namespace G4INCL {
61 
63  :theNucleus(0), maximumTime(70.0), currentTime(0.0), firstAvatar(true),
64  theLocalEnergyType(localEnergyType),
65  theLocalEnergyDeltaType(localEnergyDeltaType)
66  {
67  }
68 
70  {
71  delete theNucleus;
72  }
73 
75  {
76  return theNucleus;
77  }
78 
79  G4double StandardPropagationModel::shoot(ParticleSpecies const projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi) {
80  if(projectileSpecies.theType==Composite)
81  return shootComposite(projectileSpecies, kineticEnergy, impactParameter, phi);
82  else
83  return shootParticle(projectileSpecies.theType, kineticEnergy, impactParameter, phi);
84  }
85 
86  G4double StandardPropagationModel::shootParticle(ParticleType const type, const G4double kineticEnergy, const G4double impactParameter, const G4double phi) {
87  theNucleus->setParticleNucleusCollision();
88  currentTime = 0.0;
89 
90  // Create the projectile particle
91  const G4double projectileMass = ParticleTable::getTableParticleMass(type);
92  G4double energy = kineticEnergy + projectileMass;
93  G4double momentumZ = std::sqrt(energy*energy - projectileMass*projectileMass);
94  ThreeVector momentum(0.0, 0.0, momentumZ);
95  Particle *p= new G4INCL::Particle(type, energy, momentum, ThreeVector());
96 
97  G4double temfin = 0.0;
98  if( p->isNucleon() )
99  temfin = 29.8 * std::pow(theNucleus->getA(), 0.16);
100  else {
101  const G4double tlab = p->getEnergy() - p->getMass();
102  temfin = 30.18 * std::pow(theNucleus->getA(), 0.17*(1.0 - 5.7E-5*tlab));
103  }
104 
105  maximumTime = temfin;
106 
107  // If the incoming particle is slow, use a larger stopping time
108  const G4double rMax = theNucleus->getUniverseRadius();
109  const G4double distance = 2.*rMax;
110  const G4double projectileVelocity = p->boostVector().mag();
111  const G4double traversalTime = distance / projectileVelocity;
112  if(maximumTime < traversalTime)
113  maximumTime = traversalTime;
114  DEBUG("Cascade stopping time is " << maximumTime << std::endl);
115 
116  // If Coulomb is activated, do not process events with impact
117  // parameter larger than the maximum impact parameter, taking into
118  // account Coulomb distortion.
119  if(impactParameter>CoulombDistortion::maxImpactParameter(p->getSpecies(), kineticEnergy, theNucleus)) {
120  DEBUG("impactParameter>CoulombDistortion::maxImpactParameter" << std::endl);
121  delete p;
122  return -1.;
123  }
124 
125  ThreeVector position(impactParameter * std::cos(phi),
126  impactParameter * std::sin(phi),
127  0.);
128  p->setPosition(position);
129 
130  // Fill in the relevant kinematic variables
132  theNucleus->setIncomingMomentum(p->getMomentum());
133  theNucleus->setInitialEnergy(p->getEnergy()
134  + ParticleTable::getTableMass(theNucleus->getA(),theNucleus->getZ()));
135 
136  // Reset the particle kinematics to the INCL values
137  p->setINCLMass();
138  p->setEnergy(p->getMass() + kineticEnergy);
140 
143  firstAvatar = false;
144 
145  // Get the entry avatars from Coulomb and put them in the Store
146  ParticleEntryAvatar *theEntryAvatar = CoulombDistortion::bringToSurface(p, theNucleus);
147  if(theEntryAvatar) {
148  theNucleus->getStore()->addParticleEntryAvatar(theEntryAvatar);
149 
150  theNucleus->setProjectileChargeNumber(p->getZ());
151  theNucleus->setProjectileMassNumber(p->getA());
152 
153  return p->getTransversePosition().mag();
154  } else {
155  delete p;
156  return -1.;
157  }
158  }
159 
160  G4double StandardPropagationModel::shootComposite(ParticleSpecies const species, const G4double kineticEnergy, const G4double impactParameter, const G4double phi) {
161  theNucleus->setNucleusNucleusCollision();
162  currentTime = 0.0;
163 
164  // Create the ProjectileRemnant object
165  ProjectileRemnant *pr = new ProjectileRemnant(species, kineticEnergy);
166 
167  // Same stopping time as for nucleon-nucleus
168  maximumTime = 29.8 * std::pow(theNucleus->getA(), 0.16);
169 
170  // If the incoming cluster is slow, use a larger stopping time
171  const G4double rms = ParticleTable::getNuclearRadius(pr->getA(), pr->getZ());
172  const G4double rMax = theNucleus->getUniverseRadius();
173  const G4double distance = 2.*rMax + 2.725*rms;
174  const G4double projectileVelocity = pr->boostVector().mag();
175  const G4double traversalTime = distance / projectileVelocity;
176  if(maximumTime < traversalTime)
177  maximumTime = traversalTime;
178  DEBUG("Cascade stopping time is " << maximumTime << std::endl);
179 
180  // If Coulomb is activated, do not process events with impact
181  // parameter larger than the maximum impact parameter, taking into
182  // account Coulomb distortion.
183  if(impactParameter>CoulombDistortion::maxImpactParameter(pr,theNucleus)) {
184  pr->deleteParticles();
185  DEBUG("impactParameter>CoulombDistortion::maxImpactParameter" << std::endl);
186  delete pr;
187  return -1.;
188  }
189 
190  // Position the cluster at the right impact parameter
191  ThreeVector position(impactParameter * std::cos(phi),
192  impactParameter * std::sin(phi),
193  0.);
194  pr->setPosition(position);
195 
196  /* Store the internal kinematics of the projectile remnant.
197  *
198  * Note that this is at variance with the Fortran version, which stores
199  * the initial kinematics of the particles *after* putting the spectators
200  * on mass shell, but *before* removing the same energy from the
201  * participants. Due to the different code flow, doing so in the C++
202  * version leads to wrong excitation energies for the forced compound
203  * nucleus.
204  */
205  pr->storeComponents();
206 
207  // Fill in the relevant kinematic variables
209  theNucleus->setIncomingMomentum(pr->getMomentum());
210  theNucleus->setInitialEnergy(pr->getEnergy()
211  + ParticleTable::getTableMass(theNucleus->getA(),theNucleus->getZ()));
212 
214  firstAvatar = false;
215 
216  // Get the entry avatars from Coulomb
217  IAvatarList theAvatarList
218  = CoulombDistortion::bringToSurface(pr, theNucleus);
219 
220  if(theAvatarList.empty()) {
221  DEBUG("No ParticleEntryAvatar found, transparent event" << std::endl);
222  pr->deleteParticles();
223  delete pr;
224  return -1.;
225  }
226 
227  // Tell the Nucleus about the ProjectileRemnant
228  theNucleus->setProjectileRemnant(pr);
229 
230  // Set the number of projectile particles
231  theNucleus->setProjectileChargeNumber(pr->getZ());
232  theNucleus->setProjectileMassNumber(pr->getA());
233 
234  // Register the ParticleEntryAvatars
235  theNucleus->getStore()->addParticleEntryAvatars(theAvatarList);
236 
237  return pr->getTransversePosition().mag();
238  }
239 
241  return maximumTime;
242  }
243 
245 // assert(time>0.0);
246  maximumTime = time;
247  }
248 
250  return currentTime;
251  }
252 
254  {
255  theNucleus = nucleus;
256  }
257 
259  {
260  if(anAvatar) theNucleus->getStore()->add(anAvatar);
261  }
262 
264  // Is either particle a participant?
265  if(!p1->isParticipant() && !p2->isParticipant() && p1->getParticipantType()==p2->getParticipantType()) return NULL;
266 
267  // Is it a pi-resonance collision (we don't treat them)?
268  if((p1->isResonance() && p2->isPion()) || (p1->isPion() && p2->isResonance()))
269  return NULL;
270 
271  // Will the avatar take place between now and the end of the cascade?
272  G4double minDistOfApproachSquared = 0.0;
273  G4double t = getTime(p1, p2, &minDistOfApproachSquared);
274  if(t>maximumTime || t<currentTime) return NULL;
275 
276  // Local energy. Jump through some hoops to calculate the cross section
277  // at the collision point, and clean up after yourself afterwards.
278  G4bool hasLocalEnergy;
279  if(p1->isPion() || p2->isPion())
280  hasLocalEnergy = ((theLocalEnergyDeltaType == FirstCollisionLocalEnergy &&
281  theNucleus->getStore()->getBook()->getAcceptedCollisions()==0) ||
282  theLocalEnergyDeltaType == AlwaysLocalEnergy);
283  else
284  hasLocalEnergy = ((theLocalEnergyType == FirstCollisionLocalEnergy &&
285  theNucleus->getStore()->getBook()->getAcceptedCollisions()==0) ||
286  theLocalEnergyType == AlwaysLocalEnergy);
287  const G4bool p1HasLocalEnergy = (hasLocalEnergy && !p1->isPion());
288  const G4bool p2HasLocalEnergy = (hasLocalEnergy && !p2->isPion());
289 
290  Particle backupParticle1 = *p1;
291  if(p1HasLocalEnergy) {
292  p1->propagate(t - currentTime);
293  if(p1->getPosition().mag() > theNucleus->getSurfaceRadius(p1)) {
294  *p1 = backupParticle1;
295  return NULL;
296  }
298  }
299  Particle backupParticle2 = *p2;
300  if(p2HasLocalEnergy) {
301  p2->propagate(t - currentTime);
302  if(p2->getPosition().mag() > theNucleus->getSurfaceRadius(p2)) {
303  *p2 = backupParticle2;
304  if(p1HasLocalEnergy) {
305  *p1 = backupParticle1;
306  }
307  return NULL;
308  }
310  }
311 
312  // Compute the total cross section
313  const G4double totalCrossSection = CrossSections::total(p1, p2);
314  const G4double squareTotalEnergyInCM = KinematicsUtils::squareTotalEnergyInCM(p1,p2);
315 
316  // Restore particles to their state before the local-energy tweak
317  if(p1HasLocalEnergy) {
318  *p1 = backupParticle1;
319  }
320  if(p2HasLocalEnergy) {
321  *p2 = backupParticle2;
322  }
323 
324  // Is the CM energy > cutNN? (no cutNN on the first collision)
325  if(theNucleus->getStore()->getBook()->getAcceptedCollisions()>0
326  && p1->isNucleon() && p2->isNucleon()
327  && squareTotalEnergyInCM < BinaryCollisionAvatar::cutNNSquared) return NULL;
328 
329  // Do the particles come close enough to each other?
330  if(Math::tenPi*minDistOfApproachSquared > totalCrossSection) return NULL;
331 
332  // Bomb out if the two collision partners are the same particle
333 // assert(p1->getID() != p2->getID());
334 
335  // Return a new avatar, then!
336  return new G4INCL::BinaryCollisionAvatar(t, totalCrossSection, theNucleus, p1, p2);
337  }
338 
340  Intersection theIntersection(
342  aParticle->getPosition(),
343  aParticle->getPropagationVelocity(),
344  theNucleus->getSurfaceRadius(aParticle)));
345  G4double time;
346  if(theIntersection.exists) {
347  time = currentTime + theIntersection.time;
348  } else {
349  ERROR("Imaginary reflection time for particle: " << std::endl
350  << aParticle->print());
351  time = 10000.0;
352  }
353  return time;
354  }
355 
357  G4INCL::Particle const * const particleB, G4double *minDistOfApproach) const
358  {
359  G4double time;
360  G4INCL::ThreeVector t13 = particleA->getPropagationVelocity();
361  t13 -= particleB->getPropagationVelocity();
362  G4INCL::ThreeVector distance = particleA->getPosition();
363  distance -= particleB->getPosition();
364  const G4double t7 = t13.dot(distance);
365  const G4double dt = t13.mag2();
366  if(dt <= 1.0e-10) {
367  (*minDistOfApproach) = 100000.0;
368  return currentTime + 100000.0;
369  } else {
370  time = -t7/dt;
371  }
372  (*minDistOfApproach) = distance.mag2() + time * t7;
373  return currentTime + time;
374  }
375 
376  void StandardPropagationModel::generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles) {
377 
378  // Loop over all the updated particles
379  for(ParticleIter updated = updatedParticles.begin(); updated != updatedParticles.end(); ++updated)
380  {
381  // Loop over all the particles
382  for(ParticleIter particle = particles.begin(); particle != particles.end(); ++particle)
383  {
384  /* Consider the generation of a collision avatar only if (*particle)
385  * is not one of the updated particles.
386  * The criterion makes sure that you don't generate avatars between
387  * updated particles. */
388  if((*particle)->isInList(updatedParticles)) continue;
389 
390  registerAvatar(generateBinaryCollisionAvatar(*particle,*updated));
391  }
392  }
393  }
394 
396 
397  G4bool haveExcept;
398  haveExcept=(except.size()!=0);
399 
400  // Loop over all the particles
401  for(ParticleIter p1 = particles.begin(); p1 != particles.end(); ++p1)
402  {
403  // Loop over the rest of the particles
404  ParticleIter p2 = p1;
405  for(++p2; p2 != particles.end(); ++p2)
406  {
407  // Skip the collision if both particles must be excluded
408  if(haveExcept && (*p1)->isInList(except) && (*p2)->isInList(except)) continue;
409 
411  }
412  }
413 
414  }
415 
417 
418  for(ParticleIter iter = particles.begin(); iter != particles.end(); ++iter) {
419  G4double time = this->getReflectionTime(*iter);
420  if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*iter, time, theNucleus));
421  }
422  ParticleList const &p = theNucleus->getStore()->getParticles();
423  generateUpdatedCollisions(particles, p); // Predict collisions with spectators and participants
424  }
425 
427  ParticleList particles = theNucleus->getStore()->getParticles();
428  if(particles.empty()) { ERROR("No particles inside the nucleus!" << std::endl); }
429  for(ParticleIter i = particles.begin(); i != particles.end(); ++i) {
430  G4double time = this->getReflectionTime(*i);
431  if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*i, time, theNucleus));
432  }
433  ParticleList except;
434  if(excludeUpdated)
435  except = theNucleus->getUpdatedParticles();
436  generateCollisions(particles,except);
437  generateDecays(particles);
438  }
439 
441  for(ParticleIter i = particles.begin(); i != particles.end(); ++i) {
442  if((*i)->isDelta()) {
444  G4double time = currentTime + decayTime;
445  if(time <= maximumTime) {
446  registerAvatar(new DecayAvatar((*i), time, theNucleus));
447  }
448  }
449  }
450  }
451 
453  {
454  // We update only the information related to particles that were updated
455  // by the previous avatar.
456 #ifdef INCL_REGENERATE_AVATARS
457 #warning "The INCL_REGENERATE_AVATARS code has not been tested in a while. Use it at your peril."
458  if(theNucleus->getUpdatedParticles().size()!=0 || theNucleus->getCreatedParticles().size()!=0) {
459  // Regenerates the entire avatar list, skipping collisions between
460  // updated particles
461  theNucleus->getStore()->clearAvatars();
463  generateAllAvatars(true);
464  }
465 #else
466  // Deltas are created by transforming nucleon into a delta for
467  // efficiency reasons
468  Particle * const blockedDelta = theNucleus->getBlockedDelta();
469  ParticleList updatedParticles = theNucleus->getUpdatedParticles();
470  if(blockedDelta)
471  updatedParticles.push_back(blockedDelta);
472  generateDecays(updatedParticles);
473 
474  ParticleList needNewAvatars = theNucleus->getUpdatedParticles();
475  ParticleList created = theNucleus->getCreatedParticles();
476  needNewAvatars.splice(needNewAvatars.end(), created);
477  updateAvatars(needNewAvatars);
478 #endif
479 
480  G4INCL::IAvatar *theAvatar = theNucleus->getStore()->findSmallestTime();
481  if(theAvatar == 0) return 0; // Avatar list is empty
482  // theAvatar->dispose();
483 
484  if(theAvatar->getTime() < currentTime) {
485  ERROR("Avatar time = " << theAvatar->getTime() << ", currentTime = " << currentTime << std::endl);
486  return 0;
487  } else if(theAvatar->getTime() > currentTime) {
488  theNucleus->getStore()->timeStep(theAvatar->getTime() - currentTime);
489 
490  currentTime = theAvatar->getTime();
491  theNucleus->getStore()->getBook()->setCurrentTime(currentTime);
492  }
493 
494  return theAvatar;
495  }
496 
497  void StandardPropagationModel::putSpectatorsOnShell(IAvatarList const &entryAvatars, ParticleList const &spectators) {
498  G4double deltaE = 0.0;
499  for(ParticleIter p=spectators.begin(); p!=spectators.end(); ++p) {
500  // put the spectators on shell (conserving their momentum)
501  const G4double oldEnergy = (*p)->getEnergy();
502  (*p)->setTableMass();
503  (*p)->adjustEnergyFromMomentum();
504  deltaE += (*p)->getEnergy() - oldEnergy;
505  }
506 
507  deltaE /= entryAvatars.size(); // energy to remove from each participant
508 
509  for(IAvatarIter a=entryAvatars.begin(); a!=entryAvatars.end(); ++a) {
510  // remove the energy from the participant
511  Particle *p = (*a)->getParticles().front();
512  ParticleType const t = p->getType();
513  // we also need to slightly correct the participant mass
514  const G4double energy = p->getEnergy() - deltaE
516  p->setEnergy(energy);
517  const G4double newMass = std::sqrt(energy*energy - p->getMomentum().mag2());
518  p->setMass(newMass);
519  }
520  }
521 }