Geant4  10.02.p03
G4INCLNucleus.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 // Alain Boudard, CEA-Saclay, France
28 // Joseph Cugnon, University of Liege, Belgium
29 // Jean-Christophe David, CEA-Saclay, France
30 // Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31 // Sylvie Leray, CEA-Saclay, France
32 // Davide Mancusi, CEA-Saclay, France
33 //
34 #define INCLXX_IN_GEANT4_MODE 1
35 
36 #include "globals.hh"
37 
38 /*
39  * G4INCLNucleus.cc
40  *
41  * \date Jun 5, 2009
42  * \author Pekka Kaitaniemi
43  */
44 
45 #ifndef G4INCLNucleus_hh
46 #define G4INCLNucleus_hh 1
47 
48 #include "G4INCLGlobals.hh"
49 #include "G4INCLLogger.hh"
50 #include "G4INCLParticle.hh"
51 #include "G4INCLIAvatar.hh"
52 #include "G4INCLNucleus.hh"
53 #include "G4INCLKinematicsUtils.hh"
54 #include "G4INCLDecayAvatar.hh"
55 #include "G4INCLCluster.hh"
56 #include "G4INCLClusterDecay.hh"
57 #include "G4INCLDeJongSpin.hh"
58 #include <iterator>
59 #include <cstdlib>
60 #include <sstream>
61 // #include <cassert>
62 
63 namespace G4INCL {
64 
65  Nucleus::Nucleus(G4int mass, G4int charge, Config const * const conf, const G4double universeRadius)
66  : Cluster(charge,mass,true),
67  theInitialZ(charge), theInitialA(mass),
68  theNpInitial(0), theNnInitial(0),
69  initialInternalEnergy(0.),
70  incomingAngularMomentum(0.,0.,0.), incomingMomentum(0.,0.,0.),
71  initialCenterOfMass(0.,0.,0.),
72  remnant(true),
73  initialEnergy(0.),
74  tryCN(false),
75  theUniverseRadius(universeRadius),
76  isNucleusNucleus(false),
77  theProjectileRemnant(NULL),
78  theDensity(NULL),
79  thePotential(NULL)
80  {
81  PotentialType potentialType;
82  G4bool pionPotential;
83  if(conf) {
84  potentialType = conf->getPotentialType();
85  pionPotential = conf->getPionPotential();
86  } else { // By default we don't use energy dependent
87  // potential. This is convenient for some tests.
88  potentialType = IsospinPotential;
89  pionPotential = true;
90  }
91 
92  thePotential = NuclearPotential::createPotential(potentialType, theA, theZ, pionPotential);
93 
96 
98 
100  theParticleSampler->setDensity(theDensity);
101 
102  if(theUniverseRadius<0)
103  theUniverseRadius = theDensity->getMaximumRadius();
104  theStore = new Store(conf);
105  }
106 
108  delete theStore;
110  /* We don't delete the potential and the density here any more -- Factories
111  * are caching them
112  delete thePotential;
113  delete theDensity;*/
114  }
115 
117  // Reset the variables connected with the projectile remnant
118  delete theProjectileRemnant;
119  theProjectileRemnant = NULL;
120 
122  for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
124  }
126  particles.clear();
129  }
130 
132  if(!finalstate) // do nothing if no final state was returned
133  return;
134 
135  G4double totalEnergy = 0.0;
136 
137  FinalStateValidity const validity = finalstate->getValidity();
138  if(validity == ValidFS) {
139 
140  ParticleList const &created = finalstate->getCreatedParticles();
141  for(ParticleIter iter=created.begin(), e=created.end(); iter!=e; ++iter) {
142  theStore->add((*iter));
143  if(!(*iter)->isOutOfWell()) {
144  totalEnergy += (*iter)->getEnergy() - (*iter)->getPotentialEnergy();
145  }
146  }
147 
148  ParticleList const &deleted = finalstate->getDestroyedParticles();
149  for(ParticleIter iter=deleted.begin(), e=deleted.end(); iter!=e; ++iter) {
151  }
152 
153  ParticleList const &modified = finalstate->getModifiedParticles();
154  for(ParticleIter iter=modified.begin(), e=modified.end(); iter!=e; ++iter) {
156  totalEnergy += (*iter)->getEnergy() - (*iter)->getPotentialEnergy();
157  }
158 
159  ParticleList const &out = finalstate->getOutgoingParticles();
160  for(ParticleIter iter=out.begin(), e=out.end(); iter!=e; ++iter) {
161  if((*iter)->isCluster()) {
162  Cluster *clusterOut = dynamic_cast<Cluster*>((*iter));
163 // assert(clusterOut);
164 #ifdef INCLXX_IN_GEANT4_MODE
165  if(!clusterOut)
166  continue;
167 #endif
168  ParticleList const &components = clusterOut->getParticles();
169  for(ParticleIter in=components.begin(), end=components.end(); in!=end; ++in)
171  } else {
173  }
174  totalEnergy += (*iter)->getEnergy(); // No potential here because the particle is gone
175  theA -= (*iter)->getA();
176  theZ -= (*iter)->getZ();
177  theStore->addToOutgoing(*iter);
178  (*iter)->setEmissionTime(theStore->getBook().getCurrentTime());
179  }
180 
181  ParticleList const &entering = finalstate->getEnteringParticles();
182  for(ParticleIter iter=entering.begin(), e=entering.end(); iter!=e; ++iter) {
183  insertParticle(*iter);
184  totalEnergy += (*iter)->getEnergy() - (*iter)->getPotentialEnergy();
185  }
186 
187  // actually perform the removal of the scheduled avatars
189  } else if(validity == ParticleBelowFermiFS || validity == ParticleBelowZeroFS) {
190  INCL_DEBUG("A Particle is entering below the Fermi sea:" << '\n' << finalstate->print() << '\n');
191  tryCN = true;
192  ParticleList const &entering = finalstate->getEnteringParticles();
193  for(ParticleIter iter=entering.begin(), e=entering.end(); iter!=e; ++iter) {
194  insertParticle(*iter);
195  }
196  }
197 
198  if(validity==ValidFS &&
199  std::abs(totalEnergy - finalstate->getTotalEnergyBeforeInteraction()) > 0.1) {
200  INCL_ERROR("Energy nonconservation! Energy at the beginning of the event = "
201  << finalstate->getTotalEnergyBeforeInteraction()
202  <<" and after interaction = "
203  << totalEnergy << '\n'
204  << finalstate->print());
205  }
206  }
207 
209  INCL_WARN("Useless Nucleus::propagateParticles -method called." << '\n');
210  }
211 
213  G4double totalEnergy = 0.0;
214  ParticleList const &inside = theStore->getParticles();
215  for(ParticleIter p=inside.begin(), e=inside.end(); p!=e; ++p) {
216  if((*p)->isNucleon()) // Ugly: we should calculate everything using total energies!
217  totalEnergy += (*p)->getKineticEnergy() - (*p)->getPotentialEnergy();
218  else if((*p)->isResonance())
219  totalEnergy += (*p)->getEnergy() - (*p)->getPotentialEnergy() - ParticleTable::effectiveNucleonMass;
220  else
221  totalEnergy += (*p)->getEnergy() - (*p)->getPotentialEnergy();
222  }
223  return totalEnergy;
224  }
225 
227  // If the remnant consists of only one nucleon, we need to apply a special
228  // procedure to put it on mass shell.
229  if(theA==1) {
230  emitInsidePions();
232  remnant=false;
233  return;
234  }
235 
236  // Compute the recoil momentum and angular momentum
239 
240  ParticleList const &outgoing = theStore->getOutgoingParticles();
241  for(ParticleIter p=outgoing.begin(), e=outgoing.end(); p!=e; ++p) {
242  theMomentum -= (*p)->getMomentum();
243  theSpin -= (*p)->getAngularMomentum();
244  }
248  }
249 
250  // Subtract orbital angular momentum
253 
256  remnant=true;
257  }
258 
260  ThreeVector cm(0.,0.,0.);
261  G4double totalMass = 0.0;
262  ParticleList const &inside = theStore->getParticles();
263  for(ParticleIter p=inside.begin(), e=inside.end(); p!=e; ++p) {
264  const G4double mass = (*p)->getMass();
265  cm += (*p)->getPosition() * mass;
266  totalMass += mass;
267  }
268  cm /= totalMass;
269  return cm;
270  }
271 
273  const G4double totalEnergy = computeTotalEnergy();
274  const G4double separationEnergies = computeSeparationEnergyBalance();
275 
276  return totalEnergy - initialInternalEnergy - separationEnergies;
277  }
278 
279  std::string Nucleus::print()
280  {
281  std::stringstream ss;
282  ss << "Particles in the nucleus:" << '\n'
283  << "Inside:" << '\n';
284  G4int counter = 1;
285  ParticleList const &inside = theStore->getParticles();
286  for(ParticleIter p=inside.begin(), e=inside.end(); p!=e; ++p) {
287  ss << "index = " << counter << '\n'
288  << (*p)->print();
289  counter++;
290  }
291  ss <<"Outgoing:" << '\n';
292  ParticleList const &outgoing = theStore->getOutgoingParticles();
293  for(ParticleIter p=outgoing.begin(), e=outgoing.end(); p!=e; ++p)
294  ss << (*p)->print();
295 
296  return ss.str();
297  }
298 
301  ParticleList deltas;
302  for(ParticleIter i=out.begin(), e=out.end(); i!=e; ++i) {
303  if((*i)->isDelta()) deltas.push_back((*i));
304  }
305  if(deltas.empty()) return false;
306 
307  for(ParticleIter i=deltas.begin(), e=deltas.end(); i!=e; ++i) {
308  INCL_DEBUG("Decay outgoing delta particle:" << '\n'
309  << (*i)->print() << '\n');
310  const ThreeVector beta = -(*i)->boostVector();
311  const G4double deltaMass = (*i)->getMass();
312 
313  // Set the delta momentum to zero and sample the decay in the CM frame.
314  // This makes life simpler if we are using real particle masses.
315  (*i)->setMomentum(ThreeVector());
316  (*i)->setEnergy((*i)->getMass());
317 
318  // Use a DecayAvatar
319  IAvatar *decay = new DecayAvatar((*i), 0.0, NULL);
320  FinalState *fs = decay->getFinalState();
321  Particle * const pion = fs->getCreatedParticles().front();
322  Particle * const nucleon = fs->getModifiedParticles().front();
323 
324  // Adjust the decay momentum if we are using the real masses
325  const G4double decayMomentum = KinematicsUtils::momentumInCM(deltaMass,
326  nucleon->getTableMass(),
327  pion->getTableMass());
328  ThreeVector newMomentum = pion->getMomentum();
329  newMomentum *= decayMomentum / newMomentum.mag();
330 
331  pion->setTableMass();
332  pion->setMomentum(newMomentum);
333  pion->adjustEnergyFromMomentum();
334  pion->setEmissionTime(nucleon->getEmissionTime());
335  pion->boost(beta);
336 
337  nucleon->setTableMass();
338  nucleon->setMomentum(-newMomentum);
339  nucleon->adjustEnergyFromMomentum();
340  nucleon->boost(beta);
341 
342  theStore->addToOutgoing(pion);
343 
344  delete fs;
345  delete decay;
346  }
347 
348  return true;
349  }
350 
352  /* If there is a pion potential, do nothing (deltas will be counted as
353  * excitation energy).
354  * If, however, the remnant is unphysical (Z<0 or Z>A), force the deltas to
355  * decay and get rid of all the pions. In case you're wondering, you can
356  * end up with Z<0 or Z>A if the remnant contains more pi- than protons or
357  * more pi+ than neutrons, respectively.
358  */
359  const G4bool unphysicalRemnant = (theZ<0 || theZ>theA);
360  if(thePotential->hasPionPotential() && !unphysicalRemnant)
361  return false;
362 
363  // Build a list of deltas (avoid modifying the list you are iterating on).
364  ParticleList const &inside = theStore->getParticles();
365  ParticleList deltas;
366  for(ParticleIter i=inside.begin(), e=inside.end(); i!=e; ++i)
367  if((*i)->isDelta()) deltas.push_back((*i));
368 
369  // Loop over the deltas, make them decay
370  for(ParticleIter i=deltas.begin(), e=deltas.end(); i!=e; ++i) {
371  INCL_DEBUG("Decay inside delta particle:" << '\n'
372  << (*i)->print() << '\n');
373  // Create a forced-decay avatar. Note the last boolean parameter. Note
374  // also that if the remnant is unphysical we more or less explicitly give
375  // up energy conservation and CDPP by passing a NULL pointer for the
376  // nucleus.
377  IAvatar *decay;
378  if(unphysicalRemnant) {
379  INCL_WARN("Forcing delta decay inside an unphysical remnant (A=" << theA
380  << ", Z=" << theZ << "). Might lead to energy-violation warnings."
381  << '\n');
382  decay = new DecayAvatar((*i), 0.0, NULL, true);
383  } else
384  decay = new DecayAvatar((*i), 0.0, this, true);
385  FinalState *fs = decay->getFinalState();
386 
387  // The pion can be ejected only if we managed to satisfy energy
388  // conservation and if pion emission does not lead to negative excitation
389  // energies.
390  if(fs->getValidity()==ValidFS) {
391  // Apply the final state to the nucleus
392  applyFinalState(fs);
393  }
394  delete fs;
395  delete decay;
396  }
397 
398  // If the remnant is unphysical, emit all the pions
399  if(unphysicalRemnant) {
400  INCL_DEBUG("Remnant is unphysical: Z=" << theZ << ", A=" << theA << ", emitting all the pions" << '\n');
401  emitInsidePions();
402  }
403 
404  return true;
405  }
406 
409  ParticleList clusters;
410  for(ParticleIter i=out.begin(), e=out.end(); i!=e; ++i) {
411  if((*i)->isCluster()) clusters.push_back((*i));
412  }
413  if(clusters.empty()) return false;
414 
415  for(ParticleIter i=clusters.begin(), e=clusters.end(); i!=e; ++i) {
416  Cluster *cluster = dynamic_cast<Cluster*>(*i); // Can't avoid using a cast here
417 // assert(cluster);
418 #ifdef INCLXX_IN_GEANT4_MODE
419  if(!cluster)
420  continue;
421 #endif
422  cluster->deleteParticles(); // Don't need them
423  ParticleList decayProducts = ClusterDecay::decay(cluster);
424  for(ParticleIter j=decayProducts.begin(), end=decayProducts.end(); j!=end; ++j)
425  theStore->addToOutgoing(*j);
426  }
427  return true;
428  }
429 
431  // Do the phase-space decay only if Z=0 or Z=A
432  if(theA<=1 || (theZ!=0 && theA!=theZ))
433  return false;
434 
435  ParticleList decayProducts = ClusterDecay::decay(this);
436  for(ParticleIter j=decayProducts.begin(), e=decayProducts.end(); j!=e; ++j)
437  theStore->addToOutgoing(*j);
438 
439  return true;
440  }
441 
443  /* Forcing emissions of all pions in the nucleus. This probably violates
444  * energy conservation (although the computation of the recoil kinematics
445  * might sweep this under the carpet).
446  */
447  INCL_WARN("Forcing emissions of all pions in the nucleus." << '\n');
448 
449  // Emit the pions with this kinetic energy
450  const G4double tinyPionEnergy = 0.1; // MeV
451 
452  // Push out the emitted pions
453  ParticleList const &inside = theStore->getParticles();
454  ParticleList toEject;
455  for(ParticleIter i=inside.begin(), e=inside.end(); i!=e; ++i) {
456  if((*i)->isPion()) {
457  Particle * const thePion = *i;
458  INCL_DEBUG("Forcing emission of the following particle: "
459  << thePion->print() << '\n');
461  // Correction for real masses
462  const G4double theQValueCorrection = thePion->getEmissionQValueCorrection(theA,theZ);
463  const G4double kineticEnergyOutside = thePion->getKineticEnergy() - thePion->getPotentialEnergy() + theQValueCorrection;
464  thePion->setTableMass();
465  if(kineticEnergyOutside > 0.0)
466  thePion->setEnergy(thePion->getMass()+kineticEnergyOutside);
467  else
468  thePion->setEnergy(thePion->getMass()+tinyPionEnergy);
469  thePion->adjustMomentumFromEnergy();
470  thePion->setPotentialEnergy(0.);
471  theZ -= thePion->getZ();
472  toEject.push_back(thePion);
473  }
474  }
475  for(ParticleIter i=toEject.begin(), e=toEject.end(); i!=e; ++i) {
477  theStore->addToOutgoing(*i);
478  }
479  }
480 
482 
483  Book const &theBook = theStore->getBook();
484  const G4int nEventCollisions = theBook.getAcceptedCollisions();
485  const G4int nEventDecays = theBook.getAcceptedDecays();
486  const G4int nEventClusters = theBook.getEmittedClusters();
487  if(nEventCollisions==0 && nEventDecays==0 && nEventClusters==0)
488  return true;
489 
490  return false;
491 
492  }
493 
495  // We should be here only if the nucleus contains only one nucleon
496 // assert(theStore->getParticles().size()==1);
497 
498  // No excitation energy!
499  theExcitationEnergy = 0.0;
500 
501  // Move the nucleon to the outgoing list
502  Particle *remN = theStore->getParticles().front();
503  theA -= remN->getA();
504  theZ -= remN->getZ();
506  theStore->addToOutgoing(remN);
508 
509  // Treat the special case of a remaining delta
510  if(remN->isDelta()) {
511  IAvatar *decay = new DecayAvatar(remN, 0.0, NULL);
512  FinalState *fs = decay->getFinalState();
513  // Eject the pion
514  ParticleList const &created = fs->getCreatedParticles();
515  for(ParticleIter j=created.begin(), e=created.end(); j!=e; ++j)
516  theStore->addToOutgoing(*j);
517  delete fs;
518  delete decay;
519  }
520 
521  // Do different things depending on how many outgoing particles we have
522  ParticleList const &outgoing = theStore->getOutgoingParticles();
523  if(outgoing.size() == 2) {
524 
525  INCL_DEBUG("Two particles in the outgoing channel, applying exact two-body kinematics" << '\n');
526 
527  // Can apply exact 2-body kinematics here. Keep the CM emission angle of
528  // the first particle.
529  Particle *p1 = outgoing.front(), *p2 = outgoing.back();
530  const ThreeVector aBoostVector = incomingMomentum / initialEnergy;
531  // Boost to the initial CM
532  p1->boost(aBoostVector);
533  const G4double sqrts = std::sqrt(initialEnergy*initialEnergy - incomingMomentum.mag2());
534  const G4double pcm = KinematicsUtils::momentumInCM(sqrts, p1->getMass(), p2->getMass());
535  const G4double scale = pcm/(p1->getMomentum().mag());
536  // Reset the momenta
537  p1->setMomentum(p1->getMomentum()*scale);
538  p2->setMomentum(-p1->getMomentum());
540  p2->adjustEnergyFromMomentum();
541  // Unboost
542  p1->boost(-aBoostVector);
543  p2->boost(-aBoostVector);
544 
545  } else {
546 
547  INCL_DEBUG("Trying to adjust final-state momenta to achieve energy and momentum conservation" << '\n');
548 
549  const G4int maxIterations=8;
550  G4double totalEnergy, energyScale;
551  G4double val=1.E+100, oldVal=1.E+100, oldOldVal=1.E+100, oldOldOldVal;
552  ThreeVector totalMomentum, deltaP;
553  std::vector<ThreeVector> minMomenta; // use it to store the particle momenta that minimize the merit function
554 
555  // Reserve the vector size
556  minMomenta.reserve(outgoing.size());
557 
558  // Compute the initial total momentum
559  totalMomentum.setX(0.0);
560  totalMomentum.setY(0.0);
561  totalMomentum.setZ(0.0);
562  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i)
563  totalMomentum += (*i)->getMomentum();
564 
565  // Compute the initial total energy
566  totalEnergy = 0.0;
567  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i)
568  totalEnergy += (*i)->getEnergy();
569 
570  // Iterative algorithm starts here:
571  for(G4int iterations=0; iterations < maxIterations; ++iterations) {
572 
573  // Save the old merit-function values
574  oldOldOldVal = oldOldVal;
575  oldOldVal = oldVal;
576  oldVal = val;
577 
578  if(iterations%2 == 0) {
579  INCL_DEBUG("Momentum step" << '\n');
580  // Momentum step: modify all the particle momenta
581  deltaP = incomingMomentum - totalMomentum;
582  G4double pOldTot = 0.0;
583  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i)
584  pOldTot += (*i)->getMomentum().mag();
585  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i) {
586  const ThreeVector mom = (*i)->getMomentum();
587  (*i)->setMomentum(mom + deltaP*mom.mag()/pOldTot);
588  (*i)->adjustEnergyFromMomentum();
589  }
590  } else {
591  INCL_DEBUG("Energy step" << '\n');
592  // Energy step: modify all the particle momenta
593  energyScale = initialEnergy/totalEnergy;
594  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i) {
595  const ThreeVector mom = (*i)->getMomentum();
596  G4double pScale = ((*i)->getEnergy()*energyScale - std::pow((*i)->getMass(),2))/mom.mag2();
597  if(pScale>0) {
598  (*i)->setEnergy((*i)->getEnergy()*energyScale);
599  (*i)->adjustMomentumFromEnergy();
600  }
601  }
602  }
603 
604  // Compute the current total momentum and energy
605  totalMomentum.setX(0.0);
606  totalMomentum.setY(0.0);
607  totalMomentum.setZ(0.0);
608  totalEnergy = 0.0;
609  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i) {
610  totalMomentum += (*i)->getMomentum();
611  totalEnergy += (*i)->getEnergy();
612  }
613 
614  // Merit factor
615  val = std::pow(totalEnergy - initialEnergy,2) +
616  0.25*(totalMomentum - incomingMomentum).mag2();
617  INCL_DEBUG("Merit function: val=" << val << ", oldVal=" << oldVal << ", oldOldVal=" << oldOldVal << ", oldOldOldVal=" << oldOldOldVal << '\n');
618 
619  // Store the minimum
620  if(val < oldVal) {
621  INCL_DEBUG("New minimum found, storing the particle momenta" << '\n');
622  minMomenta.clear();
623  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i)
624  minMomenta.push_back((*i)->getMomentum());
625  }
626 
627  // Stop the algorithm if the search diverges
628  if(val > oldOldVal && oldVal > oldOldOldVal) {
629  INCL_DEBUG("Search is diverging, breaking out of the iteration loop: val=" << val << ", oldVal=" << oldVal << ", oldOldVal=" << oldOldVal << ", oldOldOldVal=" << oldOldOldVal << '\n');
630  break;
631  }
632  }
633 
634  // We should have made at least one successful iteration here
635 // assert(minMomenta.size()==outgoing.size());
636 
637  // Apply the optimal momenta
638  INCL_DEBUG("Applying the solution" << '\n');
639  std::vector<ThreeVector>::const_iterator v = minMomenta.begin();
640  for(ParticleIter i=outgoing.begin(), e=outgoing.end(); i!=e; ++i, ++v) {
641  (*i)->setMomentum(*v);
642  (*i)->adjustEnergyFromMomentum();
643  INCL_DATABLOCK((*i)->print());
644  }
645 
646  }
647 
648  }
649 
651  eventInfo->nParticles = 0;
652  G4bool isNucleonAbsorption = false;
653 
654  G4bool isPionAbsorption = false;
655  // It is possible to have pion absorption event only if the
656  // projectile is pion.
657  if(eventInfo->projectileType == PiPlus ||
658  eventInfo->projectileType == PiMinus ||
659  eventInfo->projectileType == PiZero) {
660  isPionAbsorption = true;
661  }
662 
663  // Forced CN
664  eventInfo->forcedCompoundNucleus = tryCN;
665 
666  // Outgoing particles
667  ParticleList const &outgoingParticles = getStore()->getOutgoingParticles();
668 
669  // Check if we have a nucleon absorption event: nucleon projectile
670  // and no ejected particles.
671  if(outgoingParticles.size() == 0 &&
672  (eventInfo->projectileType == Proton ||
673  eventInfo->projectileType == Neutron)) {
674  isNucleonAbsorption = true;
675  }
676 
677  // Reset the remnant counter
678  eventInfo->nRemnants = 0;
679  eventInfo->history.clear();
680 
681  for(ParticleIter i=outgoingParticles.begin(), e=outgoingParticles.end(); i!=e; ++i ) {
682  // We have a pion absorption event only if the projectile is
683  // pion and there are no ejected pions.
684  if(isPionAbsorption) {
685  if((*i)->isPion()) {
686  isPionAbsorption = false;
687  }
688  }
689 
690  eventInfo->A[eventInfo->nParticles] = (*i)->getA();
691  eventInfo->Z[eventInfo->nParticles] = (*i)->getZ();
692  eventInfo->emissionTime[eventInfo->nParticles] = (*i)->getEmissionTime();
693  eventInfo->EKin[eventInfo->nParticles] = (*i)->getKineticEnergy();
694  ThreeVector mom = (*i)->getMomentum();
695  eventInfo->px[eventInfo->nParticles] = mom.getX();
696  eventInfo->py[eventInfo->nParticles] = mom.getY();
697  eventInfo->pz[eventInfo->nParticles] = mom.getZ();
698  eventInfo->theta[eventInfo->nParticles] = Math::toDegrees(mom.theta());
699  eventInfo->phi[eventInfo->nParticles] = Math::toDegrees(mom.phi());
700  eventInfo->origin[eventInfo->nParticles] = -1;
701  eventInfo->history.push_back("");
702  eventInfo->nParticles++;
703  }
704  eventInfo->nucleonAbsorption = isNucleonAbsorption;
705  eventInfo->pionAbsorption = isPionAbsorption;
706  eventInfo->nCascadeParticles = eventInfo->nParticles;
707 
708  // Projectile-like remnant characteristics
710  eventInfo->ARem[eventInfo->nRemnants] = theProjectileRemnant->getA();
711  eventInfo->ZRem[eventInfo->nRemnants] = theProjectileRemnant->getZ();
713  if(std::abs(eStar)<1E-10)
714  eStar = 0.0; // blame rounding and set the excitation energy to zero
715  eventInfo->EStarRem[eventInfo->nRemnants] = eStar;
716  if(eventInfo->EStarRem[eventInfo->nRemnants]<0.) {
717  INCL_WARN("Negative excitation energy in projectile-like remnant! EStarRem = " << eventInfo->EStarRem[eventInfo->nRemnants] << '\n');
718  }
719  const ThreeVector &spin = theProjectileRemnant->getSpin();
720  if(eventInfo->ARem[eventInfo->nRemnants]%2==0) { // even-A nucleus
721  eventInfo->JRem[eventInfo->nRemnants] = (G4int) (spin.mag()/PhysicalConstants::hc + 0.5);
722  } else { // odd-A nucleus
723  eventInfo->JRem[eventInfo->nRemnants] = ((G4int) (spin.mag()/PhysicalConstants::hc)) + 0.5;
724  }
725  eventInfo->EKinRem[eventInfo->nRemnants] = theProjectileRemnant->getKineticEnergy();
727  eventInfo->pxRem[eventInfo->nRemnants] = mom.getX();
728  eventInfo->pyRem[eventInfo->nRemnants] = mom.getY();
729  eventInfo->pzRem[eventInfo->nRemnants] = mom.getZ();
730  eventInfo->jxRem[eventInfo->nRemnants] = spin.getX() / PhysicalConstants::hc;
731  eventInfo->jyRem[eventInfo->nRemnants] = spin.getY() / PhysicalConstants::hc;
732  eventInfo->jzRem[eventInfo->nRemnants] = spin.getZ() / PhysicalConstants::hc;
733  eventInfo->thetaRem[eventInfo->nRemnants] = Math::toDegrees(mom.theta());
734  eventInfo->phiRem[eventInfo->nRemnants] = Math::toDegrees(mom.phi());
735  eventInfo->nRemnants++;
736  }
737 
738  // Target-like remnant characteristics
739  if(hasRemnant()) {
740  eventInfo->ARem[eventInfo->nRemnants] = getA();
741  eventInfo->ZRem[eventInfo->nRemnants] = getZ();
742  eventInfo->EStarRem[eventInfo->nRemnants] = getExcitationEnergy();
743  if(eventInfo->EStarRem[eventInfo->nRemnants]<0.) {
744  INCL_WARN("Negative excitation energy in target-like remnant! EStarRem = " << eventInfo->EStarRem[eventInfo->nRemnants] << '\n');
745  }
746  const ThreeVector &spin = getSpin();
747  if(eventInfo->ARem[eventInfo->nRemnants]%2==0) { // even-A nucleus
748  eventInfo->JRem[eventInfo->nRemnants] = (G4int) (spin.mag()/PhysicalConstants::hc + 0.5);
749  } else { // odd-A nucleus
750  eventInfo->JRem[eventInfo->nRemnants] = ((G4int) (spin.mag()/PhysicalConstants::hc)) + 0.5;
751  }
752  eventInfo->EKinRem[eventInfo->nRemnants] = getKineticEnergy();
753  const ThreeVector &mom = getMomentum();
754  eventInfo->pxRem[eventInfo->nRemnants] = mom.getX();
755  eventInfo->pyRem[eventInfo->nRemnants] = mom.getY();
756  eventInfo->pzRem[eventInfo->nRemnants] = mom.getZ();
757  eventInfo->jxRem[eventInfo->nRemnants] = spin.getX() / PhysicalConstants::hc;
758  eventInfo->jyRem[eventInfo->nRemnants] = spin.getY() / PhysicalConstants::hc;
759  eventInfo->jzRem[eventInfo->nRemnants] = spin.getZ() / PhysicalConstants::hc;
760  eventInfo->thetaRem[eventInfo->nRemnants] = Math::toDegrees(mom.theta());
761  eventInfo->phiRem[eventInfo->nRemnants] = Math::toDegrees(mom.phi());
762  eventInfo->nRemnants++;
763  }
764 
765  // Global counters, flags, etc.
766  Book const &theBook = theStore->getBook();
767  eventInfo->nCollisions = theBook.getAcceptedCollisions();
768  eventInfo->nBlockedCollisions = theBook.getBlockedCollisions();
769  eventInfo->nDecays = theBook.getAcceptedDecays();
770  eventInfo->nBlockedDecays = theBook.getBlockedDecays();
771  eventInfo->firstCollisionTime = theBook.getFirstCollisionTime();
772  eventInfo->firstCollisionXSec = theBook.getFirstCollisionXSec();
776  eventInfo->nReflectionAvatars = theBook.getAvatars(SurfaceAvatarType);
777  eventInfo->nCollisionAvatars = theBook.getAvatars(CollisionAvatarType);
778  eventInfo->nDecayAvatars = theBook.getAvatars(DecayAvatarType);
780  }
781 
782  Nucleus::ConservationBalance Nucleus::getConservationBalance(const EventInfo &theEventInfo, const G4bool afterRecoil) const {
783  ConservationBalance theBalance;
784  // Initialise balance variables with the incoming values
785  theBalance.Z = theEventInfo.Zp + theEventInfo.Zt;
786  theBalance.A = theEventInfo.Ap + theEventInfo.At;
787 
788  theBalance.energy = getInitialEnergy();
789  theBalance.momentum = getIncomingMomentum();
790 
791  // Process outgoing particles
792  ParticleList const &outgoingParticles = theStore->getOutgoingParticles();
793  for(ParticleIter i=outgoingParticles.begin(), e=outgoingParticles.end(); i!=e; ++i ) {
794  theBalance.Z -= (*i)->getZ();
795  theBalance.A -= (*i)->getA();
796  // For outgoing clusters, the total energy automatically includes the
797  // excitation energy
798  theBalance.energy -= (*i)->getEnergy(); // Note that outgoing particles should have the real mass
799  theBalance.momentum -= (*i)->getMomentum();
800  }
801 
802  // Projectile-like remnant contribution, if present
804  theBalance.Z -= theProjectileRemnant->getZ();
805  theBalance.A -= theProjectileRemnant->getA();
809  theBalance.momentum -= theProjectileRemnant->getMomentum();
810  }
811 
812  // Target-like remnant contribution, if present
813  if(hasRemnant()) {
814  theBalance.Z -= getZ();
815  theBalance.A -= getA();
816  theBalance.energy -= ParticleTable::getTableMass(getA(),getZ()) +
818  if(afterRecoil)
819  theBalance.energy -= getKineticEnergy();
820  theBalance.momentum -= getMomentum();
821  }
822 
823  return theBalance;
824  }
825 
832  }
833 
834  void Nucleus::finalizeProjectileRemnant(const G4double anEmissionTime) {
835  // Deal with the projectile remnant
836  const G4int prA = theProjectileRemnant->getA();
837  if(prA>=1) {
838  // Set the mass
841 
842  // Compute the excitation energy from the invariant mass
843  const G4double anExcitationEnergy = aMass
845 
846  // Set the excitation energy
847  theProjectileRemnant->setExcitationEnergy(anExcitationEnergy);
848 
849  // No spin!
851 
852  // Set the emission time
853  theProjectileRemnant->setEmissionTime(anEmissionTime);
854  }
855  }
856 
857 }
858 
859 #endif
Short_t Zp
Charge number of the projectile nucleus.
static const double cm
Definition: G4SIunits.hh:118
void initializeParticles()
Short_t nRemnants
Number of remnants.
Bool_t nucleonAbsorption
True if the event is a nucleon absorption.
void setMass(G4double mass)
G4double phi() const
ThreeVector incomingMomentum
G4double computeTotalEnergy() const
Compute the current total energy.
Store * getStore() const
Simple class implementing De Jong&#39;s spin model for nucleus-nucleus collisions.
NuclearDensity const * createDensity(const G4int A, const G4int Z)
ThreeVector incomingAngularMomentum
Int_t nEnergyViolationInteraction
Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a so...
Float_t firstCollisionTime
Time of the first collision [fm/c].
std::string print()
G4bool pion(G4int ityp)
G4double mag2() const
G4double getExcitationEnergy() const
Get the excitation energy of the cluster.
Float_t emissionTime[maxSizeParticles]
Emission time [fm/c].
virtual ~Nucleus()
G4double getEmissionTime()
Int_t nCollisions
Number of accepted two-body collisions.
G4double getMass() const
Get the cached particle mass.
Float_t py[maxSizeParticles]
Particle momentum, y component [MeV/c].
G4bool decayOutgoingClusters()
Force the decay of unstable outgoing clusters.
Float_t firstCollisionXSec
Cross section of the first collision (mb)
ifstream in
Definition: comparison.C:7
G4bool isEventTransparent() const
Is the event transparent?
void setEmissionTime(G4double t)
Static class for carrying out cluster decays.
G4int getAcceptedDecays() const
Definition: G4INCLBook.hh:102
void computeOneNucleonRecoilKinematics()
Compute the recoil kinematics for a 1-nucleon remnant.
ParticleList const & getOutgoingParticles() const
Definition: G4INCLStore.hh:223
void applyFinalState(FinalState *)
PotentialType getPotentialType() const
Get the type of the potential for nucleons.
#define INCL_ERROR(x)
const G4double hc
[MeV*fm]
G4double getCurrentTime() const
Definition: G4INCLBook.hh:98
void boost(const ThreeVector &aBoostVector)
G4double adjustEnergyFromMomentum()
Recompute the energy to match the momentum.
Bool_t forcedCompoundNucleus
True if the event is a forced CN.
std::string print() const
G4double initialInternalEnergy
Float_t jzRem[maxSizeRemnants]
Remnant angular momentum, z component [ ].
void add(Particle *p)
Definition: G4INCLStore.cc:58
ThreeVector const & getSpin() const
Get the spin of the nucleus.
G4double computeExcitationEnergy() const
Compute the current excitation energy.
Int_t nBlockedCollisions
Number of two-body collisions blocked by Pauli or CDPP.
Struct for conservation laws.
#define INCL_WARN(x)
G4double getSeparationEnergy(const Particle *const p) const
Return the separation energy for a particle.
Bool_t pionAbsorption
True if the event is a pion absorption.
G4double toDegrees(G4double radians)
G4double mag() const
Float_t EKin[maxSizeParticles]
Particle kinetic energy [MeV].
G4double momentumInCM(Particle const *const p1, Particle const *const p2)
gives the momentum in the CM frame of two particles.
G4int getA() const
Returns the baryon number.
G4int getBlockedCollisions() const
Definition: G4INCLBook.hh:101
Float_t JRem[maxSizeRemnants]
Remnant spin [ ].
G4int getZ() const
Returns the charge number.
G4double computeSeparationEnergyBalance() const
Outgoing - incoming separation energies.
Int_t nReflectionAvatars
Number of reflection avatars.
G4double getKineticEnergy() const
Get the particle kinetic energy.
Float_t EKinRem[maxSizeRemnants]
Remnant kinetic energy [MeV].
int G4int
Definition: G4Types.hh:78
void setProtonSeparationEnergy(const G4double s)
Setter for protonSeparationEnergy.
Float_t pyRem[maxSizeRemnants]
Remnant momentum, y component [MeV/c].
Int_t nCollisionAvatars
Number of collision avatars.
Float_t pz[maxSizeParticles]
Particle momentum, z component [MeV/c].
std::string print() const
void updatePotentialEnergy(Particle *p) const
Update the particle potential energy.
Short_t At
Mass number of the target nucleus.
void setY(G4double ay)
Set the y coordinate.
G4double getInitialEnergy() const
Get the initial energy.
void setEnergy(G4double energy)
Short_t ZRem[maxSizeRemnants]
Remnant charge number.
Short_t ARem[maxSizeRemnants]
Remnant mass number.
Short_t Zt
Charge number of the target nucleus.
void removeScheduledAvatars()
Remove avatars that have been scheduled.
Definition: G4INCLStore.cc:134
Float_t theta[maxSizeParticles]
Particle momentum polar angle [radians].
void useFusionKinematics()
Adjust the kinematics for complete-fusion events.
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
G4bool decayOutgoingDeltas()
Force the decay of outgoing deltas.
Short_t Z[maxSizeParticles]
Particle charge number.
ThreeVector theSpin
ParticleList const & getOutgoingParticles() const
void setDensity(NuclearDensity const *const d)
Setter for theDensity.
G4int getEnergyViolationInteraction() const
Definition: G4INCLBook.hh:107
Book & getBook()
Definition: G4INCLStore.hh:259
Short_t nParticles
Number of particles in the final state.
G4bool nucleon(G4int ityp)
void particleHasBeenDestroyed(Particle *const)
Definition: G4INCLStore.cc:181
G4double getFirstCollisionTime() const
Definition: G4INCLBook.hh:83
Float_t jxRem[maxSizeRemnants]
Remnant angular momentum, x component [ ].
bool G4bool
Definition: G4Types.hh:79
virtual void initializeParticles()
Initialise the NuclearDensity pointer and sample the particles.
Double_t scale
ParticleList const & getDestroyedParticles() const
void particleHasBeenEjected(Particle *const)
Definition: G4INCLStore.cc:175
G4double getExcitationEnergy() const
Get the excitation energy of the nucleus.
Float_t pzRem[maxSizeRemnants]
Remnant momentum, z component [MeV/c].
ParticleList particles
G4int getAcceptedCollisions() const
Definition: G4INCLBook.hh:100
G4bool isDelta() const
Is it a Delta?
NuclearDensity const * theDensity
Pointer to the NuclearDensity object.
void setPotentialEnergy(G4double v)
Set the particle potential energy.
Float_t EStarRem[maxSizeRemnants]
Remnant excitation energy [MeV].
virtual G4double getTableMass() const
Get the tabulated particle mass.
ParticleList const & getParticles() const
Definition: G4INCLStore.hh:253
G4bool hasRemnant() const
Does the nucleus give a cascade remnant?
G4int getAvatars(AvatarType type) const
Definition: G4INCLBook.hh:104
Int_t nDecayAvatars
Number of decay avatars.
void deleteProjectileRemnant()
Delete the projectile remnant.
G4INCL::ThreeVector thePosition
G4bool decayMe()
Force the phase-space decay of the Nucleus.
Float_t jyRem[maxSizeRemnants]
Remnant angular momentum, y component [ ].
ParticleSampler * theParticleSampler
Float_t firstCollisionSpectatorPosition
Position of the spectator on the first collision (fm)
Float_t phiRem[maxSizeRemnants]
Remnant momentum azimuthal angle [radians].
Int_t nBlockedDecays
Number of decays blocked by Pauli or CDPP.
ParticleList const & getEnteringParticles() const
ParticleList const & getParticles() const
G4bool getPionPotential() const
Do we want the pion potential?
G4int getBlockedDecays() const
Definition: G4INCLBook.hh:103
Float_t pxRem[maxSizeRemnants]
Remnant momentum, x component [MeV/c].
G4double theta() const
Short_t A[maxSizeParticles]
Particle mass number.
Int_t nDecays
Number of accepted Delta decays.
Bool_t firstCollisionIsElastic
True if the first collision was elastic.
G4double theExcitationEnergy
void particleHasBeenUpdated(Particle *const)
Notify the Store about a particle update.
Definition: G4INCLStore.cc:127
ThreeVector computeCenterOfMass() const
Compute the current center-of-mass position.
void propagateParticles(G4double step)
G4bool decayInsideDeltas()
Force the decay of deltas inside the nucleus.
void finalizeProjectileRemnant(const G4double emissionTime)
Finalise the projectile remnant.
G4double getPotentialEnergy() const
Get the particle potential energy.
NuclearPotential::INuclearPotential const * thePotential
Pointer to the NuclearPotential object.
G4double getFirstCollisionSpectatorPosition() const
Definition: G4INCLBook.hh:89
void emitInsidePions()
Force emission of all pions inside the nucleus.
G4double getFirstCollisionXSec() const
Definition: G4INCLBook.hh:86
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
void setTableMass()
Set the mass of the Particle to its table mass.
G4INCL::ThreeVector getAngularMomentum() const
Get the total angular momentum (orbital + spin)
FinalStateValidity getValidity() const
G4bool hasPionPotential() const
Do we have a pion potential?
void setPotential(NuclearPotential::INuclearPotential const *const p)
Setter for thePotential.
Float_t thetaRem[maxSizeRemnants]
Remnant momentum polar angle [radians].
G4double getZ() const
FinalState * getFinalState()
void fillEventInfo(EventInfo *eventInfo)
ParticleList const & getModifiedParticles() const
G4INCL::ThreeVector theMomentum
void setX(G4double ax)
Set the x coordinate.
G4double getFirstCollisionSpectatorMomentum() const
Definition: G4INCLBook.hh:92
void insertParticle(Particle *p)
Insert a new particle (e.g. a projectile) in the nucleus.
Float_t phi[maxSizeParticles]
Particle momentum azimuthal angle [radians].
G4int getEmittedClusters() const
Definition: G4INCLBook.hh:106
double G4double
Definition: G4Types.hh:76
Nucleus(G4int mass, G4int charge, Config const *const conf, const G4double universeRadius=-1.)
virtual G4double getTableMass() const
Get the real particle mass.
void setZ(G4double az)
Set the z coordinate.
#define INCL_DEBUG(x)
const G4INCL::ThreeVector & getMomentum() const
void computeRecoilKinematics()
Compute the recoil momentum and spin of the nucleus.
Float_t firstCollisionSpectatorMomentum
Momentum of the spectator on the first collision (fm)
G4double getInvariantMass() const
Get the the particle invariant mass.
ThreeVector initialCenterOfMass
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
const ThreeVector & getIncomingMomentum() const
Get the incoming momentum vector.
void addToOutgoing(Particle *p)
add the particle to the outgoing particle list.
Definition: G4INCLStore.hh:190
ConservationBalance getConservationBalance(EventInfo const &theEventInfo, const G4bool afterRecoil) const
Compute charge, mass, energy and momentum balance.
Float_t px[maxSizeParticles]
Particle momentum, x component [MeV/c].
G4double initialEnergy
INuclearPotential const * createPotential(const PotentialType type, const G4int theA, const G4int theZ, const G4bool pionPotential)
Create an INuclearPotential object.
Short_t Ap
Mass number of the projectile nucleus.
Short_t nCascadeParticles
Number of cascade particles.
G4double getY() const
void setExcitationEnergy(const G4double e)
Set the excitation energy of the cluster.
G4bool getFirstCollisionIsElastic() const
Definition: G4INCLBook.hh:95
const G4double effectiveNucleonMass
ProjectileRemnant * theProjectileRemnant
Pointer to the quasi-projectile.
Short_t origin[maxSizeParticles]
Origin of the particle.
#define INCL_DATABLOCK(x)
void setNeutronSeparationEnergy(const G4double s)
Setter for protonSeparationEnergy.
void setSpin(const ThreeVector &j)
Set the spin of the nucleus.
ParticleList::const_iterator ParticleIter
ParticleList const & getCreatedParticles() const
G4double getEmissionQValueCorrection(const G4int AParent, const G4int ZParent) const
Computes correction on the emission Q-value.
G4double getX() const
Int_t projectileType
Projectile particle type.
virtual void setMomentum(const G4INCL::ThreeVector &momentum)
G4double getTotalEnergyBeforeInteraction() const
G4double theUniverseRadius
The radius of the universe.
std::vector< std::string > history
History of the particle.