Geant4-11
G4IntraNucleiCascader.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// 20100307 M. Kelsey -- Bug fix: momentum_out[0] should be momentum_out.e()
29// 20100309 M. Kelsey -- Eliminate some unnecessary std::pow()
30// 20100407 M. Kelsey -- Pass "all_particles" as argument to initializeCascad,
31// following recent change to G4NucleiModel.
32// 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
33// 20100517 M. Kelsey -- Inherit from common base class, make other colliders
34// simple data members
35// 20100616 M. Kelsey -- Add reporting of final residual particle
36// 20100617 M. Kelsey -- Remove "RUN" preprocessor flag and all "#else" code,
37// pass verbosity to collider. Make G4NucleiModel a data member,
38// instead of creating and deleting on every cycle.
39// 20100620 M. Kelsey -- Improved diagnostic messages. Simplify kinematics
40// of recoil nucleus.
41// 20100622 M. Kelsey -- Use local "bindingEnergy()" to call through.
42// 20100623 M. Kelsey -- Undo G4NucleiModel change from 0617. Does not work
43// properly across multiple interactions.
44// 20100627 M. Kelsey -- Protect recoil nucleus energy from floating roundoff
45// by setting small +ve or -ve values to zero.
46// 20100701 M. Kelsey -- Let excitation energy be handled by G4InuclNuclei,
47// allow for ground-state recoil (goodCase == true for Eex==0.)
48// 20100702 M. Kelsey -- Negative energy recoil should be rejected
49// 20100706 D. Wright -- Copy "abandoned" cparticles to output list, copy
50// mesonic "excitons" to output list; should be absorbed, fix up
51// diagnostic messages.
52// 20100713 M. Kelsey -- Add more diagnostics for Dennis' changes.
53// 20100714 M. Kelsey -- Switch to new G4CascadeColliderBase class, remove
54// sanity check on afin/zfin (not valid).
55// 20100715 M. Kelsey -- Add diagnostic for ekin_in vs. actual ekin; reduce
56// KE across Coulomb barrier. Rearrange end-of-loop if blocks,
57// add conservation check at end.
58// 20100716 M. Kelsey -- Eliminate inter_case; use base-class functionality.
59// Add minimum-fragment requirement for recoil, in order to
60// allow for momentum balancing
61// 20100720 M. Kelsey -- Make EPCollider pointer member
62// 20100721 M. Kelsey -- Turn on conservation checks unconditionally (override
63// new G4CASCADE_CHECK_ECONS setting
64// 20100722 M. Kelsey -- Move cascade output buffers to .hh file
65// 20100728 M. Kelsey -- Make G4NucleiModel data member for persistence,
66// delete colliders in destructor
67// 20100906 M. Kelsey -- Hide "non-physical fragment" behind verbose flag
68// 20100907 M. Kelsey -- Add makeResidualFragment function to create object
69// 20100909 M. Kelsey -- Remove all local "fragment" stuff, use RecoilMaker.
70// move goodCase() to RecoilMaker.
71// 20100910 M. Kelsey -- Use RecoilMaker::makeRecoilFragment().
72// 20100915 M. Kelsey -- Define functions to deal with trapped particles,
73// move the exciton container to a data member
74// 20100916 M. Kelsey -- Put decay photons directly onto output list
75// 20100921 M. Kelsey -- Migrate to RecoilMaker::makeRecoilNuclei().
76// 20100924 M. Kelsey -- Minor shuffling of post-cascade recoil building.
77// Create G4Fragment for recoil and store in output.
78// 20110131 M. Kelsey -- Move "momentum_in" calculation inside verbosity
79// 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
80// 20110224 M. Kelsey -- Add ::rescatter() function which takes a list of
81// pre-existing secondaries as input. Split ::collide() into
82// separate utility functions. Move cascade parameters to static
83// data members. Add setVerboseLevel().
84// 20110302 M. Kelsey -- Move G4NucleiModel::printModel() call to G4NucleiModel
85// 20110303 M. Kelsey -- Add more cascade functions to support rescattering
86// 20110304 M. Kelsey -- Get original Propagate() arguments here in rescatter()
87// and convert to particles, nuclei and G4NucleiModel state.
88// 20110308 M. Kelsey -- Don't put recoiling fragment onto output list any more
89// 20110308 M. Kelsey -- Decay unstable hadrons from pre-cascade, use daughters
90// 20110324 M. Kelsey -- Get locations of hit nuclei in ::rescatter(), pass
91// to G4NucleiModel::reset().
92// 20110404 M. Kelsey -- Reduce maximum number of retries to 100, reflection
93// cut to 50.
94// 20110721 M. Kelsey -- Put unusable pre-cascade particles directly on output,
95// do not decay.
96// 20110722 M. Kelsey -- Deprecate "output_particles" list in favor of using
97// output directly (will help with pre-cascade issues).
98// 20110801 M. Kelsey -- Use G4Inucl(Particle)::fill() functions to reduce
99// creation of temporaries. Add local target buffer for
100// rescattering, to avoid memory leak.
101// 20110808 M. Kelsey -- Pass buffer to generateParticleFate() to avoid copy
102// 20110919 M. Kelsey -- Add optional final-state clustering, controlled (for
103// now) with compiler flag G4CASCADE_DO_COALESCENCE
104// 20110922 M. Kelsey -- Follow migrations G4InuclParticle::print(ostream&)
105// and G4CascadParticle::print(ostream&); drop Q,B printing
106// 20110926 M. Kelsey -- Replace compiler flag with one-time envvar in ctor
107// for final-state clustering.
108// 20111003 M. Kelsey -- Prepare for gamma-N interactions by checking for
109// final-state tables instead of particle "isPhoton()"
110// 20120521 A. Ribon -- Specify mass when decay trapped particle.
111// 20120822 M. Kelsey -- Move envvars to G4CascadeParameters.
112// 20121205 M. Kelsey -- In processSecondary(), set generation to 1, as these
113// particles are not true projectiles, but already embedded.
114// 20130304 M. Kelsey -- Use new G4CascadeHistory to dump cascade structure
115// 20140310 M. Kelsey -- (Bug #1584) Release memory allocated by DecayIt()
116// 20140409 M. Kelsey -- Use const G4ParticleDefinition* everywhere
117// 20141204 M. Kelsey -- Add function to test for non-interacting particles,
118// move those directly to output without propagating
119// 20150608 M. Kelsey -- Label all while loops as terminating.
120// 20150619 M. Kelsey -- Replace std::exp with G4Exp
121
122#include <algorithm>
123
125#include "G4SystemOfUnits.hh"
128#include "G4CascadeHistory.hh"
129#include "G4CascadeParameters.hh"
131#include "G4CascadParticle.hh"
132#include "G4CollisionOutput.hh"
133#include "G4DecayProducts.hh"
134#include "G4DecayTable.hh"
137#include "G4Exp.hh"
138#include "G4HadTmpUtil.hh"
140#include "G4InuclNuclei.hh"
143#include "G4KineticTrack.hh"
145#include "G4LorentzConvertor.hh"
146#include "G4Neutron.hh"
147#include "G4NucleiModel.hh"
149#include "G4Proton.hh"
150#include "G4V3DNucleus.hh"
151#include "Randomize.hh"
152
153using namespace G4InuclParticleNames;
154using namespace G4InuclSpecialFunctions;
155
156
157// Configuration parameters for cascade production
162
163
164typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
165
167 : G4CascadeColliderBase("G4IntraNucleiCascader"), model(new G4NucleiModel),
168 theElementaryParticleCollider(new G4ElementaryParticleCollider),
169 theRecoilMaker(new G4CascadeRecoilMaker), theClusterMaker(0),
170 theCascadeHistory(0), tnuclei(0), bnuclei(0), bparticle(0),
171 minimum_recoil_A(0.), coulombBarrier(0.),
172 nucleusTarget(new G4InuclNuclei),
173 protonTarget(new G4InuclElementaryParticle) {
176
179}
180
182 delete model;
184 delete theRecoilMaker;
185 delete theClusterMaker;
186 delete theCascadeHistory;
187 delete nucleusTarget;
188 delete protonTarget;
189}
190
193 model->setVerboseLevel(verbose);
196
197 // Optional functionality
200}
201
202
203
205 G4InuclParticle* target,
206 G4CollisionOutput& globalOutput) {
207 if (verboseLevel) G4cout << " >>> G4IntraNucleiCascader::collide " << G4endl;
208
209 if (!initialize(bullet, target)) return; // Load buffers and drivers
210
211 G4int itry = 0;
212 do { /* Loop checking 08.06.2015 MHK */
213 newCascade(++itry);
214 setupCascade();
216 } while (!finishCascade() && itry<itry_max);
217
218 // Report full structure of final cascade if requested
220
221 finalize(itry, bullet, target, globalOutput);
222}
223
224// For use with Propagate to preload a set of secondaries
225// FIXME: So far, we don't have any bullet information from Propagate!
226
228 G4KineticTrackVector* theSecondaries,
229 G4V3DNucleus* theNucleus,
230 G4CollisionOutput& globalOutput) {
231 if (verboseLevel)
232 G4cout << " >>> G4IntraNucleiCascader::rescatter " << G4endl;
233
234 G4InuclParticle* target = createTarget(theNucleus);
235 if (!initialize(bullet, target)) return; // Load buffers and drivers
236
237 G4int itry = 0;
238 do { /* Loop checking 08.06.2015 MHK */
239 newCascade(++itry);
240 preloadCascade(theNucleus, theSecondaries);
242 } while (!finishCascade() && itry<itry_max);
243
244 // Report full structure of final cascade if requested
246
247 finalize(itry, bullet, target, globalOutput);
248}
249
250
252 G4InuclParticle* target) {
253 if (verboseLevel>1)
254 G4cout << " >>> G4IntraNucleiCascader::initialize " << G4endl;
255
256 // Configure processing modules
258
259 interCase.set(bullet,target); // Classify collision type
260
261 if (verboseLevel > 3) {
263 << *interCase.getTarget() << G4endl;
264 }
265
266 // Bullet may be nucleus or simple particle
267 bnuclei = dynamic_cast<G4InuclNuclei*>(interCase.getBullet());
269
270 if (!bnuclei && !bparticle) {
271 G4cerr << " G4IntraNucleiCascader: projectile is not a valid particle."
272 << G4endl;
273 return false;
274 }
275
276 // Target _must_ be nucleus
277 tnuclei = dynamic_cast<G4InuclNuclei*>(interCase.getTarget());
278 if (!tnuclei) {
279 if (verboseLevel)
280 G4cerr << " Target is not a nucleus. Abandoning." << G4endl;
281 return false;
282 }
283
285 coulombBarrier = 0.00126*tnuclei->getZ() / (1.+G4cbrt(tnuclei->getA()));
286
287 // Energy/momentum conservation usually requires a recoiling nuclear fragment
288 // This cut will be increased on each "itry" if momentum could not balance.
289 minimum_recoil_A = 0.;
290
291 if (verboseLevel > 3) {
292 G4LorentzVector momentum_in = bullet->getMomentum() + target->getMomentum();
293 G4cout << " intitial momentum E " << momentum_in.e() << " Px "
294 << momentum_in.x() << " Py " << momentum_in.y() << " Pz "
295 << momentum_in.z() << G4endl;
296 }
297
298 return true;
299}
300
301// Re-initialize buffers for new attempt at cascade
302
304 if (verboseLevel > 1) {
305 G4cout << " IntraNucleiCascader itry " << itry << " inter_case "
306 << interCase.code() << G4endl;
307 }
308
309 model->reset(); // Start new cascade process
310 output.reset();
311 new_cascad_particles.clear();
313 cascad_particles.clear(); // List of initial secondaries
314
316}
317
318
319// Load initial cascade using nuclear-model calculations
320
322 if (verboseLevel > 1)
323 G4cout << " >>> G4IntraNucleiCascader::setupCascade" << G4endl;
324
325 if (interCase.hadNucleus()) { // particle with nuclei
326 if (verboseLevel > 3)
327 G4cout << " bparticle charge " << bparticle->getCharge()
328 << " baryon number " << bparticle->baryon() << G4endl;
329
331 } else { // nuclei with nuclei
332 G4int ab = bnuclei->getA();
333 G4int zb = bnuclei->getZ();
334
335 G4NucleiModel::modelLists all_particles; // Buffer to receive lists
336 model->initializeCascad(bnuclei, tnuclei, all_particles);
337
338 cascad_particles = all_particles.first;
339 output.addOutgoingParticles(all_particles.second);
340
341 if (cascad_particles.size() == 0) { // compound nuclei
342 G4int i;
343
344 for (i = 0; i < ab; i++) {
345 G4int knd = i < zb ? 1 : 2;
347 };
348
349 G4int ihn = G4int(2 * (ab-zb) * inuclRndm() + 0.5);
350 G4int ihz = G4int(2 * zb * inuclRndm() + 0.5);
351
352 for (i = 0; i < ihn; i++) theExitonConfiguration.incrementHoles(2);
353 for (i = 0; i < ihz; i++) theExitonConfiguration.incrementHoles(1);
354 }
355 } // if (interCase ...
356}
357
358
359// Generate one possible cascade (all secondaries, etc.)
360
362 if (verboseLevel>1) G4cout << " generateCascade " << G4endl;
363
364 /* Loop checking 08.06.2015 MHK */
365 G4int iloop = 0;
366 while (!cascad_particles.empty() && !model->empty()) {
367 iloop++;
368
369 if (verboseLevel > 2) {
370 G4cout << " Iteration " << iloop << ": Number of cparticles "
371 << cascad_particles.size() << " last one: \n"
372 << cascad_particles.back() << G4endl;
373 }
374
375 // Record incident particle first, to get history ID
376 if (theCascadeHistory) {
378 if (verboseLevel > 2) {
379 G4cout << " active cparticle got history ID "
380 << cascad_particles.back().getHistoryId() << G4endl;
381 }
382 }
383
384 // If non-interacting particle, move directly to output
386 if (verboseLevel > 2)
387 G4cout << " particle is non-interacting; moving to output" << G4endl;
388
389 output.addOutgoingParticle(cascad_particles.back().getParticle());
390 cascad_particles.pop_back();
391 continue;
392 }
393
394 // Generate interaction with nucleon
398
399 // Record interaction for later reporting (if desired)
402
403 if (verboseLevel > 2) {
404 G4cout << " After generate fate: New particles "
405 << new_cascad_particles.size() << G4endl
406 << " Discarding last cparticle from list " << G4endl;
407 }
408
409 cascad_particles.pop_back();
410
411 // handle the result of a new step
412
413 if (new_cascad_particles.size() == 1) { // last particle goes without interaction
414 const G4CascadParticle& currentCParticle = new_cascad_particles[0];
415
416 if (model->stillInside(currentCParticle)) {
417 if (verboseLevel > 3)
418 G4cout << " particle still inside nucleus " << G4endl;
419
420 if (currentCParticle.getNumberOfReflections() < reflection_cut &&
421 model->worthToPropagate(currentCParticle)) {
422 if (verboseLevel > 3) G4cout << " continue reflections " << G4endl;
423 cascad_particles.push_back(currentCParticle);
424 } else {
425 processTrappedParticle(currentCParticle);
426 } // reflection or exciton
427
428 } else { // particle about to leave nucleus - check for Coulomb barrier
429 if (verboseLevel > 3) G4cout << " possible escape " << G4endl;
430
431 const G4InuclElementaryParticle& currentParticle =
432 currentCParticle.getParticle();
433
434 G4double KE = currentParticle.getKineticEnergy();
435 G4double mass = currentParticle.getMass();
436 G4double Q = currentParticle.getCharge();
437
438 if (verboseLevel > 3)
439 G4cout << " KE " << KE << " barrier " << Q*coulombBarrier << G4endl;
440
441 if (KE < Q*coulombBarrier) {
442 // Calculate barrier penetration
443 G4double CBP = 0.0;
444
445 // if (KE > 0.0001) CBP = std::exp(-0.00126*tnuclei->getZ()*0.25*
446 // (1./KE - 1./coulombBarrier));
447 if (KE > 0.0001) CBP = G4Exp(-0.0181*0.5*tnuclei->getZ()*
448 (1./KE - 1./coulombBarrier)*
449 std::sqrt(mass*(coulombBarrier-KE)) );
450
451 if (G4UniformRand() < CBP) {
452 if (verboseLevel > 3)
453 G4cout << " tunneled\n" << currentParticle << G4endl;
454
455 // Tunnelling through barrier leaves KE unchanged
456 output.addOutgoingParticle(currentParticle);
457 } else {
458 processTrappedParticle(currentCParticle);
459 }
460 } else {
461 output.addOutgoingParticle(currentParticle);
462
463 if (verboseLevel > 3)
464 G4cout << " Goes out\n" << output.getOutgoingParticles().back()
465 << G4endl;
466 }
467 }
468 } else { // interaction
469 if (verboseLevel > 3)
470 G4cout << " interacted, adding new to list " << G4endl;
471
473 new_cascad_particles.begin(),
475
476 std::pair<G4int, G4int> holes = model->getTypesOfNucleonsInvolved();
477 if (verboseLevel > 3)
478 G4cout << " adding new exciton holes " << holes.first << ","
479 << holes.second << G4endl;
480
482
483 if (holes.second > 0)
485 } // if (new_cascad_particles ...
486
487 // Evaluate nuclear residue
490
492 if (verboseLevel > 2) {
493 G4cout << " cparticles remaining " << cascad_particles.size()
494 << " nucleus (model) has "
495 << model->getNumberOfNeutrons() << " n, "
496 << model->getNumberOfProtons() << " p "
497 << " residual fragment A " << aresid << G4endl;
498 }
499
500 if (aresid <= minimum_recoil_A) return; // Must have minimum fragment
501 } // while cascade-list and model
502}
503
504
505// Conslidate results of cascade and evaluate success
506
508 if (verboseLevel > 1)
509 G4cout << " >>> G4IntraNucleiCascader::finishCascade ?" << G4endl;
510
511 // Add left-over cascade particles to output
513 cascad_particles.clear();
514
515 // Cascade is finished. Check if it's OK.
516 if (verboseLevel>3) {
517 G4cout << " G4IntraNucleiCascader finished" << G4endl;
519 }
520
521 // Apply cluster coalesence model to produce light ions
522 if (theClusterMaker) {
525
526 // Update recoil fragment after generating light ions
527 if (verboseLevel>3) G4cout << " Recomputing recoil fragment" << G4endl;
529 output);
530 if (verboseLevel>3) {
531 G4cout << " After cluster coalescence" << G4endl;
533 }
534 }
535
536 // Use last created recoil fragment instead of re-constructing
539
540 // FIXME: Should we deal with unbalanced (0,0) case before rejecting?
541
542 // Sanity check before proceeding
544 if (verboseLevel > 1)
545 G4cerr << " Recoil nucleus is not physical: A=" << afin << " Z="
546 << zfin << G4endl;
547 return false; // Discard event and try again
548 }
549
551
552 if (verboseLevel > 1) {
553 G4cout << " afin " << afin << " zfin " << zfin << G4endl;
554 }
555
556 if (afin == 0) return true; // Whole event fragmented, exit
557
558 if (afin == 1) { // Add bare nucleon to particle list
559 G4int last_type = (zfin==1) ? 1 : 2; // proton=1, neutron=2
560
562 G4double mres = presid.m();
563
564 // Check for sensible kinematics
565 if (mres-mass < -small_ekin) { // Insufficient recoil energy
566 if (verboseLevel > 2) G4cerr << " unphysical recoil nucleon" << G4endl;
567 return false;
568 }
569
570 if (mres-mass > small_ekin) { // Too much extra energy
571 if (verboseLevel > 2)
572 G4cerr << " extra energy with recoil nucleon" << G4endl;
573
574 // FIXME: For now, we add the nucleon as unbalanced, and let
575 // "SetOnShell" fudge things. This should be abandoned.
576 }
577
578 G4InuclElementaryParticle last_particle(presid, last_type,
580
581 if (verboseLevel > 3) {
582 G4cout << " adding recoiling nucleon to output list\n"
583 << last_particle << G4endl;
584 }
585
586 output.addOutgoingParticle(last_particle);
587
588 // Update recoil to include residual nucleon
590 output);
591 }
592
593 // Process recoil fragment for consistency, exit or reject
596 if (std::abs(Eex) < quasielast_cut) {
597 if (verboseLevel > 3) {
598 G4cout << " quasi-elastic scatter with " << Eex << " MeV recoil"
599 << G4endl;
600 }
601
603 if (verboseLevel > 3) {
604 G4cout << " Eex reset to " << theRecoilMaker->getRecoilExcitation()
605 << G4endl;
606 }
607 }
608 }
609
612
614 if (!recoilFrag) {
615 G4cerr << "Got null pointer for recoil fragment!" << G4endl;
616 return false;
617 }
618
619 if (verboseLevel > 2)
620 G4cout << " adding recoil fragment to output list" << G4endl;
621
622 output.addRecoilFragment(*recoilFrag);
623 }
624
625 // Put final-state particles in "leading order" for return
626 std::vector<G4InuclElementaryParticle>& opart = output.getOutgoingParticles();
627 std::sort(opart.begin(), opart.end(), G4ParticleLargerEkin());
628
629 // Adjust final state to balance momentum and energy if necessary
634
635 if (output.acceptable()) return true;
636 else if (verboseLevel>2) G4cerr << " Cascade setOnShell failed." << G4endl;
637 }
638
639 // Cascade not physically reasonable
640 if (afin <= minimum_recoil_A && minimum_recoil_A < tnuclei->getA()) {
642 if (verboseLevel > 3) {
643 G4cout << " minimum recoil fragment increased to A " << minimum_recoil_A
644 << G4endl;
645 }
646 }
647
648 if (verboseLevel>2) G4cerr << " Cascade failed. Retrying..." << G4endl;
649 return false;
650}
651
652
653// Transfer finished cascade to return buffer
654
655void
657 G4InuclParticle* target,
658 G4CollisionOutput& globalOutput) {
659 if (itry >= itry_max) {
660 if (verboseLevel) {
661 G4cout << " IntraNucleiCascader-> no inelastic interaction after "
662 << itry << " attempts " << G4endl;
663 }
664
665 output.trivialise(bullet, target);
666 } else if (verboseLevel) {
667 G4cout << " IntraNucleiCascader output after trials " << itry << G4endl;
668 }
669
670 // Copy final generated cascade to output buffer for return
671 globalOutput.add(output);
672}
673
674
675// Create simple nucleus from rescattering target
676
679 G4int theNucleusA = theNucleus->GetMassNumber();
680 G4int theNucleusZ = theNucleus->GetCharge();
681
682 if (theNucleusA > 1) {
683 if (!nucleusTarget) nucleusTarget = new G4InuclNuclei; // Just in case
684 nucleusTarget->fill(0., theNucleusA, theNucleusZ, 0.);
685 return nucleusTarget;
686 } else {
688 protonTarget->fill(0., (theNucleusZ==1)?proton:neutron);
689 return protonTarget;
690 }
691
692 return 0; // Can never actually get here
693}
694
695// Copy existing (rescattering) cascade for propagation
696
697void
699 G4KineticTrackVector* theSecondaries) {
700 if (verboseLevel > 1)
701 G4cout << " >>> G4IntraNucleiCascader::preloadCascade" << G4endl;
702
703 copyWoundedNucleus(theNucleus); // Update interacted nucleon counts
704 copySecondaries(theSecondaries); // Copy original to internal list
705}
706
708 if (verboseLevel > 1)
709 G4cout << " >>> G4IntraNucleiCascader::copyWoundedNucleus" << G4endl;
710
711 // Loop over nucleons and count hits as exciton holes
713 hitNucleons.clear();
714 if (theNucleus->StartLoop()) {
715 G4Nucleon* nucl = 0;
716 G4int nuclType = 0;
717 /* Loop checking 08.06.2015 MHK */
718 while ((nucl = theNucleus->GetNextNucleon())) {
719 if (nucl->AreYouHit()) { // Found previously interacted nucleon
722 hitNucleons.push_back(nucl->GetPosition());
723 }
724 }
725 }
726
727 if (verboseLevel > 3)
728 G4cout << " nucleus has " << theExitonConfiguration.neutronHoles
729 << " neutrons hit, " << theExitonConfiguration.protonHoles
730 << " protons hit" << G4endl;
731
732 // Preload nuclear model with confirmed hits, including locations
735}
736
737void
739 if (verboseLevel > 1)
740 G4cout << " >>> G4IntraNucleiCascader::copySecondaries" << G4endl;
741
742 for (size_t i=0; i<secondaries->size(); i++) {
743 if (verboseLevel > 3) G4cout << " processing secondary " << i << G4endl;
744
745 processSecondary((*secondaries)[i]); // Copy to cascade or to output
746 }
747
748 // Sort list of secondaries to put leading particle first
749 std::sort(cascad_particles.begin(), cascad_particles.end(),
751
752 if (verboseLevel > 2) {
753 G4cout << " Original list of " << secondaries->size() << " secondaries"
754 << " produced " << cascad_particles.size() << " cascade, "
755 << output.numberOfOutgoingParticles() << " released particles, "
756 << output.numberOfOutgoingNuclei() << " fragments" << G4endl;
757 }
758}
759
760
761// Convert from pre-cascade secondary to local version
762
764 if (!ktrack) return; // Sanity check
765
766 // Get particle type to determine whether to keep or release
767 const G4ParticleDefinition* kpd = ktrack->GetDefinition();
768 if (!kpd) return;
769
771 if (!ktype) {
772 releaseSecondary(ktrack);
773 return;
774 }
775
776 if (verboseLevel > 1) {
777 G4cout << " >>> G4IntraNucleiCascader::processSecondary "
778 << kpd->GetParticleName() << G4endl;
779 }
780
781 // Allocate next local particle in buffer and fill
782 cascad_particles.resize(cascad_particles.size()+1); // Like push_back();
783 G4CascadParticle& cpart = cascad_particles.back();
784
785 // Convert momentum to Bertini internal units
786 cpart.getParticle().fill(ktrack->Get4Momentum()/GeV, ktype);
787 cpart.setGeneration(1);
788 cpart.setMovingInsideNuclei();
789 cpart.initializePath(0);
790
791 // Convert position units to Bertini's internal scale
792 G4ThreeVector cpos = ktrack->GetPosition()/model->getRadiusUnits();
793
794 cpart.updatePosition(cpos);
795 cpart.updateZone(model->getZone(cpos.mag()));
796
797 if (verboseLevel > 2)
798 G4cout << " Created cascade particle \n" << cpart << G4endl;
799}
800
801
802// Transfer unusable pre-cascade secondaries directly to output
803
805 const G4ParticleDefinition* kpd = ktrack->GetDefinition();
806
807 if (verboseLevel > 1) {
808 G4cout << " >>> G4IntraNucleiCascader::releaseSecondary "
809 << kpd->GetParticleName() << G4endl;
810 }
811
812 // Convert light ion into nucleus on fragment list
813 if (dynamic_cast<const G4Ions*>(kpd)) {
814 // Use resize() and fill() to avoid memory churn
816 G4InuclNuclei& inucl = output.getOutgoingNuclei().back();
817
818 inucl.fill(ktrack->Get4Momentum()/GeV,
819 kpd->GetAtomicMass(), kpd->GetAtomicNumber());
820 if (verboseLevel > 2)
821 G4cout << " Created pre-cascade fragment\n" << inucl << G4endl;
822 } else {
823 // Use resize() and fill() to avoid memory churn
826
827 // SPECIAL: Use G4PartDef directly, allowing unknown type code
828 ipart.fill(ktrack->Get4Momentum()/GeV, ktrack->GetDefinition());
829 if (verboseLevel > 2)
830 G4cout << " Created invalid pre-cascade particle\n" << ipart << G4endl;
831 }
832}
833
834
835// Convert particles which cannot escape into excitons (or eject/decay them)
836
839 const G4InuclElementaryParticle& trappedP = trapped.getParticle();
840
841 G4int xtype = trappedP.type();
842 if (verboseLevel > 3) G4cout << " exciton of type " << xtype << G4endl;
843
844 if (trappedP.nucleon()) { // normal exciton (proton or neutron)
847 return;
848 }
849
850 if (trappedP.hyperon()) { // Not nucleon, so must be hyperon
851 decayTrappedParticle(trapped);
853 return;
854 }
855
856 // non-standard exciton; release it
857 // FIXME: this is a meson, so need to absorb it
858 if (verboseLevel > 3) {
859 G4cout << " non-standard should be absorbed, now released\n"
860 << trapped << G4endl;
861 }
862
863 output.addOutgoingParticle(trappedP);
864}
865
866
867// Decay unstable trapped particles, and add secondaries to processing list
868
871 if (verboseLevel > 3)
872 G4cout << " unstable must be decayed in flight" << G4endl;
873
874 const G4InuclElementaryParticle& trappedP = trapped.getParticle();
875
876 G4DecayTable* unstable = trappedP.getDefinition()->GetDecayTable();
877 if (!unstable) { // No decay table; cannot decay!
878 if (verboseLevel > 3)
879 G4cerr << " no decay table! Releasing trapped particle" << G4endl;
880
881 output.addOutgoingParticle(trappedP);
882 return;
883 }
884
885 // Get secondaries from decay in particle's rest frame
886 G4DecayProducts* daughters = unstable->SelectADecayChannel()->DecayIt( trappedP.getDefinition()->GetPDGMass() );
887 if (!daughters) { // No final state; cannot decay!
888 if (verboseLevel > 3)
889 G4cerr << " no daughters! Releasing trapped particle" << G4endl;
890
891 output.addOutgoingParticle(trappedP);
892 return;
893 }
894
895 if (verboseLevel > 3)
896 G4cout << " " << daughters->entries() << " decay daughters" << G4endl;
897
898 // Convert secondaries to lab frame
899 G4double decayEnergy = trappedP.getEnergy();
900 G4ThreeVector decayDir = trappedP.getMomentum().vect().unit();
901 daughters->Boost(decayEnergy, decayDir);
902
903 // Put all the secondaries onto the list for propagation
904 const G4ThreeVector& decayPos = trapped.getPosition();
905 G4int zone = trapped.getCurrentZone();
906 G4int gen = trapped.getGeneration()+1;
907
908 for (G4int i=0; i<daughters->entries(); i++) {
909 G4DynamicParticle* idaug = (*daughters)[i];
910
912
913 // Propagate hadronic secondaries with known interactions (tables)
914 if (G4CascadeChannelTables::GetTable(idaugEP.type())) {
915 if (verboseLevel > 3) G4cout << " propagating " << idaugEP << G4endl;
916 cascad_particles.push_back(G4CascadParticle(idaugEP,decayPos,zone,0.,gen));
917 } else {
918 if (verboseLevel > 3) G4cout << " releasing " << idaugEP << G4endl;
920 }
921 }
922
923 delete daughters; // Clean up memory created by DecayIt()
924}
925
926
927// Test if particle is able to interact in nucleus
928
930particleCanInteract(const G4CascadParticle& cpart) const {
931 // If we have a lookup table for particle type on proton, it interacts
932 return (0 != G4CascadeChannelTables::GetTable(cpart.getParticle().type()));
933}
static const G4double ab
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
std::vector< G4InuclElementaryParticle >::iterator particleIterator
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
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
double mag() const
Hep3Vector vect() const
G4int getGeneration() const
void updateZone(G4int izone)
const G4InuclElementaryParticle & getParticle() const
void setGeneration(G4int gen)
void updatePosition(const G4ThreeVector &pos)
void setMovingInsideNuclei(G4bool isMovingIn=true)
G4int getNumberOfReflections() const
void initializePath(G4double npath)
G4int getCurrentZone() const
const G4ThreeVector & getPosition() const
static const G4CascadeChannel * GetTable(G4int initialState)
void FindClusters(G4CollisionOutput &finalState)
void setVerboseLevel(G4int verbose)
virtual void setVerboseLevel(G4int verbose=0)
void Print(std::ostream &os) const
G4int AddVertex(G4CascadParticle &cpart, std::vector< G4CascadParticle > &daug)
void DropEntry(const G4CascadParticle &cpart)
G4int AddEntry(G4CascadParticle &cpart)
void setVerboseLevel(G4int verbose=0)
static G4bool doCoalescence()
static G4bool showHistory()
void addExcitonConfiguration(const G4ExitonConfiguration exciton)
void setTolerance(G4double tolerance)
const G4LorentzVector & getRecoilMomentum() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4double getRecoilExcitation() const
void setRecoilExcitation(G4double Eexc)
G4Fragment * makeRecoilFragment()
G4int numberOfOutgoingParticles() const
void addRecoilFragment(const G4Fragment *aFragment)
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
void printCollisionOutput(std::ostream &os=G4cout) const
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void setOnShell(G4InuclParticle *bullet, G4InuclParticle *target)
void setVerboseLevel(G4int verbose)
G4int numberOfOutgoingNuclei() const
G4bool acceptable() const
void add(const G4CollisionOutput &right)
void trivialise(G4InuclParticle *bullet, G4InuclParticle *target)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
G4int entries() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:82
G4InuclParticle * getBullet() const
void set(G4InuclParticle *part1, G4InuclParticle *part2)
G4bool hadNucleus() const
G4InuclParticle * getTarget() const
void processTrappedParticle(const G4CascadParticle &trapped)
G4CascadeRecoilMaker * theRecoilMaker
void releaseSecondary(const G4KineticTrack *aSecondary)
std::vector< G4ThreeVector > hitNucleons
void copyWoundedNucleus(G4V3DNucleus *theNucleus)
static const G4double small_ekin
G4bool particleCanInteract(const G4CascadParticle &cpart) const
static const G4double quasielast_cut
std::vector< G4CascadParticle > cascad_particles
void rescatter(G4InuclParticle *bullet, G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4CollisionOutput &globalOutput)
G4CascadeCoalescence * theClusterMaker
void setVerboseLevel(G4int verbose=0)
G4CascadeHistory * theCascadeHistory
void decayTrappedParticle(const G4CascadParticle &trapped)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
std::vector< G4CascadParticle > new_cascad_particles
G4ExitonConfiguration theExitonConfiguration
G4bool initialize(G4InuclParticle *bullet, G4InuclParticle *target)
G4InuclParticle * createTarget(G4V3DNucleus *theNucleus)
void preloadCascade(G4V3DNucleus *theNucleus, G4KineticTrackVector *theSecondaries)
G4InuclElementaryParticle * protonTarget
void copySecondaries(G4KineticTrackVector *theSecondaries)
G4InuclElementaryParticle * bparticle
void processSecondary(const G4KineticTrack *aSecondary)
G4ElementaryParticleCollider * theElementaryParticleCollider
static const G4int reflection_cut
void finalize(G4int itry, G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
void fill(G4int ityp, Model model=DefaultModel)
static G4double getParticleMass(G4int type)
G4int getZ() const
void fill(G4int a, G4int z, G4double exc=0., Model model=DefaultModel)
G4int getA() const
const G4ParticleDefinition * getDefinition() const
G4double getKineticEnergy() const
G4double getMass() const
G4LorentzVector getMomentum() const
G4double getEnergy() const
G4double getCharge() const
Definition: G4Ions.hh:52
const G4ThreeVector & GetPosition() const
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4int getZone(G4double r) const
G4int getNumberOfNeutrons() const
G4CascadParticle initializeCascad(G4InuclElementaryParticle *particle)
G4bool worthToPropagate(const G4CascadParticle &cparticle) const
G4int getNumberOfProtons() const
void generateParticleFate(G4CascadParticle &cparticle, G4ElementaryParticleCollider *theEPCollider, std::vector< G4CascadParticle > &cascade)
G4bool empty() const
void reset(G4int nHitNeutrons=0, G4int nHitProtons=0, const std::vector< G4ThreeVector > *hitPoints=0)
G4bool stillInside(const G4CascadParticle &cparticle)
void generateModel(G4InuclNuclei *nuclei)
std::pair< std::vector< G4CascadParticle >, std::vector< G4InuclElementaryParticle > > modelLists
std::pair< G4int, G4int > getTypesOfNucleonsInvolved() const
G4double getRadiusUnits() const
void setVerboseLevel(G4int verbose)
const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:140
G4bool AreYouHit() const
Definition: G4Nucleon.hh:98
const G4ParticleDefinition * GetParticleType() const
Definition: G4Nucleon.hh:85
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4DecayTable * GetDecayTable() const
const G4String & GetParticleName() const
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4bool StartLoop()=0
virtual G4int GetMassNumber()=0
virtual void setVerboseLevel(G4int verbose=0)
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
static double Q[]