ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ElementaryParticleCollider.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4ElementaryParticleCollider.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 //
27 // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
28 // 20100316 D. Wright (restored by M. Kelsey) -- Replace original (incorrect)
29 // pp, pn, nn 2-body to 2-body scattering angular distributions
30 // with a new parametrization of elastic scattering data using
31 // the sum of two exponentials.
32 // 20100319 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff;
33 // eliminate some unnecessary std::pow()
34 // 20100407 M. Kelsey -- Replace std::vector<>::resize(0) with ::clear()
35 // Eliminate return-by-value std::vector<> by creating buffers
36 // buffers for particles, momenta, and particle types.
37 // The following functions now return void and are non-const:
38 // ::generateSCMfinalState()
39 // ::generateMomModules()
40 // ::generateStrangeChannelPartTypes()
41 // ::generateSCMpionAbsorption()
42 // 20100408 M. Kelsey -- Follow changes to G4*Sampler to pass particle_kinds
43 // as input buffer.
44 // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
45 // 20100428 M. Kelsey -- Use G4InuclParticleNames enum
46 // 20100429 M. Kelsey -- Change "photon()" to "isPhoton()"
47 // 20100507 M. Kelsey -- Rationalize multiplicity returns to be actual value
48 // 20100511 M. Kelsey -- Replace G4PionSampler and G4NucleonSampler with new
49 // pi-N and N-N classes, reorganize if-cascades
50 // 20100511 M. Kelsey -- Eliminate three residual random-angles blocks.
51 // 20100511 M. Kelsey -- Bug fix: pi-N two-body final states not correctly
52 // tested for charge-exchange case.
53 // 20100512 M. Kelsey -- Rationalize multiplicity returns to be actual value
54 // 20100512 M. Kelsey -- Add some additional debugging messages for 2-to-2
55 // 20100512 M. Kelsey -- Replace "if (is==)" cascades with switch blocks.
56 // Use G4CascadeInterpolator for angular distributions.
57 // 20100517 M. Kelsey -- Inherit from common base class, make arrays static
58 // 20100519 M. Kelsey -- Use G4InteractionCase to compute "is" values.
59 // 20100625 M. Kelsey -- Two bugs in n-body momentum, last particle recoil
60 // 20100713 M. Kelsey -- Bump collide start message up to verbose > 1
61 // 20100714 M. Kelsey -- Move conservation checking to base class
62 // 20100714 M. Kelsey -- Add sanity check for two-body final state, to ensure
63 // that final state total mass is below etot_scm; also compute
64 // kinematics without "rescaling" (which led to non-conservation)
65 // 20100726 M. Kelsey -- Move remaining std::vector<> buffers to .hh file
66 // 20100804 M. Kelsey -- Add printing of final-state tables, protected by
67 // G4CASCADE_DEBUG_SAMPLER preprocessor flag
68 // 20101019 M. Kelsey -- CoVerity report: check dynamic_cast<> for null
69 // 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
70 // 20110718 V. Uzhinskiy -- Drop IL,IM reset for multiplicity-three collisions
71 // 20110718 M. Kelsey -- Use enum names in switch blocks (c.f. G4NucleiModel)
72 // 20110720 M. Kelsey -- Follow interface change for cross-section tables,
73 // eliminating switch blocks.
74 // 20110806 M. Kelsey -- Pre-allocate buffers to avoid memory churn
75 // 20110922 M. Kelsey -- Follow G4InuclParticle::print(ostream&) migration
76 // 20110926 M. Kelsey -- Protect sampleCMcosFor2to2 from unphysical arguments
77 // 20111003 M. Kelsey -- Prepare for gamma-N interactions by checking for
78 // final-state tables instead of particle "isPhoton()"
79 // 20111007 M. Kelsey -- Add gamma-N final-state tables to printFinalState
80 // 20111107 M. Kelsey -- In sampleCMmomentumFor2to2(), hide message about
81 // unrecognized gamma-N initial state behind verbosity.
82 // 20111216 M. Kelsey -- Add diagnostics to generateMomModulesFor2toMany(),
83 // defer allocation of buffer in generateSCMfinalState() so that
84 // multiplicity failures return zero output, and can be trapped.
85 // 20120308 M. Kelsey -- Allow photons to interact with dibaryons (see
86 // changes in G4NucleiModel).
87 // 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
88 // 20121206 M. Kelsey -- Add Omega to printFinalStateTables(), remove line
89 // about "preliminary" gamma-N.
90 // 20130129 M. Kelsey -- Add static arrays and interpolators for two-body
91 // angular distributions (addresses MT thread-local issue)
92 // 20130221 M. Kelsey -- Move two-body angular dist classes to factory
93 // 20130306 M. Kelsey -- Use particle-name enums in if-blocks; add comments
94 // to sections of momentum-coefficient matrix; move final state
95 // table printing to G4CascadeChannelTables.
96 // 20130307 M. Kelsey -- Reverse order of dimensions for rmn array
97 // 20130307 M. Kelsey -- Use new momentum generator factory instead of rmn
98 // 20130308 M. Kelsey -- Move 3-body angle calc to G4InuclSpecialFunctions.
99 // 20130422 M. Kelsey -- Move kinematics to G4CascadeFinalStateAlgorithm;
100 // reduce nesting and replicated code in collide().
101 // 20130508 D. Wright -- Implement muon capture
102 // 20130511 M. Kelsey -- Check for neutrinos and skip them in ::collide()
103 // 20141009 M. Kelsey -- Add pion absorption by single nucleons, with
104 // nuclear recoil. Improves pi- capture performance.
105 // 20141030 M. Kelsey -- Check flag for whether to do pi-N absorption
106 // 20141201 M. Kelsey -- Check for null vector return from G4GDecay3
107 // 20141211 M. Kelsey -- Change PIN_ABSORPTION flag to double, probability;
108 // fix handling of boosts for pi-N absorption.
109 // 20150608 M. Kelsey -- Label all while loops as terminating.
110 
112 #include "G4CascadeChannel.hh"
113 #include "G4CascadeChannelTables.hh"
114 #include "G4CascadeParameters.hh"
115 #include "G4CollisionOutput.hh"
116 #include "G4GDecay3.hh"
117 #include "G4InuclParticleNames.hh"
119 #include "G4LorentzConvertor.hh"
121 #include "G4NucleiModel.hh"
122 #include "G4ParticleLargerEkin.hh"
123 #include "G4TwoBodyAngularDist.hh"
124 #include "G4VMultiBodyMomDst.hh"
125 #include "G4VTwoBodyAngDst.hh"
126 #include "Randomize.hh"
127 #include <algorithm>
128 #include <cfloat>
129 #include <vector>
130 
131 using namespace G4InuclParticleNames;
132 using namespace G4InuclSpecialFunctions;
133 
134 typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
135 
136 
138  : G4CascadeColliderBase("G4ElementaryParticleCollider"),
139  nucleusA(0), nucleusZ(0) {;}
140 
141 
142 void
145  G4CollisionOutput& output)
146 {
147  if (verboseLevel > 1)
148  G4cout << " >>> G4ElementaryParticleCollider::collide" << G4endl;
149 
150  if (!useEPCollider(bullet,target)) { // Sanity check
151  G4cerr << " ElementaryParticleCollider -> can collide only particle with particle "
152  << G4endl;
153  return;
154  }
155 
156 #ifdef G4CASCADE_DEBUG_SAMPLER
157  static G4bool doPrintTables = true; // Once and only once per job
158  if (doPrintTables) {
160  doPrintTables = false;
161  }
162 #endif
163 
164  interCase.set(bullet, target); // To identify kind of collision
165 
166  if (verboseLevel > 1) G4cout << *bullet << G4endl << *target << G4endl;
167 
168  G4InuclElementaryParticle* particle1 =
169  dynamic_cast<G4InuclElementaryParticle*>(bullet);
170  G4InuclElementaryParticle* particle2 =
171  dynamic_cast<G4InuclElementaryParticle*>(target);
172 
173  if (!particle1 || !particle2) { // Redundant with useEPCollider()
174  G4cerr << " ElementaryParticleCollider -> can only collide hadrons"
175  << G4endl;
176  return;
177  }
178 
179  if (particle1->isNeutrino() || particle2->isNeutrino()) return;
180 
181  // Check for available interaction, or pion+dibaryon special case
183  !particle1->quasi_deutron() && !particle2->quasi_deutron()) {
184  G4cerr << " ElementaryParticleCollider -> cannot collide "
185  << particle1->getDefinition()->GetParticleName() << " with "
186  << particle2->getDefinition()->GetParticleName() << G4endl;
187  return;
188  }
189 
190  G4LorentzConvertor convertToSCM; // Utility to handle frame manipulation
191  if (particle2->nucleon() || particle2->quasi_deutron()) {
192  convertToSCM.setBullet(particle1);
193  convertToSCM.setTarget(particle2);
194  } else {
195  convertToSCM.setBullet(particle2);
196  convertToSCM.setTarget(particle1);
197  }
198 
199  convertToSCM.setVerbose(verboseLevel);
200  convertToSCM.toTheCenterOfMass();
201 
202  G4double etot_scm = convertToSCM.getTotalSCMEnergy();
203 
204  // Generate any particle collision with nucleon
205  if (particle1->nucleon() || particle2->nucleon()) {
206  G4double ekin = convertToSCM.getKinEnergyInTheTRS();
207 
208  // SPECIAL: Very low energy pions may be absorbed by a nucleon
209  if (pionNucleonAbsorption(ekin)) {
210  generateSCMpionNAbsorption(etot_scm, particle1, particle2);
211  } else {
212  generateSCMfinalState(ekin, etot_scm, particle1, particle2);
213  }
214  }
215 
216  // Generate pion or photon collision with quasi-deuteron
217  if (particle1->quasi_deutron() || particle2->quasi_deutron()) {
218  if (!G4NucleiModel::useQuasiDeuteron(particle1->type(),particle2->type()) &&
219  !G4NucleiModel::useQuasiDeuteron(particle2->type(),particle1->type())) {
220  G4cerr << " ElementaryParticleCollider -> can only collide pi,mu,gamma with"
221  << " dibaryons " << G4endl;
222  return;
223  }
224 
225  if (particle1->isMuon() || particle2->isMuon()) {
226  generateSCMmuonAbsorption(etot_scm, particle1, particle2);
227  } else { // Currently, pion absoprtion also handles gammas
228  generateSCMpionAbsorption(etot_scm, particle1, particle2);
229  }
230  }
231 
232  if (particles.empty()) { // No final state possible, pass bullet through
233  if (verboseLevel) {
234  G4cerr << " ElementaryParticleCollider -> failed to collide "
235  << particle1->getMomModule() << " GeV/c "
236  << particle1->getDefinition()->GetParticleName() << " with "
237  << particle2->getDefinition()->GetParticleName() << G4endl;
238  }
239  return;
240  }
241 
242  // Convert final state back to lab frame
243  G4LorentzVector mom; // Buffer to avoid memory churn
245  for(ipart = particles.begin(); ipart != particles.end(); ipart++) {
246  mom = convertToSCM.backToTheLab(ipart->getMomentum());
247  ipart->setMomentum(mom);
248  };
249 
250  if (verboseLevel) {
251  // Check conservation in multibody final state
252  const G4ParticleDefinition* bDef = bullet->getDefinition();
253  const G4ParticleDefinition* tDef = target->getDefinition();
254 
255  G4int initBaryonNumber = bDef->GetBaryonNumber() + tDef->GetBaryonNumber();
256  G4int initCharge = bullet->getCharge() + target->getCharge();
257  G4int initStrangeness = bDef->GetQuarkContent(3) - bDef->GetAntiQuarkContent(3) +
258  tDef->GetQuarkContent(3) - tDef->GetAntiQuarkContent(3);
259 
260  G4int finalBaryonNumber = 0;
261  G4int finalCharge = 0;
262  G4int finalStrangeness = 0;
263 
264  for (ipart = particles.begin(); ipart != particles.end(); ipart++) {
265  finalBaryonNumber += ipart->baryon();
266  finalCharge += ipart->getCharge();
267  finalStrangeness += ipart->getStrangeness();
268  }
269 
270  G4int bnc = finalBaryonNumber - initBaryonNumber;
271  G4int cnc = finalCharge - initCharge;
272  G4int snc = finalStrangeness - initStrangeness;
273 
274  if (bnc != 0 || cnc != 0 || snc != 0) {
275  G4cout << " G4ElementaryParticleCollider: quantum number non-conservation " << G4endl;
276  G4cout << " Baryon number: initial = " << initBaryonNumber << ", final = "
277  << finalBaryonNumber << G4endl;
278  G4cout << " Charge: initial = " << initCharge << ", final = "
279  << finalCharge << G4endl;
280  G4cout << " Strangeness: initial = " << initStrangeness << ", final = "
281  << finalStrangeness << G4endl;
282 
283  G4cout << " bullet = " << bDef->GetParticleName() << G4endl;
284  G4cout << " target = " << tDef->GetParticleName() << G4endl;
285  G4cout << " secondaries = " ;
286  for (ipart = particles.begin(); ipart != particles.end(); ipart++) {
287  G4cout << ipart->getDefinition()->GetParticleName() << " " ;
288  }
289  G4cout << G4endl;
290  }
291  }
292 
293  if (verboseLevel && !validateOutput(bullet, target, particles)) {
294  G4cout << " incoming particles: \n" << *particle1 << G4endl
295  << *particle2 << G4endl
296  << " outgoing particles: " << G4endl;
297  for(ipart = particles.begin(); ipart != particles.end(); ipart++)
298  G4cout << *ipart << G4endl;
299 
300  G4cout << " <<< Non-conservation in G4ElementaryParticleCollider"
301  << G4endl;
302  }
303 
304  std::sort(particles.begin(), particles.end(), G4ParticleLargerEkin());
306 }
307 
308 
309 G4int
311  G4double ekin) const
312 {
313  G4int mul = 0;
314 
315  const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(is);
316 
317  if (xsecTable) mul = xsecTable->getMultiplicity(ekin);
318  else {
319  G4cerr << " G4ElementaryParticleCollider: Unknown interaction channel "
320  << is << " - multiplicity not generated " << G4endl;
321  }
322 
323  if(verboseLevel > 3){
324  G4cout << " G4ElementaryParticleCollider::generateMultiplicity: "
325  << " multiplicity = " << mul << G4endl;
326  }
327 
328  return mul;
329 }
330 
331 
332 void
334  G4double ekin)
335 {
336  particle_kinds.clear(); // Initialize buffer for generation
337 
338  const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(is);
339 
340  if (xsecTable)
341  xsecTable->getOutgoingParticleTypes(particle_kinds, mult, ekin);
342  else {
343  G4cerr << " G4ElementaryParticleCollider: Unknown interaction channel "
344  << is << " - outgoing kinds not generated " << G4endl;
345  }
346 
347  return;
348 }
349 
350 
351 void
353  G4double etot_scm,
354  G4InuclElementaryParticle* particle1,
355  G4InuclElementaryParticle* particle2) {
356  if (verboseLevel > 2) {
357  G4cout << " >>> G4ElementaryParticleCollider::generateSCMfinalState"
358  << G4endl;
359  }
360 
362 
363  const G4int itry_max = 10;
364 
365  G4int type1 = particle1->type();
366  G4int type2 = particle2->type();
367 
368  G4int is = type1 * type2;
369 
370  if (verboseLevel > 3) G4cout << " is " << is << G4endl;
371 
372  G4int multiplicity = 0;
373  G4bool generate = true;
374 
375  G4int itry = 0;
376  while (generate && itry++ < itry_max) { /* Loop checking 08.06.2015 MHK */
377  particles.clear(); // Initialize buffers for this event
378  particle_kinds.clear();
379 
380  // Generate list of final-state particles
381  multiplicity = generateMultiplicity(is, ekin);
382 
383  generateOutgoingPartTypes(is, multiplicity, ekin);
384  if (particle_kinds.empty()) {
385  if (verboseLevel > 3) {
386  G4cout << " generateOutgoingPartTypes failed mult " << multiplicity
387  << G4endl;
388  }
389  continue;
390  }
391 
392  fillOutgoingMasses(); // Fill mass buffer from particle types
393 
394  // Attempt to produce final state kinematics
395  fsGenerator.Configure(particle1, particle2, particle_kinds);
396  generate = !fsGenerator.Generate(etot_scm, masses, scm_momentums);
397  } // while (generate)
398 
399  if (itry >= itry_max) { // Unable to generate valid final state
400  if (verboseLevel > 2)
401  G4cout << " generateSCMfinalState failed " << itry << " attempts"
402  << G4endl;
403  return;
404  }
405 
406  // Store generated momenta into outgoing particles
407 
408  particles.resize(multiplicity); // Preallocate buffer
409  for (G4int i=0; i<multiplicity; i++) {
410  particles[i].fill(scm_momentums[i], particle_kinds[i],
412  }
413 
414  if (verboseLevel > 3) {
415  G4cout << " <<< G4ElementaryParticleCollider::generateSCMfinalState"
416  << G4endl;
417  }
418 
419  return; // Particles buffer filled
420 }
421 
422 
423 // Use generated list of final states to fill mass buffers
424 
426  G4int mult = particle_kinds.size();
427 
428  masses.resize(mult,0.);
429  masses2.resize(mult,0.); // Allows direct [i] setting
430 
431  for (G4int i = 0; i < mult; i++) {
433  masses2[i] = masses[i] * masses[i];
434  }
435 }
436 
437 
438 // Generate nucleon momenta for pion or photon absorption by dibaryon
439 // The nucleon distribution is assumed to be isotropic in SCM
440 
441 void
443  G4InuclElementaryParticle* particle1,
444  G4InuclElementaryParticle* particle2)
445 {
446  if (verboseLevel > 3)
447  G4cout << " >>> G4ElementaryParticleCollider::generateSCMpionAbsorption"
448  << G4endl;
449 
450  particles.clear(); // Initialize buffers for this event
451  particles.resize(2);
452 
453  particle_kinds.clear();
454 
455  G4int typeProduct = particle1->type() * particle2->type();
456 
457  if (typeProduct == pi0*diproton || typeProduct == pip*unboundPN ||
458  typeProduct == gam*diproton) {
459  particle_kinds.push_back(pro);
460  particle_kinds.push_back(pro);
461 
462  } else if (typeProduct == pim*diproton || typeProduct == pip*dineutron ||
463  typeProduct == pi0*unboundPN || typeProduct == gam*unboundPN) {
464  particle_kinds.push_back(pro);
465  particle_kinds.push_back(neu);
466 
467  } else if (typeProduct == pi0*dineutron || typeProduct == pim*unboundPN ||
468  typeProduct == gam*dineutron) {
469  particle_kinds.push_back(neu);
470  particle_kinds.push_back(neu);
471 
472  } else {
473  G4cerr << " Illegal absorption: "
474  << particle1->getDefinition()->GetParticleName() << " + "
475  << particle2->getDefinition()->GetParticleName() << " -> ?"
476  << G4endl;
477  return;
478  }
479 
481 
482  G4double a = 0.5 * (etot_scm * etot_scm - masses2[0] - masses2[1]);
483 
484  G4double pmod = std::sqrt((a*a - masses2[0]*masses2[1])
485  / (etot_scm*etot_scm) );
487  G4LorentzVector mom2;
488  mom2.setVectM(-mom1.vect(), masses[1]);
489 
492 }
493 
494 
495 // Generate nucleon momenta for muon absorption by dibaryon
496 // The nucleon distribution is assumed to be isotropic in SCM
497 
498 void
500  G4InuclElementaryParticle* particle1,
501  G4InuclElementaryParticle* particle2)
502 {
503  if (verboseLevel > 3)
504  G4cout << " >>> G4ElementaryParticleCollider::generateSCMmuonAbsorption"
505  << G4endl;
506 
507  particles.clear(); // Initialize buffers for this event
508  particles.resize(3);
509 
510  scm_momentums.clear();
511  scm_momentums.resize(3);
512 
513  particle_kinds.clear();
514 
515  G4int typeProduct = particle1->type() * particle2->type();
516 
517  if (typeProduct == mum*diproton) {
518  particle_kinds.push_back(pro);
519  particle_kinds.push_back(neu);
520  } else if (typeProduct == mum*unboundPN) {
521  particle_kinds.push_back(neu);
522  particle_kinds.push_back(neu);
523  } else {
524  G4cerr << " Illegal absorption: "
525  << particle1->getDefinition()->GetParticleName() << " + "
526  << particle2->getDefinition()->GetParticleName() << " -> ?"
527  << G4endl;
528  return;
529  }
530  particle_kinds.push_back(mnu);
531 
533 
534  // Phase space generator required for the 3-body final state
535  G4GDecay3 breakup(etot_scm, masses[0], masses[1], masses[2]);
536  std::vector<G4ThreeVector> theMomenta = breakup.GetThreeBodyMomenta();
537 
538  if (theMomenta.empty()) {
539  G4cerr << " generateSCMmuonAbsorption: GetThreeBodyMomenta() failed"
540  << " for " << particle2->type() << " dibaryon" << G4endl;
541  particle_kinds.clear();
542  masses.clear();
543  particles.clear();
544  return;
545  }
546 
547  for (size_t i=0; i<3; i++) {
548  scm_momentums[i].setVectM(theMomenta[i], masses[i]);
550  }
551 }
552 
553 
554 // generate nucleons momenta for pion absorption by single nucleon
555 
556 void
558  G4InuclElementaryParticle* particle1,
559  G4InuclElementaryParticle* particle2) {
560  if (verboseLevel > 3)
561  G4cout << " >>> G4ElementaryParticleCollider::generateSCMpionNAbsorption"
562  << G4endl;
563 
564  particles.clear(); // Initialize buffers for this event
565  particles.resize(1);
566 
567  particle_kinds.clear();
568 
569  G4int type1 = particle1->type();
570  G4int type2 = particle2->type();
571 
572  // Ensure that single-nucleon absportion is valid (charge exchangeable)
573  if ((type1*type2 != pim*pro && type1*type2 != pip*neu)) {
574  G4cerr << " pion-nucleon absorption: "
575  << particle1->getDefinition()->GetParticleName() << " + "
576  << particle2->getDefinition()->GetParticleName() << " -> ?"
577  << G4endl;
578  return;
579  }
580 
581  // Get outgoing nucleon type using charge exchange
582  // Proton code is 1, neutron code is 2, so 3-# swaps them
583  G4int ntype = (particle2->nucleon() ? type2 : type1);
584  G4int outType = 3 - ntype;
585  particle_kinds.push_back(outType);
586 
588 
589  // Get mass of residual nucleus (2-ntype = 1 for proton, 0 for neutron)
590  G4double mRecoil =
592  G4double mRecoil2 = mRecoil*mRecoil;
593 
594  // Recompute Ecm to include nucleus (for recoil kinematics)
595  G4LorentzVector piN4 = particle1->getMomentum() + particle2->getMomentum();
596  G4LorentzVector vsum(0.,0.,0.,mRecoil);
597  vsum += piN4;
598 
599  // Two-body kinematics (nucleon against nucleus) in overall CM system
600  G4double esq_scm = vsum.m2();
601  G4double a = 0.5 * (esq_scm - masses2[0] - mRecoil2);
602 
603  G4double pmod = std::sqrt((a*a - masses2[0]*mRecoil2) / esq_scm );
605 
606  if (verboseLevel > 3) {
607  G4cout << " outgoing type " << outType << " recoiling on nuclear mass "
608  << mRecoil << "\n a " << a << " p " << pmod << " Ekin "
609  << mom1.e()-masses[0] << G4endl;
610  }
611 
612  mom1.boost(-piN4.boostVector()); // Boost into CM of pi-N collision
613 
614  if (verboseLevel > 3) {
615  G4cout << " in original pi-N frame p(SCM) " << mom1.rho() << " Ekin "
616  << mom1.e()-masses[0] << G4endl;
617  }
618 
619  // Fill only the ejected nucleon
621 }
622 
623 
624 // Evaluate whether interaction is candidate for absorption on nucleon
625 G4bool
627  if (verboseLevel > 3)
628  G4cout << " >>> G4ElementaryParticleCollider::pionNucleonAbsorption ?"
629  << " ekin " << ekin << " is " << interCase.hadrons() << G4endl;
630 
631  // Absorption occurs with specified probability
633 
634  // Absorption occurs only for pi- p -> n, or pi+ n -> p
635  // Restrict to "very slow" pions, to allow for some normal scattering
636  return ((interCase.hadrons() == pim*pro ||
637  interCase.hadrons() == pip*neu)
638  && (ekin < 0.05) // 50 MeV kinetic energy or less
639  && (G4UniformRand() < absProb)
640  );
641 }
642