Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
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 // $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 // 20141211 M. Kelsey -- Change PIN_ABSORPTION flag to double, probability;
109 // fix handling of boosts for pi-N absorption.
110 // 20150608 M. Kelsey -- Label all while loops as terminating.
111 
113 #include "G4CascadeChannel.hh"
114 #include "G4CascadeChannelTables.hh"
115 #include "G4CascadeParameters.hh"
116 #include "G4CollisionOutput.hh"
117 #include "G4GDecay3.hh"
118 #include "G4InuclParticleNames.hh"
120 #include "G4LorentzConvertor.hh"
122 #include "G4NucleiModel.hh"
123 #include "G4ParticleLargerEkin.hh"
124 #include "G4TwoBodyAngularDist.hh"
125 #include "G4VMultiBodyMomDst.hh"
126 #include "G4VTwoBodyAngDst.hh"
127 #include "Randomize.hh"
128 #include <algorithm>
129 #include <cfloat>
130 #include <vector>
131 
132 using namespace G4InuclParticleNames;
133 using namespace G4InuclSpecialFunctions;
134 
135 typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
136 
137 
139  : G4CascadeColliderBase("G4ElementaryParticleCollider"),
140  nucleusA(0), nucleusZ(0) {;}
141 
142 
143 void
146  G4CollisionOutput& output)
147 {
148  if (verboseLevel > 1)
149  G4cout << " >>> G4ElementaryParticleCollider::collide" << G4endl;
150 
151  if (!useEPCollider(bullet,target)) { // Sanity check
152  G4cerr << " ElementaryParticleCollider -> can collide only particle with particle "
153  << G4endl;
154  return;
155  }
156 
157 #ifdef G4CASCADE_DEBUG_SAMPLER
158  static G4bool doPrintTables = true; // Once and only once per job
159  if (doPrintTables) {
161  doPrintTables = false;
162  }
163 #endif
164 
165  interCase.set(bullet, target); // To identify kind of collision
166 
167  if (verboseLevel > 1) G4cout << *bullet << G4endl << *target << G4endl;
168 
169  G4InuclElementaryParticle* particle1 =
170  dynamic_cast<G4InuclElementaryParticle*>(bullet);
171  G4InuclElementaryParticle* particle2 =
172  dynamic_cast<G4InuclElementaryParticle*>(target);
173 
174  if (!particle1 || !particle2) { // Redundant with useEPCollider()
175  G4cerr << " ElementaryParticleCollider -> can only collide hadrons"
176  << G4endl;
177  return;
178  }
179 
180  if (particle1->isNeutrino() || particle2->isNeutrino()) return;
181 
182  // Check for available interaction, or pion+dibaryon special case
184  !particle1->quasi_deutron() && !particle2->quasi_deutron()) {
185  G4cerr << " ElementaryParticleCollider -> cannot collide "
186  << particle1->getDefinition()->GetParticleName() << " with "
187  << particle2->getDefinition()->GetParticleName() << G4endl;
188  return;
189  }
190 
191  G4LorentzConvertor convertToSCM; // Utility to handle frame manipulation
192  if (particle2->nucleon() || particle2->quasi_deutron()) {
193  convertToSCM.setBullet(particle1);
194  convertToSCM.setTarget(particle2);
195  } else {
196  convertToSCM.setBullet(particle2);
197  convertToSCM.setTarget(particle1);
198  }
199 
200  convertToSCM.setVerbose(verboseLevel);
201  convertToSCM.toTheCenterOfMass();
202 
203  G4double etot_scm = convertToSCM.getTotalSCMEnergy();
204 
205  // Generate any particle collision with nucleon
206  if (particle1->nucleon() || particle2->nucleon()) {
207  G4double ekin = convertToSCM.getKinEnergyInTheTRS();
208 
209  // SPECIAL: Very low energy pions may be absorbed by a nucleon
210  if (pionNucleonAbsorption(ekin)) {
211  generateSCMpionNAbsorption(etot_scm, particle1, particle2);
212  } else {
213  generateSCMfinalState(ekin, etot_scm, particle1, particle2);
214  }
215  }
216 
217  // Generate pion or photon collision with quasi-deuteron
218  if (particle1->quasi_deutron() || particle2->quasi_deutron()) {
219  if (!G4NucleiModel::useQuasiDeuteron(particle1->type(),particle2->type()) &&
220  !G4NucleiModel::useQuasiDeuteron(particle2->type(),particle1->type())) {
221  G4cerr << " ElementaryParticleCollider -> can only collide pi,mu,gamma with"
222  << " dibaryons " << G4endl;
223  return;
224  }
225 
226  if (particle1->isMuon() || particle2->isMuon()) {
227  generateSCMmuonAbsorption(etot_scm, particle1, particle2);
228  } else { // Currently, pion absoprtion also handles gammas
229  generateSCMpionAbsorption(etot_scm, particle1, particle2);
230  }
231  }
232 
233  if (particles.empty()) { // No final state possible, pass bullet through
234  if (verboseLevel) {
235  G4cerr << " ElementaryParticleCollider -> failed to collide "
236  << particle1->getMomModule() << " GeV/c "
237  << particle1->getDefinition()->GetParticleName() << " with "
238  << particle2->getDefinition()->GetParticleName() << G4endl;
239  }
240  return;
241  }
242 
243  // Convert final state back to lab frame
244  G4LorentzVector mom; // Buffer to avoid memory churn
246  for(ipart = particles.begin(); ipart != particles.end(); ipart++) {
247  mom = convertToSCM.backToTheLab(ipart->getMomentum());
248  ipart->setMomentum(mom);
249  };
250 
251  if (verboseLevel) {
252  // Check conservation in multibody final state
253  const G4ParticleDefinition* bDef = bullet->getDefinition();
254  const G4ParticleDefinition* tDef = target->getDefinition();
255 
256  G4int initBaryonNumber = bDef->GetBaryonNumber() + tDef->GetBaryonNumber();
257  G4int initCharge = bullet->getCharge() + target->getCharge();
258  G4int initStrangeness = bDef->GetQuarkContent(3) - bDef->GetAntiQuarkContent(3) +
259  tDef->GetQuarkContent(3) - tDef->GetAntiQuarkContent(3);
260 
261  G4int finalBaryonNumber = 0;
262  G4int finalCharge = 0;
263  G4int finalStrangeness = 0;
264 
265  for (ipart = particles.begin(); ipart != particles.end(); ipart++) {
266  finalBaryonNumber += ipart->baryon();
267  finalCharge += ipart->getCharge();
268  finalStrangeness += ipart->getStrangeness();
269  }
270 
271  G4int bnc = finalBaryonNumber - initBaryonNumber;
272  G4int cnc = finalCharge - initCharge;
273  G4int snc = finalStrangeness - initStrangeness;
274 
275  if (bnc != 0 || cnc != 0 || snc != 0) {
276  G4cout << " G4ElementaryParticleCollider: quantum number non-conservation " << G4endl;
277  G4cout << " Baryon number: initial = " << initBaryonNumber << ", final = "
278  << finalBaryonNumber << G4endl;
279  G4cout << " Charge: initial = " << initCharge << ", final = "
280  << finalCharge << G4endl;
281  G4cout << " Strangeness: initial = " << initStrangeness << ", final = "
282  << finalStrangeness << G4endl;
283 
284  G4cout << " bullet = " << bDef->GetParticleName() << G4endl;
285  G4cout << " target = " << tDef->GetParticleName() << G4endl;
286  G4cout << " secondaries = " ;
287  for (ipart = particles.begin(); ipart != particles.end(); ipart++) {
288  G4cout << ipart->getDefinition()->GetParticleName() << " " ;
289  }
290  G4cout << G4endl;
291  }
292  }
293 
294  if (verboseLevel && !validateOutput(bullet, target, particles)) {
295  G4cout << " incoming particles: \n" << *particle1 << G4endl
296  << *particle2 << G4endl
297  << " outgoing particles: " << G4endl;
298  for(ipart = particles.begin(); ipart != particles.end(); ipart++)
299  G4cout << *ipart << G4endl;
300 
301  G4cout << " <<< Non-conservation in G4ElementaryParticleCollider"
302  << G4endl;
303  }
304 
305  std::sort(particles.begin(), particles.end(), G4ParticleLargerEkin());
307 }
308 
309 
310 G4int
312  G4double ekin) const
313 {
314  G4int mul = 0;
315 
316  const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(is);
317 
318  if (xsecTable) mul = xsecTable->getMultiplicity(ekin);
319  else {
320  G4cerr << " G4ElementaryParticleCollider: Unknown interaction channel "
321  << is << " - multiplicity not generated " << G4endl;
322  }
323 
324  if(verboseLevel > 3){
325  G4cout << " G4ElementaryParticleCollider::generateMultiplicity: "
326  << " multiplicity = " << mul << G4endl;
327  }
328 
329  return mul;
330 }
331 
332 
333 void
335  G4double ekin)
336 {
337  particle_kinds.clear(); // Initialize buffer for generation
338 
339  const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(is);
340 
341  if (xsecTable)
342  xsecTable->getOutgoingParticleTypes(particle_kinds, mult, ekin);
343  else {
344  G4cerr << " G4ElementaryParticleCollider: Unknown interaction channel "
345  << is << " - outgoing kinds not generated " << G4endl;
346  }
347 
348  return;
349 }
350 
351 
352 void
354  G4double etot_scm,
355  G4InuclElementaryParticle* particle1,
356  G4InuclElementaryParticle* particle2) {
357  if (verboseLevel > 2) {
358  G4cout << " >>> G4ElementaryParticleCollider::generateSCMfinalState"
359  << G4endl;
360  }
361 
363 
364  const G4int itry_max = 10;
365 
366  G4int type1 = particle1->type();
367  G4int type2 = particle2->type();
368 
369  G4int is = type1 * type2;
370 
371  if (verboseLevel > 3) G4cout << " is " << is << G4endl;
372 
373  G4int multiplicity = 0;
374  G4bool generate = true;
375 
376  G4int itry = 0;
377  while (generate && itry++ < itry_max) { /* Loop checking 08.06.2015 MHK */
378  particles.clear(); // Initialize buffers for this event
379  particle_kinds.clear();
380 
381  // Generate list of final-state particles
382  multiplicity = generateMultiplicity(is, ekin);
383 
384  generateOutgoingPartTypes(is, multiplicity, ekin);
385  if (particle_kinds.empty()) {
386  if (verboseLevel > 3) {
387  G4cout << " generateOutgoingPartTypes failed mult " << multiplicity
388  << G4endl;
389  }
390  continue;
391  }
392 
393  fillOutgoingMasses(); // Fill mass buffer from particle types
394 
395  // Attempt to produce final state kinematics
396  fsGenerator.Configure(particle1, particle2, particle_kinds);
397  generate = !fsGenerator.Generate(etot_scm, masses, scm_momentums);
398  } // while (generate)
399 
400  if (itry >= itry_max) { // Unable to generate valid final state
401  if (verboseLevel > 2)
402  G4cout << " generateSCMfinalState failed " << itry << " attempts"
403  << G4endl;
404  return;
405  }
406 
407  // Store generated momenta into outgoing particles
408 
409  particles.resize(multiplicity); // Preallocate buffer
410  for (G4int i=0; i<multiplicity; i++) {
411  particles[i].fill(scm_momentums[i], particle_kinds[i],
413  }
414 
415  if (verboseLevel > 3) {
416  G4cout << " <<< G4ElementaryParticleCollider::generateSCMfinalState"
417  << G4endl;
418  }
419 
420  return; // Particles buffer filled
421 }
422 
423 
424 // Use generated list of final states to fill mass buffers
425 
427  G4int mult = particle_kinds.size();
428 
429  masses.resize(mult,0.);
430  masses2.resize(mult,0.); // Allows direct [i] setting
431 
432  for (G4int i = 0; i < mult; i++) {
434  masses2[i] = masses[i] * masses[i];
435  }
436 }
437 
438 
439 // Generate nucleon momenta for pion or photon absorption by dibaryon
440 // The nucleon distribution is assumed to be isotropic in SCM
441 
442 void
444  G4InuclElementaryParticle* particle1,
445  G4InuclElementaryParticle* particle2)
446 {
447  if (verboseLevel > 3)
448  G4cout << " >>> G4ElementaryParticleCollider::generateSCMpionAbsorption"
449  << G4endl;
450 
451  particles.clear(); // Initialize buffers for this event
452  particles.resize(2);
453 
454  particle_kinds.clear();
455 
456  G4int typeProduct = particle1->type() * particle2->type();
457 
458  if (typeProduct == pi0*diproton || typeProduct == pip*unboundPN ||
459  typeProduct == gam*diproton) {
460  particle_kinds.push_back(pro);
461  particle_kinds.push_back(pro);
462 
463  } else if (typeProduct == pim*diproton || typeProduct == pip*dineutron ||
464  typeProduct == pi0*unboundPN || typeProduct == gam*unboundPN) {
465  particle_kinds.push_back(pro);
466  particle_kinds.push_back(neu);
467 
468  } else if (typeProduct == pi0*dineutron || typeProduct == pim*unboundPN ||
469  typeProduct == gam*dineutron) {
470  particle_kinds.push_back(neu);
471  particle_kinds.push_back(neu);
472 
473  } else {
474  G4cerr << " Illegal absorption: "
475  << particle1->getDefinition()->GetParticleName() << " + "
476  << particle2->getDefinition()->GetParticleName() << " -> ?"
477  << G4endl;
478  return;
479  }
480 
482 
483  G4double a = 0.5 * (etot_scm * etot_scm - masses2[0] - masses2[1]);
484 
485  G4double pmod = std::sqrt((a*a - masses2[0]*masses2[1])
486  / (etot_scm*etot_scm) );
488  G4LorentzVector mom2;
489  mom2.setVectM(-mom1.vect(), masses[1]);
490 
493 }
494 
495 
496 // Generate nucleon momenta for muon absorption by dibaryon
497 // The nucleon distribution is assumed to be isotropic in SCM
498 
499 void
501  G4InuclElementaryParticle* particle1,
502  G4InuclElementaryParticle* particle2)
503 {
504  if (verboseLevel > 3)
505  G4cout << " >>> G4ElementaryParticleCollider::generateSCMmuonAbsorption"
506  << G4endl;
507 
508  particles.clear(); // Initialize buffers for this event
509  particles.resize(3);
510 
511  scm_momentums.clear();
512  scm_momentums.resize(3);
513 
514  particle_kinds.clear();
515 
516  G4int typeProduct = particle1->type() * particle2->type();
517 
518  if (typeProduct == mum*diproton) {
519  particle_kinds.push_back(pro);
520  particle_kinds.push_back(neu);
521  } else if (typeProduct == mum*unboundPN) {
522  particle_kinds.push_back(neu);
523  particle_kinds.push_back(neu);
524  } else {
525  G4cerr << " Illegal absorption: "
526  << particle1->getDefinition()->GetParticleName() << " + "
527  << particle2->getDefinition()->GetParticleName() << " -> ?"
528  << G4endl;
529  return;
530  }
531  particle_kinds.push_back(mnu);
532 
534 
535  // Phase space generator required for the 3-body final state
536  G4GDecay3 breakup(etot_scm, masses[0], masses[1], masses[2]);
537  std::vector<G4ThreeVector> theMomenta = breakup.GetThreeBodyMomenta();
538 
539  if (theMomenta.empty()) {
540  G4cerr << " generateSCMmuonAbsorption: GetThreeBodyMomenta() failed"
541  << " for " << particle2->type() << " dibaryon" << G4endl;
542  particle_kinds.clear();
543  masses.clear();
544  particles.clear();
545  return;
546  }
547 
548  for (size_t i=0; i<3; i++) {
549  scm_momentums[i].setVectM(theMomenta[i], masses[i]);
551  }
552 }
553 
554 
555 // generate nucleons momenta for pion absorption by single nucleon
556 
557 void
559  G4InuclElementaryParticle* particle1,
560  G4InuclElementaryParticle* particle2) {
561  if (verboseLevel > 3)
562  G4cout << " >>> G4ElementaryParticleCollider::generateSCMpionNAbsorption"
563  << G4endl;
564 
565  particles.clear(); // Initialize buffers for this event
566  particles.resize(1);
567 
568  particle_kinds.clear();
569 
570  G4int type1 = particle1->type();
571  G4int type2 = particle2->type();
572 
573  // Ensure that single-nucleon absportion is valid (charge exchangeable)
574  if ((type1*type2 != pim*pro && type1*type2 != pip*neu)) {
575  G4cerr << " pion-nucleon absorption: "
576  << particle1->getDefinition()->GetParticleName() << " + "
577  << particle2->getDefinition()->GetParticleName() << " -> ?"
578  << G4endl;
579  return;
580  }
581 
582  // Get outgoing nucleon type using charge exchange
583  // Proton code is 1, neutron code is 2, so 3-# swaps them
584  G4int ntype = (particle2->nucleon() ? type2 : type1);
585  G4int outType = 3 - ntype;
586  particle_kinds.push_back(outType);
587 
589 
590  // Get mass of residual nucleus (2-ntype = 1 for proton, 0 for neutron)
591  G4double mRecoil =
593  G4double mRecoil2 = mRecoil*mRecoil;
594 
595  // Recompute Ecm to include nucleus (for recoil kinematics)
596  G4LorentzVector piN4 = particle1->getMomentum() + particle2->getMomentum();
597  G4LorentzVector vsum(0.,0.,0.,mRecoil);
598  vsum += piN4;
599 
600  // Two-body kinematics (nucleon against nucleus) in overall CM system
601  G4double esq_scm = vsum.m2();
602  G4double a = 0.5 * (esq_scm - masses2[0] - mRecoil2);
603 
604  G4double pmod = std::sqrt((a*a - masses2[0]*mRecoil2) / esq_scm );
606 
607  if (verboseLevel > 3) {
608  G4cout << " outgoing type " << outType << " recoiling on nuclear mass "
609  << mRecoil << "\n a " << a << " p " << pmod << " Ekin "
610  << mom1.e()-masses[0] << G4endl;
611  }
612 
613  mom1.boost(-piN4.boostVector()); // Boost into CM of pi-N collision
614 
615  if (verboseLevel > 3) {
616  G4cout << " in original pi-N frame p(SCM) " << mom1.rho() << " Ekin "
617  << mom1.e()-masses[0] << G4endl;
618  }
619 
620  // Fill only the ejected nucleon
622 }
623 
624 
625 // Evaluate whether interaction is candidate for absorption on nucleon
626 G4bool
628  if (verboseLevel > 3)
629  G4cout << " >>> G4ElementaryParticleCollider::pionNucleonAbsorption ?"
630  << " ekin " << ekin << " is " << interCase.hadrons() << G4endl;
631 
632  // Absorption occurs with specified probability
634 
635  // Absorption occurs only for pi- p -> n, or pi+ n -> p
636  // Restrict to "very slow" pions, to allow for some normal scattering
637  return ((interCase.hadrons() == pim*pro ||
638  interCase.hadrons() == pip*neu)
639  && (ekin < 0.05) // 50 MeV kinetic energy or less
640  && (G4UniformRand() < absProb)
641  );
642 }
643 
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static G4double getParticleMass(G4int type)
const XML_Char * target
Definition: expat.h:268
Int_t ipart
Definition: Style.C:10
void SetVerboseLevel(G4int verbose)
G4int hadrons() const
void setBullet(const G4InuclParticle *bullet)
virtual void getOutgoingParticleTypes(std::vector< G4int > &kinds, G4int mult, G4double ke) const =0
#define G4endl
Definition: G4ios.hh:61
virtual G4int getMultiplicity(G4double ke) const =0
const G4String & GetParticleName() const
static void Print(std::ostream &os=G4cout)
void setVerbose(G4int vb=0)
G4LorentzVector generateWithRandomAngles(G4double p, G4double mass=0.)
void generateOutgoingPartTypes(G4int is, G4int mult, G4double ekin)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
const G4ParticleDefinition * getDefinition() const
void set(G4InuclParticle *part1, G4InuclParticle *part2)
void generate(const G4double sqrtS, ParticleList &particles)
Generate an event in the CM system.
G4CascadeFinalStateGenerator fsGenerator
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
static const G4CascadeChannel * GetTable(G4int initialState)
void generateSCMmuonAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4InuclElementaryParticle > particles
virtual G4bool useEPCollider(G4InuclParticle *bullet, G4InuclParticle *target) const
void generateSCMpionNAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
#define G4UniformRand()
Definition: Randomize.hh:53
G4int GetQuarkContent(G4int flavor) const
double m2() const
static G4bool useQuasiDeuteron(G4int ptype, G4int qdtype=0)
G4GLOB_DLL std::ostream G4cerr
static G4double piNAbsorption()
double rho() const
std::vector< G4ThreeVector > GetThreeBodyMomenta()
Definition: G4GDecay3.cc:100
void setTarget(const G4InuclParticle *target)
int G4int
Definition: G4Types.hh:78
G4double getNucleiMass() const
void generateSCMpionAbsorption(G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
G4bool pionNucleonAbsorption(G4double ekin) const
G4bool Generate(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
G4double getCharge() const
G4LorentzVector getMomentum() const
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4double getKinEnergyInTheTRS() const
G4GLOB_DLL std::ostream G4cout
Hep3Vector boostVector() const
Hep3Vector vect() const
G4double getTotalSCMEnergy() const
G4int generateMultiplicity(G4int is, G4double ekin) const
void generateSCMfinalState(G4double ekin, G4double etot_scm, G4InuclElementaryParticle *particle1, G4InuclElementaryParticle *particle2)
G4int GetAntiQuarkContent(G4int flavor) const
std::vector< G4LorentzVector > scm_momentums
G4double getMomModule() const
virtual G4bool validateOutput(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
HepLorentzVector & boost(double, double, double)
void setVectM(const Hep3Vector &spatial, double mass)
void Configure(G4InuclElementaryParticle *bullet, G4InuclElementaryParticle *target, const std::vector< G4int > &particle_kinds)