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