Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4CascadeInterface.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 // $Id$
27 //
28 // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29 // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
30 // 20100414 M. Kelsey -- Check for K0L/K0S before using G4InuclElemPart::type
31 // 20100418 M. Kelsey -- Reference output particle lists via const-ref, use
32 // const_iterator for both.
33 // 20100428 M. Kelsey -- Use G4InuclParticleNames enum
34 // 20100429 M. Kelsey -- Change "case gamma:" to "case photon:"
35 // 20100517 M. Kelsey -- Follow new ctors for G4*Collider family.
36 // 20100520 M. Kelsey -- Simplify collision loop, move momentum rotations to
37 // G4CollisionOutput, copy G4DynamicParticle directly from
38 // G4InuclParticle, no switch-block required.
39 // 20100615 M. Kelsey -- Bug fix: For K0's need ekin in GEANT4 units
40 // 20100617 M. Kelsey -- Rename "debug_" preprocessor flag to G4CASCADE_DEBUG,
41 // and "BERTDEV" to "G4CASCADE_COULOMB_DEV"
42 // 20100617 M. Kelsey -- Make G4InuclCollider a local data member
43 // 20100618 M. Kelsey -- Deploy energy-conservation test on final state, with
44 // preprocessor flag G4CASCADE_SKIP_ECONS to skip test.
45 // 20100620 M. Kelsey -- Use new energy-conservation pseudo-collider
46 // 20100621 M. Kelsey -- Fix compiler warning from GCC 4.5
47 // 20100624 M. Kelsey -- Fix cascade loop to check nTries every time (had
48 // allowed for infinite loop on E-violation); dump event data
49 // to output if E-violation exceeds maxTries; use CheckBalance
50 // for baryon and charge conservation.
51 // 20100701 M. Kelsey -- Pass verbosity through to G4CollisionOutput
52 // 20100714 M. Kelsey -- Report number of iterations before success
53 // 20100720 M. Kelsey -- Use G4CASCADE_SKIP_ECONS flag for reporting
54 // 20100723 M. Kelsey -- Move G4CollisionOutput to .hh file for reuse
55 // 20100914 M. Kelsey -- Migrate to integer A and Z
56 // 20100916 M. Kelsey -- Simplify ApplyYourself() by encapsulating code blocks
57 // into numerous functions; make data-member colliders pointers;
58 // provide support for projectile nucleus
59 // 20100919 M. Kelsey -- Fix incorrect logic in retryInelasticNucleus()
60 // 20100922 M. Kelsey -- Add functions to select de-excitation method
61 // 20100924 M. Kelsey -- Migrate to "OutgoingNuclei" names in CollisionOutput
62 // 20110224 M. Kelsey -- Add createTarget() for use with Propagate(); split
63 // conservation law messages to separate function; begin to add
64 // infrastructure code to Propagate. Move verbose
65 // setting to .cc file, and apply to all member objects.
66 // 20110301 M. Kelsey -- Add copyPreviousCascade() for use with Propagate()
67 // along with new buffers and related particle-conversion
68 // functions. Encapulate buffer deletion in clear(). Add some
69 // diagnostic messages.
70 // 20110302 M. Kelsey -- Redo diagnostics inside G4CASCADE_DEBUG_INTERFACE
71 // 20110304 M. Kelsey -- Drop conversion of Propagate() arguments; pass
72 // directly to collider for processing. Rename makeReactionProduct
73 // to makeDynamicParticle.
74 // 20110316 M. Kelsey -- Move kaon-mixing for type-code into G4InuclParticle;
75 // add placeholders to capture and use bullet in Propagte
76 // 20110327 G. Folger -- Set up for E/p checking by G4HadronicProcess in ctor
77 // 20110328 M. Kelsey -- Modify balance() initialization to match Gunter's
78 // 20110404 M. Kelsey -- Get primary projectile from base class (ref-03)
79 // 20110502 M. Kelsey -- Add interface to capture random seeds for user
80 // 20110719 M. Kelsey -- Use trivialise() in case maximum retries are reached
81 // 20110720 M. Kelsey -- Discard elastic-cut array (no longer needed),
82 // discard local "theFinalState" (in base as "theParticleChange"),
83 // Modify createBullet() to set null pointer if bullet unusable,
84 // return empty final-state on failures.
85 // Fix charge violation report before throwing exception.
86 // 20110722 M. Kelsey -- In makeDynamicParticle(), allow invalid type codes
87 // in order to process, e.g., resonances from Propagate() input
88 // 20110728 M. Kelsey -- Per V.Ivantchenko, change NoInteraction to return
89 // zero particles, but set kinetic energy from projectile.
90 // 20110801 M. Kelsey -- Make bullet, target point to local buffers, no delete
91 // 20110802 M. Kelsey -- Use new decay handler for Propagate interface
92 // 20110922 M. Kelsey -- Follow migration of G4InuclParticle::print(), use
93 // G4ExceptionDescription for reporting before throwing exception
94 // 20120125 M. Kelsey -- In retryInelasticProton() check for empty output
95 // 20120525 M. Kelsey -- In NoInteraction, check for Ekin<0., set to zero;
96 // use SetEnergyChange(0.) explicitly for good final states.
97 // 20120822 M. Kelsey -- Move envvars to G4CascadeParameters.
98 
99 #include <cmath>
100 #include <iostream>
101 
102 #include "G4CascadeInterface.hh"
103 #include "globals.hh"
104 #include "G4SystemOfUnits.hh"
105 #include "G4CascadeChannelTables.hh"
106 #include "G4CascadeCheckBalance.hh"
107 #include "G4CascadeParameters.hh"
108 #include "G4CollisionOutput.hh"
109 #include "G4DecayKineticTracks.hh"
110 #include "G4DynamicParticle.hh"
111 #include "G4HadronicException.hh"
112 #include "G4InuclCollider.hh"
114 #include "G4InuclNuclei.hh"
115 #include "G4InuclParticle.hh"
116 #include "G4InuclParticleNames.hh"
117 #include "G4KaonZeroLong.hh"
118 #include "G4KaonZeroShort.hh"
119 #include "G4KineticTrack.hh"
120 #include "G4KineticTrackVector.hh"
121 #include "G4Nucleus.hh"
122 #include "G4ParticleDefinition.hh"
124 #include "G4Track.hh"
125 #include "G4V3DNucleus.hh"
126 
127 using namespace G4InuclParticleNames;
128 
129 typedef std::vector<G4InuclElementaryParticle>::const_iterator particleIterator;
130 typedef std::vector<G4InuclNuclei>::const_iterator nucleiIterator;
131 
132 
133 // Filename to capture random seed, if specified by user at runtime
134 
135 const G4String G4CascadeInterface::randomFile = G4CascadeParameters::randomFile();
136 
137 // Maximum number of iterations allowed for inelastic collision attempts
138 
139 const G4int G4CascadeInterface::maximumTries = 20;
140 
141 
142 // Constructor and destrutor
143 
145  : G4VIntraNuclearTransportModel(name), numberOfTries(0),
146  collider(new G4InuclCollider), balance(new G4CascadeCheckBalance(name)),
147  bullet(0), target(0), output(new G4CollisionOutput) {
149  balance->setLimits(5*perCent, 10*MeV/GeV); // Bertini internal units
150  // Description(); // Model description
151 
154 }
155 
156 
158  clear();
159  delete collider; collider=0;
160  delete balance; balance=0;
161  delete output; output=0;
162 }
163 
165 {
166  outFile << "The Bertini-style cascade implements the inelastic scattering\n"
167  << "of hadrons by nuclei. Nucleons, pions, kaons and hyperons\n"
168  << "from 0 to 15 GeV may be used as projectiles in this model.\n"
169  << "Final state hadrons are produced by a classical cascade\n"
170  << "consisting of individual hadron-nucleon scatterings which use\n"
171  << "free-space partial cross sections, corrected for various\n"
172  << "nuclear medium effects. The target nucleus is modeled as a\n"
173  << "set of 1, 3 or 6 spherical shells, in which scattered hadrons\n"
174  << "travel in straight lines until they are reflected from or\n"
175  << "transmitted through shell boundaries.\n";
176 }
177 
180 }
181 
182 
184  bullet=0;
185  target=0;
186 }
187 
188 
189 // Select post-cascade processing (default will be CascadeDeexcitation)
190 // NOTE: Currently just calls through to Collider, in future will do something
191 
193  collider->useCascadeDeexcitation();
194 }
195 
197  collider->usePreCompoundDeexcitation();
198 }
199 
200 
201 // Apply verbosity to all member objects (overrides base class version)
202 
205  collider->setVerboseLevel(verbose);
206  balance->setVerboseLevel(verbose);
207  output->setVerboseLevel(verbose);
208 }
209 
210 
211 // Test whether inputs are valid for this model
212 
214  G4Nucleus& /* theNucleus */) {
215  return IsApplicable(aTrack.GetDefinition());
216 }
217 
219  if (aPD->GetAtomicMass() > 1) return true; // Nuclei are okay
220 
221  // Valid particle and have interactions available
223  return (type>0 && G4CascadeChannelTables::GetTable(type));
224 }
225 
226 
227 // Main Actions
228 
231  G4Nucleus& theNucleus) {
232  if (verboseLevel)
233  G4cout << " >>> G4CascadeInterface::ApplyYourself" << G4endl;
234 
235  if (aTrack.GetKineticEnergy() < 0.) {
236  G4cerr << " >>> G4CascadeInterface got negative-energy track: "
237  << aTrack.GetDefinition()->GetParticleName() << " Ekin = "
238  << aTrack.GetKineticEnergy() << G4endl;
239  }
240 
241 #ifdef G4CASCADE_DEBUG_INTERFACE
242  static G4int counter(0);
243  counter++;
244  G4cerr << "Reaction number "<< counter << " "
245  << aTrack.GetDefinition()->GetParticleName() << " "
246  << aTrack.GetKineticEnergy() << " MeV" << G4endl;
247 #endif
248 
249  if (!randomFile.empty()) { // User requested random-seed capture
250  if (verboseLevel>1)
251  G4cout << " Saving random engine state to " << randomFile << G4endl;
253  }
254 
256  clear();
257 
258  // Abort processing if no interaction is possible
259  if (!IsApplicable(aTrack, theNucleus)) {
260  if (verboseLevel) G4cerr << " No interaction possible " << G4endl;
261  return NoInteraction(aTrack, theNucleus);
262  }
263 
264  // Make conversion between native Geant4 and Bertini cascade classes.
265  if (!createBullet(aTrack)) {
266  if (verboseLevel) G4cerr << " Unable to create usable bullet" << G4endl;
267  return NoInteraction(aTrack, theNucleus);
268  }
269 
270  if (!createTarget(theNucleus)) {
271  if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
272  return NoInteraction(aTrack, theNucleus);
273  }
274 
275  // Different retry conditions for proton target vs. nucleus
276  const G4bool isHydrogen = (theNucleus.GetA_asInt() == 1);
277 
278  numberOfTries = 0;
279  do { // we try to create inelastic interaction
280  if (verboseLevel > 1)
281  G4cout << " Generating cascade attempt " << numberOfTries << G4endl;
282 
283  output->reset();
284  collider->collide(bullet, target, *output);
285  balance->collide(bullet, target, *output);
286 
287  numberOfTries++;
288  } while ( isHydrogen ? retryInelasticProton() : retryInelasticNucleus() );
289 
290  // Null event if unsuccessful
291  if (numberOfTries >= maximumTries) {
292  if (verboseLevel)
293  G4cout << " Cascade aborted after trials " << numberOfTries << G4endl;
294  return NoInteraction(aTrack, theNucleus);
295  }
296 
297  // Abort job if energy or momentum are not conserved
298  if (!balance->okay()) {
300  return NoInteraction(aTrack, theNucleus);
301  }
302 
303  // Successful cascade -- clean up and return
304  if (verboseLevel) {
305  G4cout << " Cascade output after trials " << numberOfTries << G4endl;
306  if (verboseLevel > 1) output->printCollisionOutput();
307  }
308 
309  // Rotate event to put Z axis along original projectile direction
310  output->rotateEvent(bulletInLabFrame);
311 
313 
314  // Report violations of conservation laws in original frame
316 
317  // Clean up and return final result;
318  clear();
319  return &theParticleChange;
320 }
321 
324  G4V3DNucleus* theNucleus) {
325  if (verboseLevel) G4cout << " >>> G4CascadeInterface::Propagate" << G4endl;
326 
327 #ifdef G4CASCADE_DEBUG_INTERFACE
328  if (verboseLevel>1) {
329  G4cout << " G4V3DNucleus A " << theNucleus->GetMassNumber()
330  << " Z " << theNucleus->GetCharge()
331  << "\n " << theSecondaries->size() << " secondaries:" << G4endl;
332  for (size_t i=0; i<theSecondaries->size(); i++) {
333  G4KineticTrack* kt = (*theSecondaries)[i];
334  G4cout << " " << i << ": " << kt->GetDefinition()->GetParticleName()
335  << " p " << kt->Get4Momentum() << " @ " << kt->GetPosition()
336  << " t " << kt->GetFormationTime() << G4endl;
337  }
338  }
339 #endif
340 
341  if (!randomFile.empty()) { // User requested random-seed capture
342  if (verboseLevel>1)
343  G4cout << " Saving random engine state to " << randomFile << G4endl;
345  }
346 
348  clear();
349 
350  // Process input secondaries list to eliminate resonances
351  G4DecayKineticTracks decay(theSecondaries);
352 
353  // NOTE: Requires 9.4-ref-03 mods to base class and G4TheoFSGenerator
354  const G4HadProjectile* projectile = GetPrimaryProjectile();
355  if (projectile) createBullet(*projectile);
356 
357  if (!createTarget(theNucleus)) {
358  if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
359  return 0; // FIXME: This will cause a segfault later
360  }
361 
362  numberOfTries = 0;
363  do {
364  if (verboseLevel > 1)
365  G4cout << " Generating rescatter attempt " << numberOfTries << G4endl;
366 
367  output->reset();
368  collider->rescatter(bullet, theSecondaries, theNucleus, *output);
369  balance->collide(bullet, target, *output);
370 
371  numberOfTries++;
372  // FIXME: retry checks will SEGFAULT until we can define the bullet!
373  } while (retryInelasticNucleus());
374 
375  // Check whether repeated attempts have all failed; report and exit
376  if (numberOfTries >= maximumTries && !balance->okay()) {
377  throwNonConservationFailure(); // This terminates the job
378  }
379 
380  // Successful cascade -- clean up and return
381  if (verboseLevel) {
382  G4cout << " Cascade rescatter after trials " << numberOfTries << G4endl;
383  if (verboseLevel > 1) output->printCollisionOutput();
384  }
385 
386  // Does calling code take ownership? I hope so!
388 
389  // Clean up and and return final result
390  clear();
391  return propResult;
392 }
393 
394 
395 // Replicate input particles onto output
396 
399  G4Nucleus& /*theNucleus*/) {
400  if (verboseLevel)
401  G4cout << " >>> G4CascadeInterface::NoInteraction" << G4endl;
402 
405 
406  G4double ekin = aTrack.GetKineticEnergy()>0. ? aTrack.GetKineticEnergy() : 0.;
407  theParticleChange.SetEnergyChange(ekin); // Protect against rounding
408 
409  return &theParticleChange;
410 }
411 
412 
413 // Convert input projectile to Bertini internal object
414 
416  const G4ParticleDefinition* trkDef = aTrack.GetDefinition();
417 
418  G4int bulletType = 0; // For elementary particles
419  G4int bulletA = 0, bulletZ = 0; // For nucleus projectile
420 
421  if (trkDef->GetAtomicMass() <= 1) {
422  bulletType = G4InuclElementaryParticle::type(trkDef);
423  } else {
424  bulletA = trkDef->GetAtomicMass();
425  bulletZ = trkDef->GetAtomicNumber();
426  }
427 
428  if (0 == bulletType && 0 == bulletA*bulletZ) {
429  if (verboseLevel) {
430  G4cerr << " G4CascadeInterface: " << trkDef->GetParticleName()
431  << " not usable as bullet." << G4endl;
432  }
433  bullet = 0;
434  return false;
435  }
436 
437  // Code momentum and energy -- Bertini wants z-axis and GeV units
438  G4LorentzVector projectileMomentum = aTrack.Get4Momentum()/GeV;
439 
440  // Rotation/boost to get from z-axis back to original frame
441  bulletInLabFrame = G4LorentzRotation::IDENTITY; // Initialize
442  bulletInLabFrame.rotateZ(-projectileMomentum.phi());
443  bulletInLabFrame.rotateY(-projectileMomentum.theta());
444  bulletInLabFrame.invert();
445 
446  G4LorentzVector momentumBullet(0., 0., projectileMomentum.rho(),
447  projectileMomentum.e());
448 
449  if (bulletType > 0) {
450  hadronBullet.fill(momentumBullet, bulletType);
451  bullet = &hadronBullet;
452  } else {
453  nucleusBullet.fill(momentumBullet, bulletA, bulletZ);
454  bullet = &nucleusBullet;
455  }
456 
457  if (verboseLevel > 2) G4cout << "Bullet: \n" << *bullet << G4endl;
458 
459  return true;
460 }
461 
462 
463 // Convert input nuclear target to Bertini internal object
464 
466  return createTarget(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt());
467 }
468 
470  return createTarget(theNucleus->GetMassNumber(), theNucleus->GetCharge());
471 }
472 
474  if (A > 1) {
475  nucleusTarget.fill(A, Z);
476  target = &nucleusTarget;
477  } else {
478  hadronTarget.fill(0., (Z=1?proton:neutron));
479  target = &hadronTarget;
480  }
481 
482  if (verboseLevel > 2) G4cout << "Target: \n" << *target << G4endl;
483 
484  return true; // Right now, target never fails
485 }
486 
487 
488 // Convert Bertini particle to output (G4DynamicParticle)
489 
492  G4int outgoingType = iep.type();
493 
494  if (iep.quasi_deutron()) {
495  G4cerr << " ERROR: G4CascadeInterface incompatible particle type "
496  << outgoingType << G4endl;
497  return 0;
498  }
499 
500  // Copy local G4DynPart to public output (handle kaon mixing specially)
501  if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
502  G4ThreeVector momDir = iep.getMomentum().vect().unit();
503  G4double ekin = iep.getKineticEnergy()*GeV; // Bertini -> G4 units
504 
506  if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
507 
508  return new G4DynamicParticle(pd, momDir, ekin);
509  } else {
510  return new G4DynamicParticle(iep.getDynamicParticle());
511  }
512 
513  return 0; // Should never get here!
514 }
515 
518  if (verboseLevel > 2) {
519  G4cout << " Nuclei fragment: \n" << inuc << G4endl;
520  }
521 
522  // Copy local G4DynPart to public output
523  return new G4DynamicParticle(inuc.getDynamicParticle());
524 }
525 
526 
527 // Transfer Bertini internal final state to hadronics interface
528 
530  if (verboseLevel > 1)
531  G4cout << " >>> G4CascadeInterface::copyOutputToHadronicResult" << G4endl;
532 
533  const std::vector<G4InuclNuclei>& outgoingNuclei = output->getOutgoingNuclei();
534  const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
535 
538 
539  // Get outcoming particles
540  if (!particles.empty()) {
541  particleIterator ipart = particles.begin();
542  for (; ipart != particles.end(); ipart++) {
544  }
545  }
546 
547  // get nuclei fragments
548  if (!outgoingNuclei.empty()) {
549  nucleiIterator ifrag = outgoingNuclei.begin();
550  for (; ifrag != outgoingNuclei.end(); ifrag++) {
552  }
553  }
554 }
555 
557  if (verboseLevel > 1)
558  G4cout << " >>> G4CascadeInterface::copyOutputToReactionProducts" << G4endl;
559 
560  const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
561  const std::vector<G4InuclNuclei>& fragments = output->getOutgoingNuclei();
562 
564 
565  G4ReactionProduct* rp = 0; // Buffers to create outgoing tracks
566  G4DynamicParticle* dp = 0;
567 
568  // Get outcoming particles
569  if (!particles.empty()) {
570  particleIterator ipart = particles.begin();
571  for (; ipart != particles.end(); ipart++) {
572  rp = new G4ReactionProduct;
573  dp = makeDynamicParticle(*ipart);
574  (*rp) = (*dp); // This does all the necessary copying
575  propResult->push_back(rp);
576  delete dp;
577  }
578  }
579 
580  // get nuclei fragments
581  if (!fragments.empty()) {
582  nucleiIterator ifrag = fragments.begin();
583  for (; ifrag != fragments.end(); ifrag++) {
584  rp = new G4ReactionProduct;
585  dp = makeDynamicParticle(*ifrag);
586  (*rp) = (*dp); // This does all the necessary copying
587  propResult->push_back(rp);
588  delete dp;
589  }
590  }
591 
592  return propResult;
593 }
594 
595 
596 // Report violations of conservation laws in original frame
597 
599  balance->collide(bullet, target, *output);
600 
601  if (verboseLevel > 2) {
602  if (!balance->baryonOkay()) {
603  G4cerr << "ERROR: no baryon number conservation, sum of baryons = "
604  << balance->deltaB() << G4endl;
605  }
606 
607  if (!balance->chargeOkay()) {
608  G4cerr << "ERROR: no charge conservation, sum of charges = "
609  << balance->deltaQ() << G4endl;
610  }
611 
612  if (std::abs(balance->deltaKE()) > 0.01 ) { // GeV
613  G4cerr << "Kinetic energy conservation violated by "
614  << balance->deltaKE() << " GeV" << G4endl;
615  }
616 
617  G4double eInit = bullet->getEnergy() + target->getEnergy();
618  G4double eFinal = eInit + balance->deltaE();
619 
620  G4cout << "Initial energy " << eInit << " final energy " << eFinal
621  << "\nTotal energy conservation at level "
622  << balance->deltaE() * GeV << " MeV" << G4endl;
623 
624  if (balance->deltaKE() > 5.0e-5 ) { // 0.05 MeV
625  G4cerr << "FATAL ERROR: kinetic energy created "
626  << balance->deltaKE() * GeV << " MeV" << G4endl;
627  }
628  }
629 }
630 
631 
632 // Evaluate whether any outgoing particles penetrated Coulomb barrier
633 
635  G4bool violated = false; // by default coulomb analysis is OK
636 
637  const G4double coulumbBarrier = 8.7 * MeV/GeV; // Bertini uses GeV
638 
639  const std::vector<G4InuclElementaryParticle>& p =
640  output->getOutgoingParticles();
641 
642  for (particleIterator ipart=p.begin(); ipart != p.end(); ipart++) {
643  if (ipart->type() == proton) {
644  violated |= (ipart->getKineticEnergy() < coulumbBarrier);
645  }
646  }
647 
648  return violated;
649 }
650 
651 // Check whether inelastic collision on proton failed
652 
654  const std::vector<G4InuclElementaryParticle>& out =
655  output->getOutgoingParticles();
656 
657 #ifdef G4CASCADE_DEBUG_INTERFACE
658  // Report on all retry conditions, in order of return logic
659  G4cout << " retryInelasticProton: number of Tries "
660  << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
661  << "\n retryInelasticProton: AND collision type ";
662  if (out.empty()) G4cout << "FAILED" << G4endl;
663  else {
664  G4cout << (out.size() == 2 ? "ELASTIC (t)" : "INELASTIC (f)")
665  << "\n retryInelasticProton: AND Leading particles bullet "
666  << (out.size() >= 2 &&
667  (out[0].getDefinition() == bullet->getDefinition() ||
668  out[1].getDefinition() == bullet->getDefinition())
669  ? "YES (t)" : "NO (f)")
670  << G4endl;
671  }
672 #endif
673 
674  return ( (numberOfTries < maximumTries) &&
675  (out.empty() ||
676  (out.size() == 2 &&
677  (out[0].getDefinition() == bullet->getDefinition() ||
678  out[1].getDefinition() == bullet->getDefinition())))
679  );
680 }
681 
682 // Check whether generic inelastic collision failed
683 // NOTE: some conditions are set by compiler flags
684 
686  // Quantities necessary for conditional block below
688  G4int nfrag = output->numberOfOutgoingNuclei();
689 
690  const G4ParticleDefinition* firstOut = (npart == 0) ? 0 :
691  output->getOutgoingParticles().begin()->getDefinition();
692 
693 #ifdef G4CASCADE_DEBUG_INTERFACE
694  // Report on all retry conditions, in order of return logic
695  G4cout << " retryInelasticNucleus: numberOfTries "
696  << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
697  << "\n retryInelasticNucleus: AND outputParticles "
698  << ((npart != 0) ? "NON-ZERO (t)" : "EMPTY (f)")
699 #ifdef G4CASCADE_COULOMB_DEV
700  << "\n retryInelasticNucleus: AND coulombBarrier (COULOMB_DEV) "
701  << (coulombBarrierViolation() ? "VIOLATED (t)" : "PASSED (f)")
702  << "\n retryInelasticNucleus: AND collision type (COULOMB_DEV) "
703  << ((npart+nfrag > 2) ? "INELASTIC (t)" : "ELASTIC (f)")
704 #else
705  << "\n retryInelasticNucleus: AND collsion type "
706  << ((npart+nfrag < 3) ? "ELASTIC (t)" : "INELASTIC (f)")
707  << "\n retryInelasticNucleus: AND Leading particle bullet "
708  << ((firstOut == bullet->getDefinition()) ? "YES (t)" : "NO (f)")
709 #endif
710  << "\n retryInelasticNucleus: OR conservation "
711  << (!balance->okay() ? "FAILED (t)" : "PASSED (f)")
712  << G4endl;
713 #endif
714 
715  return ( ((numberOfTries < maximumTries) &&
716  (npart != 0) &&
717 #ifdef G4CASCADE_COULOMB_DEV
718  (coulombBarrierViolation() && npart+nfrag > 2)
719 #else
720  (npart+nfrag < 3 && firstOut == bullet->getDefinition())
721 #endif
722  )
723 #ifndef G4CASCADE_SKIP_ECONS
724  || (!balance->okay())
725 #endif
726  );
727 }
728 
729 
730 // Terminate job in case of persistent non-conservation
731 // FIXME: Need to migrate to G4ExceptionDescription
732 
734  // NOTE: Once G4HadronicException is changed, use the following line!
735  // G4ExceptionDescription errInfo;
736  std::ostream& errInfo = G4cerr;
737 
738  errInfo << " >>> G4CascadeInterface has non-conserving"
739  << " cascade after " << numberOfTries << " attempts." << G4endl;
740 
741  G4String throwMsg = "G4CascadeInterface - ";
742  if (!balance->energyOkay()) {
743  throwMsg += "Energy";
744  errInfo << " Energy conservation violated by " << balance->deltaE()
745  << " GeV (" << balance->relativeE() << ")" << G4endl;
746  }
747 
748  if (!balance->momentumOkay()) {
749  throwMsg += "Momentum";
750  errInfo << " Momentum conservation violated by " << balance->deltaP()
751  << " GeV/c (" << balance->relativeP() << ")" << G4endl;
752  }
753 
754  if (!balance->baryonOkay()) {
755  throwMsg += "Baryon number";
756  errInfo << " Baryon number violated by " << balance->deltaB() << G4endl;
757  }
758 
759  if (!balance->chargeOkay()) {
760  throwMsg += "Charge";
761  errInfo << " Charge conservation violated by " << balance->deltaQ()
762  << G4endl;
763  }
764 
765  errInfo << " Final event output, for debugging:\n"
766  << " Bullet: \n" << *bullet << G4endl
767  << " Target: \n" << *target << G4endl;
768  output->printCollisionOutput(errInfo);
769 
770  throwMsg += " non-conservation. More info in output.";
771  throw G4HadronicException(__FILE__, __LINE__, throwMsg); // Job ends here!
772 }