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