Geant4  10.01.p01
G4ElementaryParticleCollider.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4ElementaryParticleCollider.cc 71940 2013-06-28 19:04:44Z mkelsey $
27 //
28 // 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29 // 20100316 D. Wright (restored by M. Kelsey) -- Replace original (incorrect)
30 // pp, pn, nn 2-body to 2-body scattering angular distributions
31 // with a new parametrization of elastic scattering data using
32 // the sum of two exponentials.
33 // 20100319 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff;
34 // eliminate some unnecessary std::pow()
35 // 20100407 M. Kelsey -- Replace std::vector<>::resize(0) with ::clear()
36 // Eliminate return-by-value std::vector<> by creating buffers
37 // buffers for particles, momenta, and particle types.
38 // The following functions now return void and are non-const:
39 // ::generateSCMfinalState()
40 // ::generateMomModules()
41 // ::generateStrangeChannelPartTypes()
42 // ::generateSCMpionAbsorption()
43 // 20100408 M. Kelsey -- Follow changes to G4*Sampler to pass particle_kinds
44 // as input buffer.
45 // 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
46 // 20100428 M. Kelsey -- Use G4InuclParticleNames enum
47 // 20100429 M. Kelsey -- Change "photon()" to "isPhoton()"
48 // 20100507 M. Kelsey -- Rationalize multiplicity returns to be actual value
49 // 20100511 M. Kelsey -- Replace G4PionSampler and G4NucleonSampler with new
50 // pi-N and N-N classes, reorganize if-cascades
51 // 20100511 M. Kelsey -- Eliminate three residual random-angles blocks.
52 // 20100511 M. Kelsey -- Bug fix: pi-N two-body final states not correctly
53 // tested for charge-exchange case.
54 // 20100512 M. Kelsey -- Rationalize multiplicity returns to be actual value
55 // 20100512 M. Kelsey -- Add some additional debugging messages for 2-to-2
56 // 20100512 M. Kelsey -- Replace "if (is==)" cascades with switch blocks.
57 // Use G4CascadeInterpolator for angular distributions.
58 // 20100517 M. Kelsey -- Inherit from common base class, make arrays static
59 // 20100519 M. Kelsey -- Use G4InteractionCase to compute "is" values.
60 // 20100625 M. Kelsey -- Two bugs in n-body momentum, last particle recoil
61 // 20100713 M. Kelsey -- Bump collide start message up to verbose > 1
62 // 20100714 M. Kelsey -- Move conservation checking to base class
63 // 20100714 M. Kelsey -- Add sanity check for two-body final state, to ensure
64 // that final state total mass is below etot_scm; also compute
65 // kinematics without "rescaling" (which led to non-conservation)
66 // 20100726 M. Kelsey -- Move remaining std::vector<> buffers to .hh file
67 // 20100804 M. Kelsey -- Add printing of final-state tables, protected by
68 // G4CASCADE_DEBUG_SAMPLER preprocessor flag
69 // 20101019 M. Kelsey -- CoVerity report: check dynamic_cast<> for null
70 // 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
71 // 20110718 V. Uzhinskiy -- Drop IL,IM reset for multiplicity-three collisions
72 // 20110718 M. Kelsey -- Use enum names in switch blocks (c.f. G4NucleiModel)
73 // 20110720 M. Kelsey -- Follow interface change for cross-section tables,
74 // eliminating switch blocks.
75 // 20110806 M. Kelsey -- Pre-allocate buffers to avoid memory churn
76 // 20110922 M. Kelsey -- Follow G4InuclParticle::print(ostream&) migration
77 // 20110926 M. Kelsey -- Protect sampleCMcosFor2to2 from unphysical arguments
78 // 20111003 M. Kelsey -- Prepare for gamma-N interactions by checking for
79 // final-state tables instead of particle "isPhoton()"
80 // 20111007 M. Kelsey -- Add gamma-N final-state tables to printFinalState
81 // 20111107 M. Kelsey -- In sampleCMmomentumFor2to2(), hide message about
82 // unrecognized gamma-N initial state behind verbosity.
83 // 20111216 M. Kelsey -- Add diagnostics to generateMomModulesFor2toMany(),
84 // defer allocation of buffer in generateSCMfinalState() so that
85 // multiplicity failures return zero output, and can be trapped.
86 // 20120308 M. Kelsey -- Allow photons to interact with dibaryons (see
87 // changes in G4NucleiModel).
88 // 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
89 // 20121206 M. Kelsey -- Add Omega to printFinalStateTables(), remove line
90 // about "preliminary" gamma-N.
91 // 20130129 M. Kelsey -- Add static arrays and interpolators for two-body
92 // angular distributions (addresses MT thread-local issue)
93 // 20130221 M. Kelsey -- Move two-body angular dist classes to factory
94 // 20130306 M. Kelsey -- Use particle-name enums in if-blocks; add comments
95 // to sections of momentum-coefficient matrix; move final state
96 // table printing to G4CascadeChannelTables.
97 // 20130307 M. Kelsey -- Reverse order of dimensions for rmn array
98 // 20130307 M. Kelsey -- Use new momentum generator factory instead of rmn
99 // 20130308 M. Kelsey -- Move 3-body angle calc to G4InuclSpecialFunctions.
100 // 20130422 M. Kelsey -- Move kinematics to G4CascadeFinalStateAlgorithm;
101 // reduce nesting and replicated code in collide().
102 // 20130508 D. Wright -- Implement muon capture
103 // 20130511 M. Kelsey -- Check for neutrinos and skip them in ::collide()
104 // 20141009 M. Kelsey -- Add pion absorption by single nucleons, with
105 // nuclear recoil. Improves pi- capture performance.
106 // 20141030 M. Kelsey -- Check flag for whether to do pi-N absorption
107 // 20141201 M. Kelsey -- Check for null vector return from G4GDecay3
108 
110 #include "G4CascadeChannel.hh"
111 #include "G4CascadeChannelTables.hh"
112 #include "G4CascadeParameters.hh"
113 #include "G4CollisionOutput.hh"
114 #include "G4GDecay3.hh"
115 #include "G4InuclParticleNames.hh"
117 #include "G4LorentzConvertor.hh"
119 #include "G4NucleiModel.hh"
120 #include "G4ParticleLargerEkin.hh"
121 #include "G4TwoBodyAngularDist.hh"
122 #include "G4VMultiBodyMomDst.hh"
123 #include "G4VTwoBodyAngDst.hh"
124 #include "Randomize.hh"
125 #include <algorithm>
126 #include <cfloat>
127 #include <vector>
128 
129 using namespace G4InuclParticleNames;
130 using namespace G4InuclSpecialFunctions;
131 
132 typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
133 
134 
136  : G4CascadeColliderBase("G4ElementaryParticleCollider"),
137  nucleusA(0), nucleusZ(0) {;}
138 
139 
140 void
142  G4InuclParticle* target,
143  G4CollisionOutput& output)
144 {
145  if (verboseLevel > 1)
146  G4cout << " >>> G4ElementaryParticleCollider::collide" << G4endl;
147 
148  if (!useEPCollider(bullet,target)) { // Sanity check
149  G4cerr << " ElementaryParticleCollider -> can collide only particle with particle "
150  << G4endl;
151  return;
152  }
153 
154 #ifdef G4CASCADE_DEBUG_SAMPLER
155  static G4bool doPrintTables = true; // Once and only once per job
156  if (doPrintTables) {
158  doPrintTables = false;
159  }
160 #endif
161 
162  interCase.set(bullet, target); // To identify kind of collision
163 
164  if (verboseLevel > 1) G4cout << *bullet << G4endl << *target << G4endl;
165 
166  G4InuclElementaryParticle* particle1 =
167  dynamic_cast<G4InuclElementaryParticle*>(bullet);
168  G4InuclElementaryParticle* particle2 =
169  dynamic_cast<G4InuclElementaryParticle*>(target);
170 
171  if (!particle1 || !particle2) { // Redundant with useEPCollider()
172  G4cerr << " ElementaryParticleCollider -> can only collide hadrons"
173  << G4endl;
174  return;
175  }
176 
177  if (particle1->isNeutrino() || particle2->isNeutrino()) return;
178 
179  // Check for available interaction, or pion+dibaryon special case
181  !particle1->quasi_deutron() && !particle2->quasi_deutron()) {
182  G4cerr << " ElementaryParticleCollider -> cannot collide "
183  << particle1->getDefinition()->GetParticleName() << " with "
184  << particle2->getDefinition()->GetParticleName() << G4endl;
185  return;
186  }
187 
188  G4LorentzConvertor convertToSCM; // Utility to handle frame manipulation
189  if (particle2->nucleon() || particle2->quasi_deutron()) {
190  convertToSCM.setBullet(particle1);
191  convertToSCM.setTarget(particle2);
192  } else {
193  convertToSCM.setBullet(particle2);
194  convertToSCM.setTarget(particle1);
195  }
196 
197  convertToSCM.setVerbose(verboseLevel);
198  convertToSCM.toTheCenterOfMass();
199 
200  // Generate any particle collision with nucleon
201  if (particle1->nucleon() || particle2->nucleon()) {
202  G4double ekin = convertToSCM.getKinEnergyInTheTRS();
203  G4double etot_scm = convertToSCM.getTotalSCMEnergy();
204 
205  // SPECIAL: Very low energy pions may be absorbed by a nucleon
206  if (pionNucleonAbsorption(ekin)) {
207  generateSCMpionNAbsorption(etot_scm, particle1, particle2);
208  } else {
209  generateSCMfinalState(ekin, etot_scm, particle1, particle2);
210  }
211  }
212 
213  // Generate pion or photon collision with quasi-deuteron
214  if (particle1->quasi_deutron() || particle2->quasi_deutron()) {
215  if (!G4NucleiModel::useQuasiDeuteron(particle1->type(),particle2->type()) &&
216  !G4NucleiModel::useQuasiDeuteron(particle2->type(),particle1->type())) {
217  G4cerr << " ElementaryParticleCollider -> can only collide pi,mu,gamma with"
218  << " dibaryons " << G4endl;
219  return;
220  }
221 
222  G4double etot_scm = convertToSCM.getTotalSCMEnergy();
223 
224  if (particle1->isMuon() || particle2->isMuon()) {
225  generateSCMmuonAbsorption(etot_scm, particle1, particle2);
226  } else { // Currently, pion absoprtion also handles gammas
227  generateSCMpionAbsorption(etot_scm, particle1, particle2);
228  }
229  }
230 
231  if (particles.empty()) { // No final state possible, pass bullet through
232  if (verboseLevel) {
233  G4cerr << " ElementaryParticleCollider -> failed to collide "
234  << particle1->getMomModule() << " GeV/c "
235  << particle1->getDefinition()->GetParticleName() << " with "
236  << particle2->getDefinition()->GetParticleName() << G4endl;
237  }
238  return;
239  }
240 
241  // Convert final state back to lab frame
242  G4LorentzVector mom; // Buffer to avoid memory churn
243  particleIterator ipart;
244  for(ipart = particles.begin(); ipart != particles.end(); ipart++) {
245  mom = convertToSCM.backToTheLab(ipart->getMomentum());
246  ipart->setMomentum(mom);
247  };
248 
249  // Check conservation in multibody final state
250  if (verboseLevel && !validateOutput(bullet, target, particles)) {
251  G4cout << " incoming particles: \n" << *particle1 << G4endl
252  << *particle2 << G4endl
253  << " outgoing particles: " << G4endl;
254  for(ipart = particles.begin(); ipart != particles.end(); ipart++)
255  G4cout << *ipart << G4endl;
256 
257  G4cout << " <<< Non-conservation in G4ElementaryParticleCollider"
258  << G4endl;
259  }
260 
261  std::sort(particles.begin(), particles.end(), G4ParticleLargerEkin());
263 }
264 
265 
266 G4int
268  G4double ekin) const
269 {
270  G4int mul = 0;
271 
272  const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(is);
273 
274  if (xsecTable) mul = xsecTable->getMultiplicity(ekin);
275  else {
276  G4cerr << " G4ElementaryParticleCollider: Unknown interaction channel "
277  << is << " - multiplicity not generated " << G4endl;
278  }
279 
280  if(verboseLevel > 3){
281  G4cout << " G4ElementaryParticleCollider::generateMultiplicity: "
282  << " multiplicity = " << mul << G4endl;
283  }
284 
285  return mul;
286 }
287 
288 
289 void
291  G4double ekin)
292 {
293  particle_kinds.clear(); // Initialize buffer for generation
294 
295  const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(is);
296 
297  if (xsecTable)
298  xsecTable->getOutgoingParticleTypes(particle_kinds, mult, ekin);
299  else {
300  G4cerr << " G4ElementaryParticleCollider: Unknown interaction channel "
301  << is << " - outgoing kinds not generated " << G4endl;
302  }
303 
304  return;
305 }
306 
307 
308 void
310  G4double etot_scm,
311  G4InuclElementaryParticle* particle1,
312  G4InuclElementaryParticle* particle2) {
313  if (verboseLevel > 2) {
314  G4cout << " >>> G4ElementaryParticleCollider::generateSCMfinalState"
315  << G4endl;
316  }
317 
319 
320  const G4int itry_max = 10;
321 
322  G4int type1 = particle1->type();
323  G4int type2 = particle2->type();
324 
325  G4int is = type1 * type2;
326 
327  if (verboseLevel > 3) G4cout << " is " << is << G4endl;
328 
329  G4int multiplicity = 0;
330  G4bool generate = true;
331 
332  G4int itry = 0;
333  while (generate && itry++ < itry_max) {
334  particles.clear(); // Initialize buffers for this event
335  particle_kinds.clear();
336 
337  // Generate list of final-state particles
338  multiplicity = generateMultiplicity(is, ekin);
339 
340  generateOutgoingPartTypes(is, multiplicity, ekin);
341  if (particle_kinds.empty()) {
342  if (verboseLevel > 3) {
343  G4cout << " generateOutgoingPartTypes failed mult " << multiplicity
344  << G4endl;
345  }
346  continue;
347  }
348 
349  fillOutgoingMasses(); // Fill mass buffer from particle types
350 
351  // Attempt to produce final state kinematics
352  fsGenerator.Configure(particle1, particle2, particle_kinds);
353  generate = !fsGenerator.Generate(etot_scm, masses, scm_momentums);
354  } // while (generate)
355 
356  if (itry >= itry_max) { // Unable to generate valid final state
357  if (verboseLevel > 2)
358  G4cout << " generateSCMfinalState failed " << itry << " attempts"
359  << G4endl;
360  return;
361  }
362 
363  // Store generated momenta into outgoing particles
364 
365  particles.resize(multiplicity); // Preallocate buffer
366  for (G4int i=0; i<multiplicity; i++) {
367  particles[i].fill(scm_momentums[i], particle_kinds[i],
369  }
370 
371  if (verboseLevel > 3) {
372  G4cout << " <<< G4ElementaryParticleCollider::generateSCMfinalState"
373  << G4endl;
374  }
375 
376  return; // Particles buffer filled
377 }
378 
379 
380 // Use generated list of final states to fill mass buffers
381 
383  G4int mult = particle_kinds.size();
384 
385  masses.resize(mult,0.);
386  masses2.resize(mult,0.); // Allows direct [i] setting
387 
388  for (G4int i = 0; i < mult; i++) {
390  masses2[i] = masses[i] * masses[i];
391  }
392 }
393 
394 
395 // generate nucleons momenta for pion or photon absorption by dibaryon
396 // the nucleon distribution assumed to be isotropic in SCM
397 
398 void
400  G4InuclElementaryParticle* particle1,
401  G4InuclElementaryParticle* particle2) {
402  if (verboseLevel > 3)
403  G4cout << " >>> G4ElementaryParticleCollider::generateSCMpionAbsorption"
404  << G4endl;
405 
406  particles.clear(); // Initialize buffers for this event
407  particles.resize(2);
408 
409  particle_kinds.clear();
410 
411  G4int type1 = particle1->type();
412  G4int type2 = particle2->type();
413 
414  // Ensure that absportion is valid (charge conservable)
415  if (!G4NucleiModel::useQuasiDeuteron(type1, type2)) {
416  G4cerr << " pion absorption: "
417  << particle1->getDefinition()->GetParticleName() << " + "
418  << particle2->getDefinition()->GetParticleName() << " -> ?"
419  << G4endl;
420  return;
421  }
422 
423  if (!splitQuasiDeuteron(type2)) return; // Get constituents of [NN]
424 
426 
427  G4double a = 0.5 * (etot_scm * etot_scm - masses2[0] - masses2[1]);
428 
429  G4double pmod = std::sqrt((a*a - masses2[0]*masses2[1])
430  / (masses2[0] + masses2[1] + 2.0*a));
432  G4LorentzVector mom2;
433  mom2.setVectM(-mom1.vect(), masses[1]);
434 
437 }
438 
439 
440 // generate nucleons momenta for muon absorption by dibaryon
441 // the nucleon distribution assumed to be isotropic in SCM
442 
443 void
445  G4InuclElementaryParticle* particle1,
446  G4InuclElementaryParticle* particle2)
447 {
448  if (verboseLevel > 3)
449  G4cout << " >>> G4ElementaryParticleCollider::generateSCMmuonAbsorption"
450  << G4endl;
451 
452  // A phase space generator is required for the 3-body final state
453 
454  particles.clear(); // Initialize buffers for this event
455  particles.resize(3);
456 
457  scm_momentums.clear();
458  scm_momentums.resize(3);
459 
460  particle_kinds.clear();
461 
462  G4int type1 = particle1->type();
463  G4int type2 = particle2->type();
464 
465  if (type1 != muonMinus) return; // Sanity check, only mu- absorption
466 
467  // Ensure that absportion is valid (charge conservable)
468  if (!G4NucleiModel::useQuasiDeuteron(type1, type2)) {
469  G4cerr << " mu- absorption: "
470  << particle1->getDefinition()->GetParticleName() << " + "
471  << particle2->getDefinition()->GetParticleName() << " -> ?"
472  << G4endl;
473  return;
474  }
475 
476  if (!splitQuasiDeuteron(type2)) return; // Get constituents of [NN]
477  particle_kinds.push_back(mnu);
478 
480 
481  G4GDecay3 breakup(etot_scm, masses[0], masses[1], masses[2]);
482  std::vector<G4ThreeVector> theMomenta = breakup.GetThreeBodyMomenta();
483 
484  if (theMomenta.empty()) {
485  G4cerr << " generateSCMmuonAbsorption: GetThreeBodyMomenta() failed"
486  << " for " << type2 << " dibaryon" << G4endl;
487  particle_kinds.clear();
488  masses.clear();
489  return;
490  }
491 
492  for (size_t i=0; i<3; i++) {
493  scm_momentums[i].setVectM(theMomenta[i], masses[i]);
495  }
496 }
497 
498 
499 // generate nucleons momenta for pion absorption by single nucleon
500 
501 void
503  G4InuclElementaryParticle* particle1,
504  G4InuclElementaryParticle* particle2) {
505  if (verboseLevel > 3)
506  G4cout << " >>> G4ElementaryParticleCollider::generateSCMpionNAbsorption"
507  << G4endl;
508 
509  particles.clear(); // Initialize buffers for this event
510  particles.resize(1);
511 
512  particle_kinds.clear();
513 
514  G4int type1 = particle1->type();
515  G4int type2 = particle2->type();
516 
517  // Ensure that single-nucleon absportion is valid (charge exchangeable)
518  if ((type1*type2 != pim*pro && type1*type2 != pip*neu)) {
519  G4cerr << " pion-nucleon absorption: "
520  << particle1->getDefinition()->GetParticleName() << " + "
521  << particle2->getDefinition()->GetParticleName() << " -> ?"
522  << G4endl;
523  return;
524  }
525 
526  // Get outgoing nucleon type using charge exchange
527  // Proton code is 1, neutron code is 2, so 3-# swaps them
528  G4int ntype = (particle2->nucleon() ? type2 : type1);
529  G4int outType = 3 - ntype;
530  particle_kinds.push_back(outType);
531 
533 
534  // Get mass of residual nucleus (2-ntype = 1 for proton, 0 for neutron)
535  G4double mRecoil =
537  G4double mRecoil2 = mRecoil*mRecoil;
538 
539  // Recompute Ecm to include nucleus (for recoil kinematics)
540  G4LorentzVector vsum(0.,0.,0.,mRecoil);
541  vsum += particle1->getMomentum() + particle2->getMomentum();
542  etot_scm = vsum.m();
543 
544  // Two-body kinematics (nucleon against nucleus)
545  G4double a = 0.5 * (etot_scm*etot_scm - masses2[0] - mRecoil2);
546 
547  G4double pmod = std::sqrt((a*a - masses2[0]*mRecoil2)
548  / (masses2[0] + mRecoil2 + 2.0*a));
550 
551  if (verboseLevel > 3) {
552  G4cout << " outgoing type " << outType << " recoiling on nuclear mass "
553  << mRecoil << " a " << a << " p(SCM) " << pmod << " Ekin "
554  << mom1.e()-masses[0] << G4endl;
555  }
556 
557  // Nuclear recoil four momentum, for boosting back
558  G4LorentzVector mom2;
559  mom2.setVectM(-mom1.vect(), mRecoil);
560  mom1.boost(-mom2.boostVector());
561 
562  if (verboseLevel > 3) {
563  G4cout << " after nuclear recoil p " << mom1.rho() << " Ekin "
564  << mom1.e()-masses[0] << G4endl;
565  }
566 
567  // Fill only the ejected nucleon
569 }
570 
571 
572 // Evaluate whether interaction is candidate for absorption on nucleon
573 
574 G4bool
576  if (!G4CascadeParameters::piNAbsorption()) return false;
577 
578  if (verboseLevel > 3)
579  G4cout << " >>> G4ElementaryParticleCollider::pionNucleonAbsorption ?"
580  << " ekin " << ekin << " is " << interCase.hadrons() << G4endl;
581 
583 
584  // Absorption occurs only for pi- p -> n, or pi+ n -> p
585  // Restrict to "very slow" pions, to allow for some normal scattering
586  G4bool isAbsorbable = ((interCase.hadrons() == pim*pro ||
587  interCase.hadrons() == pip*neu)
588  && ekin < 0.5*mpi0);
589 
590  // FIXME: Should have an effective cross-section to throw random here
591  return isAbsorbable;
592 }
593 
594 
595 // generate constituents of dibaryon for "explosion"
596 
598  if (qdtype != diproton && qdtype != unboundPN && qdtype != dineutron) {
599  G4cerr << " type " << qdtype << " not dibaryon!" << G4endl;
600  return false;
601  }
602 
603  G4int b2 = qdtype % 10; // Dibaryon codes are 1ab (a=1,2; b=1,2)
604  G4int b1 = (qdtype/10) % 10;
605 
606  particle_kinds.push_back(b1);
607  particle_kinds.push_back(b2);
608 
609  return true;
610 }
void generateSCMpionNAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
void generateSCMpionAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
std::vector< G4ThreeVector > GetThreeBodyMomenta()
Definition: G4GDecay3.cc:99
static const G4CascadeChannel * GetTable(G4int initialState)
G4bool Generate(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
G4CascadeFinalStateGenerator fsGenerator
G4LorentzVector getMomentum() const
virtual G4bool useEPCollider(G4InuclParticle *bullet, G4InuclParticle *target) const
virtual void getOutgoingParticleTypes(std::vector< G4int > &kinds, G4int mult, G4double ke) const =0
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
const G4ParticleDefinition * getDefinition() const
G4int hadrons() const
G4double a
Definition: TRTMaterials.hh:39
static G4bool useQuasiDeuteron(G4int ptype, G4int qdtype=0)
static G4double getParticleMass(G4int type)
void setBullet(const G4InuclParticle *bullet)
G4double getTotalSCMEnergy() const
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
int G4int
Definition: G4Types.hh:78
void generateSCMmuonAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
const G4String & GetParticleName() const
std::vector< G4LorentzVector > scm_momentums
void generateSCMfinalState(G4double ekin, G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
G4GLOB_DLL std::ostream G4cout
void setVerbose(G4int vb=0)
std::vector< G4InuclElementaryParticle > particles
bool G4bool
Definition: G4Types.hh:79
G4double getNucleiMass() const
void Configure(G4InuclElementaryParticle *bullet, G4InuclElementaryParticle *target, const std::vector< G4int > &particle_kinds)
void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin)
G4LorentzVector generateWithRandomAngles(G4double p, G4double mass=0.)
static const G4double b2
G4double getKinEnergyInTheTRS() const
G4int generateMultiplicity(G4int is, G4double ekin) const
static void Print(std::ostream &os=G4cout)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
void set(G4InuclParticle *part1, G4InuclParticle *part2)
void generate(const G4double sqrtS, ParticleList &particles)
Generate an event in the CM system.
void SetVerboseLevel(G4int verbose)
virtual G4bool validateOutput(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
static const G4double b1
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4bool piNAbsorption()
virtual G4int getMultiplicity(G4double ke) const =0
std::vector< G4InuclElementaryParticle >::iterator particleIterator
G4bool pionNucleonAbsorption(G4double ekin) const
G4double getMomModule() const
G4GLOB_DLL std::ostream G4cerr
void setTarget(const G4InuclParticle *target)
CLHEP::HepLorentzVector G4LorentzVector