Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4CascadeInterface.cc
이 파일의 문서화 페이지로 가기
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: G4CascadeInterface.cc 71719 2013-06-21 00:01:54Z mkelsey $
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 // 20130508 D. Wright -- Add support for muon capture
99 // 20130804 M. Kelsey -- Fix bug #1513 -- "(Z=1)" in boolean expression
100 // 20140116 M. Kelsey -- Move statics to const data members to avoid weird
101 // interactions with MT.
102 // 20140929 M. Kelsey -- Explicitly call useCascadeDeexcitation() in ctor
103 // 20150506 M. Kelsey -- Call Initialize() in ctor for master thread only
104 // 20150608 M. Kelsey -- Label all while loops as terminating.
105 
106 #include <cmath>
107 #include <iostream>
108 
109 #include "G4CascadeInterface.hh"
110 #include "globals.hh"
111 #include "G4SystemOfUnits.hh"
112 #include "G4CascadeChannelTables.hh"
113 #include "G4CascadeCheckBalance.hh"
114 #include "G4CascadeParameters.hh"
115 #include "G4CollisionOutput.hh"
116 #include "G4DecayKineticTracks.hh"
117 #include "G4DynamicParticle.hh"
118 #include "G4HadronicException.hh"
119 #include "G4InuclCollider.hh"
121 #include "G4InuclNuclei.hh"
122 #include "G4InuclParticle.hh"
123 #include "G4InuclParticleNames.hh"
124 #include "G4KaonZeroLong.hh"
125 #include "G4KaonZeroShort.hh"
126 #include "G4KineticTrack.hh"
127 #include "G4KineticTrackVector.hh"
128 #include "G4Nucleus.hh"
129 #include "G4ParticleDefinition.hh"
131 #include "G4Threading.hh"
132 #include "G4Track.hh"
133 #include "G4V3DNucleus.hh"
134 #include "G4UnboundPN.hh"
135 #include "G4Dineutron.hh"
136 #include "G4Diproton.hh"
137 
138 using namespace G4InuclParticleNames;
139 
140 typedef std::vector<G4InuclElementaryParticle>::const_iterator particleIterator;
141 typedef std::vector<G4InuclNuclei>::const_iterator nucleiIterator;
142 
143 
144 // Constructor and destrutor
145 
148  randomFile(G4CascadeParameters::randomFile()),
149  maximumTries(20), numberOfTries(0),
150  collider(new G4InuclCollider), balance(new G4CascadeCheckBalance(name)),
151  bullet(0), target(0), output(new G4CollisionOutput) {
152  // Set up global objects for master thread or sequential build
154 
156  balance->setLimits(5*perCent, 10*MeV/GeV); // Bertini internal units
158 
161  else
163 }
164 
166  clear();
167  delete collider; collider=0;
168  delete balance; balance=0;
169  delete output; output=0;
170 }
171 
172 void G4CascadeInterface::ModelDescription(std::ostream& outFile) const
173 {
174  outFile << "The Bertini-style cascade implements the inelastic scattering\n"
175  << "of hadrons by nuclei. Nucleons, pions, kaons and hyperons\n"
176  << "from 0 to 15 GeV may be used as projectiles in this model.\n"
177  << "Final state hadrons are produced by a classical cascade\n"
178  << "consisting of individual hadron-nucleon scatterings which use\n"
179  << "free-space partial cross sections, corrected for various\n"
180  << "nuclear medium effects. The target nucleus is modeled as a\n"
181  << "set of 1, 3 or 6 spherical shells, in which scattered hadrons\n"
182  << "travel in straight lines until they are reflected from or\n"
183  << "transmitted through shell boundaries.\n";
184 }
185 
186 void G4CascadeInterface::DumpConfiguration(std::ostream& outFile) const {
188 }
189 
191  bullet=0;
192  target=0;
193 }
194 
195 
196 // Initialize shared objects for use across multiple threads
197 
203  if (!ch || !pn || !nn || !pp) return; // Avoid "unused variables"
204 }
205 
206 
207 // Select post-cascade processing (default will be CascadeDeexcitation)
208 // NOTE: Currently just calls through to Collider, in future will do something
209 
212 }
213 
216 }
217 
218 
219 // Apply verbosity to all member objects (overrides base class version)
220 
223  collider->setVerboseLevel(verbose);
224  balance->setVerboseLevel(verbose);
225  output->setVerboseLevel(verbose);
226 }
227 
228 
229 // Test whether inputs are valid for this model
230 
232  G4Nucleus& /* theNucleus */) {
233  return IsApplicable(aTrack.GetDefinition());
234 }
235 
237  if (aPD->GetAtomicMass() > 1) return true; // Nuclei are okay
238 
239  // Valid particle and have interactions available
241  return (G4CascadeChannelTables::GetTable(type));
242 }
243 
244 
245 // Main Actions
246 
249  G4Nucleus& theNucleus) {
250  if (verboseLevel)
251  G4cout << " >>> G4CascadeInterface::ApplyYourself" << G4endl;
252 
253  if (aTrack.GetKineticEnergy() < 0.) {
254  G4cerr << " >>> G4CascadeInterface got negative-energy track: "
255  << aTrack.GetDefinition()->GetParticleName() << " Ekin = "
256  << aTrack.GetKineticEnergy() << G4endl;
257  }
258 
259 #ifdef G4CASCADE_DEBUG_INTERFACE
260  static G4int counter(0);
261  counter++;
262  G4cerr << "Reaction number "<< counter << " "
263  << aTrack.GetDefinition()->GetParticleName() << " "
264  << aTrack.GetKineticEnergy() << " MeV" << G4endl;
265 #endif
266 
267  if (!randomFile.empty()) { // User requested random-seed capture
268  if (verboseLevel>1)
269  G4cout << " Saving random engine state to " << randomFile << G4endl;
271  }
272 
274  clear();
275 
276  // Abort processing if no interaction is possible
277  if (!IsApplicable(aTrack, theNucleus)) {
278  if (verboseLevel) G4cerr << " No interaction possible " << G4endl;
279  return NoInteraction(aTrack, theNucleus);
280  }
281 
282  // Make conversion between native Geant4 and Bertini cascade classes.
283  if (!createBullet(aTrack)) {
284  if (verboseLevel) G4cerr << " Unable to create usable bullet" << G4endl;
285  return NoInteraction(aTrack, theNucleus);
286  }
287 
288  if (!createTarget(theNucleus)) {
289  if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
290  return NoInteraction(aTrack, theNucleus);
291  }
292 
293  // Different retry conditions for proton target vs. nucleus
294  const G4bool isHydrogen = (theNucleus.GetA_asInt() == 1);
295 
296  numberOfTries = 0;
297  do { // we try to create inelastic interaction
298  if (verboseLevel > 1)
299  G4cout << " Generating cascade attempt " << numberOfTries << G4endl;
300 
301  output->reset();
304 
305  numberOfTries++;
306  /* Loop checking 08.06.2015 MHK */
307  } while ( isHydrogen ? retryInelasticProton() : retryInelasticNucleus() );
308 
309  // Null event if unsuccessful
310  if (numberOfTries >= maximumTries) {
311  if (verboseLevel)
312  G4cout << " Cascade aborted after trials " << numberOfTries << G4endl;
313  return NoInteraction(aTrack, theNucleus);
314  }
315 
316  // Abort job if energy or momentum are not conserved
317  if (!balance->okay()) {
319  return NoInteraction(aTrack, theNucleus);
320  }
321 
322  // Successful cascade -- clean up and return
323  if (verboseLevel) {
324  G4cout << " Cascade output after trials " << numberOfTries << G4endl;
326  }
327 
328  // Rotate event to put Z axis along original projectile direction
329  // Removed by DHW to fix bug #1990
330  // output->rotateEvent(bulletInLabFrame);
331 
333 
334  // Report violations of conservation laws in original frame
336 
337  // Clean up and return final result;
338  clear();
339 /*
340  G4int nSec = theParticleChange.GetNumberOfSecondaries();
341  for (G4int i = 0; i < nSec; i++) {
342  G4HadSecondary* sec = theParticleChange.GetSecondary(i);
343  G4DynamicParticle* dp = sec->GetParticle();
344  if (dp->GetDefinition()->GetParticleName() == "neutron")
345  G4cout << dp->GetDefinition()->GetParticleName() << " has "
346  << dp->GetKineticEnergy()/MeV << " MeV " << G4endl;
347  }
348 */
349  return &theParticleChange;
350 }
351 
354  G4V3DNucleus* theNucleus) {
355  if (verboseLevel) G4cout << " >>> G4CascadeInterface::Propagate" << G4endl;
356 
357 #ifdef G4CASCADE_DEBUG_INTERFACE
358  if (verboseLevel>1) {
359  G4cout << " G4V3DNucleus A " << theNucleus->GetMassNumber()
360  << " Z " << theNucleus->GetCharge()
361  << "\n " << theSecondaries->size() << " secondaries:" << G4endl;
362  for (size_t i=0; i<theSecondaries->size(); i++) {
363  G4KineticTrack* kt = (*theSecondaries)[i];
364  G4cout << " " << i << ": " << kt->GetDefinition()->GetParticleName()
365  << " p " << kt->Get4Momentum() << " @ " << kt->GetPosition()
366  << " t " << kt->GetFormationTime() << G4endl;
367  }
368  }
369 #endif
370 
371  if (!randomFile.empty()) { // User requested random-seed capture
372  if (verboseLevel>1)
373  G4cout << " Saving random engine state to " << randomFile << G4endl;
375  }
376 
378  clear();
379 
380  // Process input secondaries list to eliminate resonances
381  G4DecayKineticTracks decay(theSecondaries);
382 
383  // NOTE: Requires 9.4-ref-03 mods to base class and G4TheoFSGenerator
384  const G4HadProjectile* projectile = GetPrimaryProjectile();
385  if (projectile) createBullet(*projectile);
386 
387  if (!createTarget(theNucleus)) {
388  if (verboseLevel) G4cerr << " Unable to create usable target" << G4endl;
389  return 0; // FIXME: This will cause a segfault later
390  }
391 
392  numberOfTries = 0;
393  do {
394  if (verboseLevel > 1)
395  G4cout << " Generating rescatter attempt " << numberOfTries << G4endl;
396 
397  output->reset();
398  collider->rescatter(bullet, theSecondaries, theNucleus, *output);
400 
401  numberOfTries++;
402  // FIXME: retry checks will SEGFAULT until we can define the bullet!
403  } while (retryInelasticNucleus()); /* Loop checking 08.06.2015 MHK */
404 
405  // Check whether repeated attempts have all failed; report and exit
406  if (numberOfTries >= maximumTries && !balance->okay()) {
407  throwNonConservationFailure(); // This terminates the job
408  }
409 
410  // Successful cascade -- clean up and return
411  if (verboseLevel) {
412  G4cout << " Cascade rescatter after trials " << numberOfTries << G4endl;
414  }
415 
416  // Does calling code take ownership? I hope so!
418 
419  // Clean up and and return final result
420  clear();
421  return propResult;
422 }
423 
424 
425 // Replicate input particles onto output
426 
429  G4Nucleus& /*theNucleus*/) {
430  if (verboseLevel)
431  G4cout << " >>> G4CascadeInterface::NoInteraction" << G4endl;
432 
435 
436  G4double ekin = aTrack.GetKineticEnergy()>0. ? aTrack.GetKineticEnergy() : 0.;
437  theParticleChange.SetEnergyChange(ekin); // Protect against rounding
438 
439  return &theParticleChange;
440 }
441 
442 
443 // Convert input projectile to Bertini internal object
444 
446  const G4ParticleDefinition* trkDef = aTrack.GetDefinition();
447  G4int bulletType = 0; // For elementary particles
448  G4int bulletA = 0, bulletZ = 0; // For nucleus projectile
449 
450  if (trkDef->GetAtomicMass() <= 1) {
451  bulletType = G4InuclElementaryParticle::type(trkDef);
452  } else {
453  bulletA = trkDef->GetAtomicMass();
454  bulletZ = trkDef->GetAtomicNumber();
455  }
456 
457  if (0 == bulletType && 0 == bulletA*bulletZ) {
458  if (verboseLevel) {
459  G4cerr << " G4CascadeInterface: " << trkDef->GetParticleName()
460  << " not usable as bullet." << G4endl;
461  }
462  bullet = 0;
463  return false;
464  }
465 
466  // Code momentum and energy -- Bertini wants z-axis and GeV units
467  G4LorentzVector projectileMomentum = aTrack.Get4Momentum()/GeV;
468 
469  // Rotation/boost to get from z-axis back to original frame
470  // According to bug report #1990 this rotation is unnecessary and causes
471  // irreproducibility. Verifed and fixed by DHW 27 Nov 2017
472  // bulletInLabFrame = G4LorentzRotation::IDENTITY; // Initialize
473  // bulletInLabFrame.rotateZ(-projectileMomentum.phi());
474  // bulletInLabFrame.rotateY(-projectileMomentum.theta());
475  // bulletInLabFrame.invert();
476 
477  G4LorentzVector momentumBullet(0., 0., projectileMomentum.rho(),
478  projectileMomentum.e());
479 
480  if (G4InuclElementaryParticle::valid(bulletType)) {
481  hadronBullet.fill(momentumBullet, bulletType);
482  bullet = &hadronBullet;
483  } else {
484  nucleusBullet.fill(momentumBullet, bulletA, bulletZ);
486  }
487 
488  if (verboseLevel > 2) G4cout << "Bullet: \n" << *bullet << G4endl;
489 
490  return true;
491 }
492 
493 
494 // Convert input nuclear target to Bertini internal object
495 
497  return createTarget(theNucleus.GetA_asInt(), theNucleus.GetZ_asInt());
498 }
499 
501  return createTarget(theNucleus->GetMassNumber(), theNucleus->GetCharge());
502 }
503 
505  if (A > 1) {
506  nucleusTarget.fill(A, Z);
508  } else {
509  hadronTarget.fill(0., (Z==1?proton:neutron));
510  target = &hadronTarget;
511  }
512 
513  if (verboseLevel > 2) G4cout << "Target: \n" << *target << G4endl;
514 
515  return true; // Right now, target never fails
516 }
517 
518 
519 // Convert Bertini particle to output (G4DynamicParticle)
520 
523  G4int outgoingType = iep.type();
524 
525  if (iep.quasi_deutron()) {
526  G4cerr << " ERROR: G4CascadeInterface incompatible particle type "
527  << outgoingType << G4endl;
528  return 0;
529  }
530 
531  // Copy local G4DynPart to public output (handle kaon mixing specially)
532  if (outgoingType == kaonZero || outgoingType == kaonZeroBar) {
533  G4ThreeVector momDir = iep.getMomentum().vect().unit();
534  G4double ekin = iep.getKineticEnergy()*GeV; // Bertini -> G4 units
535 
537  if (G4UniformRand() > 0.5) pd = G4KaonZeroLong::Definition();
538 
539  return new G4DynamicParticle(pd, momDir, ekin);
540  } else {
541  return new G4DynamicParticle(iep.getDynamicParticle());
542  }
543 
544  return 0; // Should never get here!
545 }
546 
549  if (verboseLevel > 2) {
550  G4cout << " Nuclei fragment: \n" << inuc << G4endl;
551  }
552 
553  // Copy local G4DynPart to public output
554  return new G4DynamicParticle(inuc.getDynamicParticle());
555 }
556 
557 
558 // Transfer Bertini internal final state to hadronics interface
559 
561  if (verboseLevel > 1)
562  G4cout << " >>> G4CascadeInterface::copyOutputToHadronicResult" << G4endl;
563 
564  const std::vector<G4InuclNuclei>& outgoingNuclei = output->getOutgoingNuclei();
565  const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
566 
569 
570  // Get outcoming particles
571  if (!particles.empty()) {
572  particleIterator ipart = particles.begin();
573  for (; ipart != particles.end(); ipart++) {
575  }
576  }
577 
578  // get nuclei fragments
579  if (!outgoingNuclei.empty()) {
580  nucleiIterator ifrag = outgoingNuclei.begin();
581  for (; ifrag != outgoingNuclei.end(); ifrag++) {
583  }
584  }
585 }
586 
588  if (verboseLevel > 1)
589  G4cout << " >>> G4CascadeInterface::copyOutputToReactionProducts" << G4endl;
590 
591  const std::vector<G4InuclElementaryParticle>& particles = output->getOutgoingParticles();
592  const std::vector<G4InuclNuclei>& fragments = output->getOutgoingNuclei();
593 
595 
596  G4ReactionProduct* rp = 0; // Buffers to create outgoing tracks
597  G4DynamicParticle* dp = 0;
598 
599  // Get outcoming particles
600  if (!particles.empty()) {
601  particleIterator ipart = particles.begin();
602  for (; ipart != particles.end(); ipart++) {
603  rp = new G4ReactionProduct;
604  dp = makeDynamicParticle(*ipart);
605  (*rp) = (*dp); // This does all the necessary copying
606  propResult->push_back(rp);
607  delete dp;
608  }
609  }
610 
611  // get nuclei fragments
612  if (!fragments.empty()) {
613  nucleiIterator ifrag = fragments.begin();
614  for (; ifrag != fragments.end(); ifrag++) {
615  rp = new G4ReactionProduct;
616  dp = makeDynamicParticle(*ifrag);
617  (*rp) = (*dp); // This does all the necessary copying
618  propResult->push_back(rp);
619  delete dp;
620  }
621  }
622 
623  return propResult;
624 }
625 
626 
627 // Report violations of conservation laws in original frame
628 
631 
632  if (verboseLevel > 2) {
633  if (!balance->baryonOkay()) {
634  G4cerr << "ERROR: no baryon number conservation, sum of baryons = "
635  << balance->deltaB() << G4endl;
636  }
637 
638  if (!balance->chargeOkay()) {
639  G4cerr << "ERROR: no charge conservation, sum of charges = "
640  << balance->deltaQ() << G4endl;
641  }
642 
643  if (std::abs(balance->deltaKE()) > 0.01 ) { // GeV
644  G4cerr << "Kinetic energy conservation violated by "
645  << balance->deltaKE() << " GeV" << G4endl;
646  }
647 
648  G4double eInit = bullet->getEnergy() + target->getEnergy();
649  G4double eFinal = eInit + balance->deltaE();
650 
651  G4cout << "Initial energy " << eInit << " final energy " << eFinal
652  << "\nTotal energy conservation at level "
653  << balance->deltaE() * GeV << " MeV" << G4endl;
654 
655  if (balance->deltaKE() > 5.0e-5 ) { // 0.05 MeV
656  G4cerr << "FATAL ERROR: kinetic energy created "
657  << balance->deltaKE() * GeV << " MeV" << G4endl;
658  }
659  }
660 }
661 
662 
663 // Evaluate whether any outgoing particles penetrated Coulomb barrier
664 
666  G4bool violated = false; // by default coulomb analysis is OK
667 
668  const G4double coulumbBarrier = 8.7 * MeV/GeV; // Bertini uses GeV
669 
670  const std::vector<G4InuclElementaryParticle>& p =
672 
673  for (particleIterator ipart=p.begin(); ipart != p.end(); ipart++) {
674  if (ipart->type() == proton) {
675  violated |= (ipart->getKineticEnergy() < coulumbBarrier);
676  }
677  }
678 
679  return violated;
680 }
681 
682 // Check whether inelastic collision on proton failed
683 
685  const std::vector<G4InuclElementaryParticle>& out =
687 
688 #ifdef G4CASCADE_DEBUG_INTERFACE
689  // Report on all retry conditions, in order of return logic
690  G4cout << " retryInelasticProton: number of Tries "
691  << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
692  << "\n retryInelasticProton: AND collision type ";
693  if (out.empty()) G4cout << "FAILED" << G4endl;
694  else {
695  G4cout << (out.size() == 2 ? "ELASTIC (t)" : "INELASTIC (f)")
696  << "\n retryInelasticProton: AND Leading particles bullet "
697  << (out.size() >= 2 &&
698  (out[0].getDefinition() == bullet->getDefinition() ||
699  out[1].getDefinition() == bullet->getDefinition())
700  ? "YES (t)" : "NO (f)")
701  << G4endl;
702  }
703 #endif
704 
705  return ( (numberOfTries < maximumTries) &&
706  (out.empty() ||
707  (out.size() == 2 &&
708  (out[0].getDefinition() == bullet->getDefinition() ||
709  out[1].getDefinition() == bullet->getDefinition())))
710  );
711 }
712 
713 // Check whether generic inelastic collision failed
714 // NOTE: some conditions are set by compiler flags
715 
717  // Quantities necessary for conditional block below
720 
721  const G4ParticleDefinition* firstOut = (npart == 0) ? 0 :
722  output->getOutgoingParticles().begin()->getDefinition();
723 
724 #ifdef G4CASCADE_DEBUG_INTERFACE
725  // Report on all retry conditions, in order of return logic
726  G4cout << " retryInelasticNucleus: numberOfTries "
727  << ((numberOfTries < maximumTries) ? "RETRY (t)" : "EXIT (f)")
728  << "\n retryInelasticNucleus: AND outputParticles "
729  << ((npart != 0) ? "NON-ZERO (t)" : "EMPTY (f)")
730 #ifdef G4CASCADE_COULOMB_DEV
731  << "\n retryInelasticNucleus: AND coulombBarrier (COULOMB_DEV) "
732  << (coulombBarrierViolation() ? "VIOLATED (t)" : "PASSED (f)")
733  << "\n retryInelasticNucleus: AND collision type (COULOMB_DEV) "
734  << ((npart+nfrag > 2) ? "INELASTIC (t)" : "ELASTIC (f)")
735 #else
736  << "\n retryInelasticNucleus: AND collsion type "
737  << ((npart+nfrag < 3) ? "ELASTIC (t)" : "INELASTIC (f)")
738  << "\n retryInelasticNucleus: AND Leading particle bullet "
739  << ((firstOut == bullet->getDefinition()) ? "YES (t)" : "NO (f)")
740 #endif
741  << "\n retryInelasticNucleus: OR conservation "
742  << (!balance->okay() ? "FAILED (t)" : "PASSED (f)")
743  << G4endl;
744 #endif
745 
746  return ( (numberOfTries < maximumTries) &&
747  ( ((npart != 0) &&
748 #ifdef G4CASCADE_COULOMB_DEV
749  (coulombBarrierViolation() && npart+nfrag > 2)
750 #else
751  (npart+nfrag < 3 && firstOut == bullet->getDefinition())
752 #endif
753  )
754 #ifndef G4CASCADE_SKIP_ECONS
755  || (!balance->okay())
756 #endif
757  )
758  );
759 }
760 
761 
762 // Terminate job in case of persistent non-conservation
763 // FIXME: Need to migrate to G4ExceptionDescription
764 
766  // NOTE: Once G4HadronicException is changed, use the following line!
767  // G4ExceptionDescription errInfo;
768  std::ostream& errInfo = G4cerr;
769 
770  errInfo << " >>> G4CascadeInterface has non-conserving"
771  << " cascade after " << numberOfTries << " attempts." << G4endl;
772 
773  G4String throwMsg = "G4CascadeInterface - ";
774  if (!balance->energyOkay()) {
775  throwMsg += "Energy";
776  errInfo << " Energy conservation violated by " << balance->deltaE()
777  << " GeV (" << balance->relativeE() << ")" << G4endl;
778  }
779 
780  if (!balance->momentumOkay()) {
781  throwMsg += "Momentum";
782  errInfo << " Momentum conservation violated by " << balance->deltaP()
783  << " GeV/c (" << balance->relativeP() << ")" << G4endl;
784  }
785 
786  if (!balance->baryonOkay()) {
787  throwMsg += "Baryon number";
788  errInfo << " Baryon number violated by " << balance->deltaB() << G4endl;
789  }
790 
791  if (!balance->chargeOkay()) {
792  throwMsg += "Charge";
793  errInfo << " Charge conservation violated by " << balance->deltaQ()
794  << G4endl;
795  }
796 
797  errInfo << " Final event output, for debugging:\n"
798  << " Bullet: \n" << *bullet << G4endl
799  << " Target: \n" << *target << G4endl;
800  output->printCollisionOutput(errInfo);
801 
802  throwMsg += " non-conservation. More info in output.";
803  throw G4HadronicException(__FILE__, __LINE__, throwMsg); // Job ends here!
804 }
static G4Dineutron * Definition()
Definition: G4Dineutron.cc:68
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
const XML_Char * name
Definition: expat.h:151
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
G4DynamicParticle * makeDynamicParticle(const G4InuclElementaryParticle &iep) const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4int GetAtomicNumber() const
const G4LorentzVector & Get4Momentum() const
const XML_Char * target
Definition: expat.h:268
static constexpr double MeV
Definition: G4SIunits.hh:214
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
Int_t ipart
Definition: Style.C:10
void fill(G4int ityp, Model model=DefaultModel)
static G4Diproton * Definition()
Definition: G4Diproton.cc:68
const G4HadProjectile * GetPrimaryProjectile() const
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
static constexpr double perCent
Definition: G4SIunits.hh:332
void rescatter(G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
const G4String & GetParticleName() const
void usePreCompoundDeexcitation()
const G4DynamicParticle & getDynamicParticle() const
const G4ParticleDefinition * GetDefinition() const
G4bool IsMasterThread()
Definition: G4Threading.cc:130
static G4KaonZeroLong * Definition()
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
static void saveEngineStatus(const char filename[]="Config.conf")
Definition: Random.cc:275
static void DumpConfiguration(std::ostream &os)
const G4ParticleDefinition * getDefinition() const
G4InuclElementaryParticle hadronBullet
virtual void DumpConfiguration(std::ostream &outFile) const
G4HadFinalState * NoInteraction(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4bool createTarget(G4Nucleus &theNucleus)
Float_t Z
void setVerboseLevel(G4int verbose=0)
void SetEnergyChange(G4double anEnergy)
void useCascadeDeexcitation()
static const G4CascadeChannel * GetTable(G4int initialState)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
Int_t npart
Definition: Style.C:7
G4double getKineticEnergy() const
G4CascadeInterface(const G4String &name="BertiniCascade")
void SetEnergyMomentumCheckLevels(G4double relativeLevel, G4double absoluteLevel)
void SetVerboseLevel(G4int value)
void setLimits(G4double relative, G4double absolute)
const G4String randomFile
#define G4UniformRand()
Definition: Randomize.hh:53
G4int numberOfOutgoingNuclei() const
double A(double temperature)
void setVerboseLevel(G4int verbose)
G4bool retryInelasticProton() const
virtual G4int GetCharge()=0
G4InuclNuclei nucleusTarget
G4CascadeCheckBalance * balance
G4bool retryInelasticNucleus() const
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4GLOB_DLL std::ostream G4cerr
G4int numberOfOutgoingParticles() const
const G4LorentzVector & Get4Momentum() const
virtual G4int GetMassNumber()=0
double rho() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
const G4ParticleDefinition * GetDefinition() const
G4CollisionOutput * output
G4double GetFormationTime() const
int G4int
Definition: G4Types.hh:78
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
static G4bool usePreCompound()
G4ReactionProductVector * copyOutputToReactionProducts()
static G4KaonZeroShort * Definition()
G4int GetAtomicMass() const
G4LorentzVector getMomentum() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4GLOB_DLL std::ostream G4cout
virtual void setVerboseLevel(G4int verbose=0)
G4InuclParticle * bullet
G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
G4bool createBullet(const G4HadProjectile &aTrack)
Hep3Vector vect() const
virtual void ModelDescription(std::ostream &outFile) const
const G4ThreeVector & GetPosition() const
void SetVerboseLevel(G4int verbose)
void printCollisionOutput(std::ostream &os=G4cout) const
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
G4bool coulombBarrierViolation() const
G4double getEnergy() const
static G4UnboundPN * Definition()
Definition: G4UnboundPN.cc:67
static constexpr double GeV
Definition: G4SIunits.hh:217
G4InuclElementaryParticle hadronTarget
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4InuclCollider * collider
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
void SetStatusChange(G4HadFinalStateStatus aS)
G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
G4InuclNuclei nucleusBullet