Geant4-11
G4CollisionOutput.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// 20100309 M. Kelsey -- Introduced bug checking i3 for valid tuning pair
29// 20100409 M. Kelsey -- Move non-inlinable code here out of .hh file, replace
30// loop over push_back() with block insert().
31// 20100418 M. Kelsey -- Add function to boost output lists to lab frame
32// 20100520 M. Kelsey -- Add function to rotate Z axis, from G4Casc.Interface
33// 20100620 M. Kelsey -- Add some diagnostics in setOnShell, simplify if's
34// 20100630 M. Kelsey -- Use "getExcitationEnergyInGeV()" instead of ".001*"
35// 20100701 M. Kelsey -- G4InuclNuclei now includes excitation energy as part
36// of the reported mass and four-vector.
37// 20100714 M. Kelsey -- Modify setOnShell() to avoid creating particles
38// with negative kinetic energy.
39// 20100715 M. Kelsey -- Add total charge and baryon number functions, and a
40// combined "add()" function to put two of these together
41// 20100924 M. Kelsey -- Use "OutgoingNuclei" name consistently, replacing
42// old "TargetFragment". Add new (reusable) G4Fragment buffer
43// and access functions for initial post-cascade processing.
44// Move implementation of add() to .cc file.
45// 20101019 M. Kelsey -- CoVerity report: unitialized constructor
46// 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
47// 20110225 M. Kelsey -- Add non-const functions to remove list elements
48// 20110302 M. Kelsey -- Fix message in setOnShell() by moving ini_mom calc
49// 20110307 M. Kelsey -- Need to discard existing ouput lists in trivialize()
50// 20110311 M. Kelsey -- Include nuclear fragment in setOnShell balancing,
51// including calculation of final-state momentum
52// 20110519 M. Kelsey -- Drop unused "p2" variable from selectPairToTune()
53// 20110801 M. Kelsey -- Use resize to avoid temporaries when copying from
54// G4ReactionProductVector
55// 20110922 M. Kelsey -- Follow G4InuclParticle::print(ostream&) migration,
56// Add optional stream argument to printCollisionOutput
57// 20121002 M. Kelsey -- Add strangeness calculation
58// 20130628 M. Kelsey -- Support multiple recoil fragments (for G4Fissioner)
59// 20141208 M. Kelsey -- Split function to do pair-wise "hard" tuning
60
61#include <algorithm>
62
63#include "G4CollisionOutput.hh"
64#include "G4SystemOfUnits.hh"
66#include "G4CascadParticle.hh"
68#include "G4Electron.hh"
70#include "G4LorentzConvertor.hh"
71#include "G4LorentzRotation.hh"
72#include "G4LorentzVector.hh"
74#include "G4ReactionProduct.hh"
75#include "G4ThreeVector.hh"
76
77typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
78typedef std::vector<G4InuclNuclei>::iterator nucleiIterator;
79typedef std::vector<G4Fragment>::iterator fragmentIterator;
80
82
83
85 : verboseLevel(0), eex_rest(0), on_shell(false) {
86 if (verboseLevel > 1)
87 G4cout << " >>> G4CollisionOutput::G4CollisionOutput" << G4endl;
88}
89
90
92{
93 if (this != &right) {
98 eex_rest = right.eex_rest;
99 on_shell = right.on_shell;
100 }
101 return *this;
102}
103
105 outgoingNuclei.clear();
106 outgoingParticles.clear();
107 recoilFragments.clear();
108 eex_rest = 0.;
109 on_shell = false;
110}
111
112
113// Get requested recoil fragment from list, or return dummy version
114
116 return ( (index >= 0 && index < numberOfFragments())
117 ? recoilFragments[index] : emptyFragment);
118}
119
120
121// Merge two complete objects
122
126 recoilFragments = right.recoilFragments; // Replace, don't combine
127 eex_rest = 0.;
128 on_shell = false;
129}
130
131
132// Append to lists
133
134void G4CollisionOutput::addOutgoingParticles(const std::vector<G4InuclElementaryParticle>& particles) {
136 particles.begin(), particles.end());
137}
138
139void G4CollisionOutput::addOutgoingNuclei(const std::vector<G4InuclNuclei>& nuclea) {
140 outgoingNuclei.insert(outgoingNuclei.end(), nuclea.begin(), nuclea.end());
141}
142
143// These are primarily for G4IntraNucleiCascader internal checks
144
146 addOutgoingParticle(cparticle.getParticle());
147}
148
149void G4CollisionOutput::addOutgoingParticles(const std::vector<G4CascadParticle>& cparticles) {
150 for (unsigned i=0; i<cparticles.size(); i++)
151 addOutgoingParticle(cparticles[i]);
152}
153
154// This comes from PreCompound de-excitation, both particles and nuclei
155
157 if (!rproducts) return; // Sanity check, no error if null
158
159 if (verboseLevel) {
160 G4cout << " >>> G4CollisionOutput::addOutgoingParticles(G4RPVector)"
161 << G4endl;
162 }
163
164 G4ReactionProductVector::const_iterator j;
165 for (j=rproducts->begin(); j!=rproducts->end(); ++j) {
166 const G4ParticleDefinition* pd = (*j)->GetDefinition();
168
169 // FIXME: Momentum returned by value; extra copying here!
170 G4LorentzVector mom((*j)->GetMomentum(), (*j)->GetTotalEnergy());
171 mom /= GeV; // Convert from GEANT4 to Bertini units
172
173 if (verboseLevel>1)
174 G4cout << " Processing " << pd->GetParticleName() << " (" << type
175 << "), momentum " << mom << " GeV" << G4endl;
176
177 // Nucleons and nuclei are jumbled together in the list
178 // NOTE: Resize and set/fill avoid unnecessary temporary copies
179 if (type) {
182
183 if (verboseLevel>1) G4cout << outgoingParticles.back() << G4endl;
184 } else {
186 outgoingNuclei.back().fill(mom,pd->GetAtomicMass(),pd->GetAtomicNumber(),
188
189 if (verboseLevel>1) G4cout << outgoingNuclei.back() << G4endl;
190 }
191 }
192}
193
194
195// Removing elements from lists by index
196
198 if (index >= 0 && index < numberOfOutgoingParticles())
199 outgoingParticles.erase(outgoingParticles.begin()+(size_t)index);
200}
201
203 if (index >= 0 && index < numberOfOutgoingNuclei())
204 outgoingNuclei.erase(outgoingNuclei.begin()+(size_t)index);
205}
206
207// Remove elements from list by reference, or by value comparison
208
211 std::find(outgoingParticles.begin(), outgoingParticles.end(), particle);
212 if (pos != outgoingParticles.end()) outgoingParticles.erase(pos);
213}
214
217 std::find(outgoingNuclei.begin(), outgoingNuclei.end(), nuclei);
218 if (pos != outgoingNuclei.end()) outgoingNuclei.erase(pos);
219}
220
221// Remove specified recoil fragment(s) from buffer
222
224 if (index < 0) recoilFragments.clear();
225 else if (index < numberOfFragments())
226 recoilFragments.erase(recoilFragments.begin()+(size_t)index);
227}
228
229
230// Kinematics of final state, for recoil and conservation checks
231
233 if (verboseLevel > 1)
234 G4cout << " >>> G4CollisionOutput::getTotalOutputMomentum" << G4endl;
235
236 G4LorentzVector tot_mom;
237 G4int i(0);
238 for(i=0; i < numberOfOutgoingParticles(); i++) {
239 tot_mom += outgoingParticles[i].getMomentum();
240 }
241 for(i=0; i < numberOfOutgoingNuclei(); i++) {
242 tot_mom += outgoingNuclei[i].getMomentum();
243 }
244 for(i=0; i < numberOfFragments(); i++) {
245 tot_mom += recoilFragments[i].GetMomentum()/GeV; // Need Bertini units!
246 }
247
248 return tot_mom;
249}
250
252 if (verboseLevel > 1)
253 G4cout << " >>> G4CollisionOutput::getTotalCharge" << G4endl;
254
255 G4int charge = 0;
256 G4int i(0);
257
258 for (i = 0; i < numberOfOutgoingParticles(); i++)
259 charge += G4int(outgoingParticles[i].getCharge());
260
261 for (i = 0; i < numberOfOutgoingNuclei(); i++)
262 charge += G4int(outgoingNuclei[i].getCharge());
263
264 for (i = 0; i < numberOfFragments(); i++)
265 charge += recoilFragments[i].GetZ_asInt();
266
267 return charge;
268}
269
270
272 if (verboseLevel > 1)
273 G4cout << " >>> G4CollisionOutput::getTotalBaryonNumber" << G4endl;
274
275 G4int baryon = 0;
276 G4int i(0);
277 for(i=0; i < numberOfOutgoingParticles(); i++) {
278 baryon += outgoingParticles[i].baryon();
279 }
280 for(i=0; i < numberOfOutgoingNuclei(); i++) {
281 baryon += G4int(outgoingNuclei[i].getA());
282 }
283 for(i=0; i < numberOfFragments(); i++) {
284 baryon += recoilFragments[i].GetA_asInt();
285 }
286
287 return baryon;
288}
289
291 if (verboseLevel > 1)
292 G4cout << " >>> G4CollisionOutput::getTotalStrangeness" << G4endl;
293
294 G4int strange = 0;
295 G4int i(0);
296 for(i=0; i < numberOfOutgoingParticles(); i++) {
297 strange += outgoingParticles[i].getStrangeness();
298 }
299
300 return strange;
301}
302
303
304void G4CollisionOutput::printCollisionOutput(std::ostream& os) const {
305 os << " Output: " << G4endl
306 << " Outgoing Particles: " << numberOfOutgoingParticles() << G4endl;
307
308 G4int i(0);
309 for(i=0; i < numberOfOutgoingParticles(); i++)
310 os << outgoingParticles[i] << G4endl;
311
312 os << " Outgoing Nuclei: " << numberOfOutgoingNuclei() << G4endl;
313 for(i=0; i < numberOfOutgoingNuclei(); i++)
314 os << outgoingNuclei[i] << G4endl;
315 for(i=0; i < numberOfFragments(); i++)
316 os << recoilFragments[i] << G4endl;
317}
318
319
320// Boost particles and fragment to LAB -- "convertor" must already be configured
321
323 if (verboseLevel > 1)
324 G4cout << " >>> G4CollisionOutput::boostToLabFrame" << G4endl;
325
326 particleIterator ipart = outgoingParticles.begin();
327 for(; ipart != outgoingParticles.end(); ipart++) {
328 ipart->setMomentum(boostToLabFrame(ipart->getMomentum(), convertor));
329 }
330
331 std::sort(outgoingParticles.begin(), outgoingParticles.end(),
333
334 nucleiIterator inuc = outgoingNuclei.begin();
335 for (; inuc != outgoingNuclei.end(); inuc++) {
336 inuc->setMomentum(boostToLabFrame(inuc->getMomentum(), convertor));
337 }
338
339 // NOTE: Fragment momentum must be converted to and from Bertini units
340 G4LorentzVector fmom;
341 fragmentIterator ifrag = recoilFragments.begin();
342 for (; ifrag != recoilFragments.end(); ifrag++) {
343 fmom = ifrag->GetMomentum() / GeV;
344 ifrag->SetMomentum(boostToLabFrame(fmom, convertor)*GeV);
345 }
346}
347
348// Pass by value to allow direct (internal) manipulation
351 const G4LorentzConvertor& convertor) const {
352 if (convertor.reflectionNeeded()) mom.setZ(-mom.z());
353 mom = convertor.rotate(mom);
354 mom = convertor.backToTheLab(mom);
355
356 return mom;
357}
358
359// Apply LorentzRotation to all particles in event
360
362 if (verboseLevel > 1)
363 G4cout << " >>> G4CollisionOutput::rotateEvent" << G4endl;
364
365 particleIterator ipart = outgoingParticles.begin();
366 for(; ipart != outgoingParticles.end(); ipart++)
367 ipart->setMomentum(ipart->getMomentum()*=rotate);
368
369 nucleiIterator inuc = outgoingNuclei.begin();
370 for (; inuc != outgoingNuclei.end(); inuc++)
371 inuc->setMomentum(inuc->getMomentum()*=rotate);
372
373 fragmentIterator ifrag = recoilFragments.begin();
374 for (; ifrag != recoilFragments.end(); ifrag++) {
375 G4LorentzVector mom = ifrag->GetMomentum(); // Need copy
376 ifrag->SetMomentum(mom*=rotate);
377 }
378}
379
380
382 G4InuclParticle* target) {
383 if (verboseLevel > 1)
384 G4cout << " >>> G4CollisionOutput::trivialize" << G4endl;
385
386 reset(); // Discard existing output, replace with bullet/target
387
388 if (G4InuclNuclei* nuclei_target = dynamic_cast<G4InuclNuclei*>(target)) {
389 outgoingNuclei.push_back(*nuclei_target);
390 } else {
391 G4InuclElementaryParticle* particle =
392 dynamic_cast<G4InuclElementaryParticle*>(target);
393 outgoingParticles.push_back(*particle);
394 }
395
396 if (G4InuclNuclei* nuclei_bullet = dynamic_cast<G4InuclNuclei*>(bullet)) {
397 outgoingNuclei.push_back(*nuclei_bullet);
398 } else {
399 G4InuclElementaryParticle* particle =
400 dynamic_cast<G4InuclElementaryParticle*>(bullet);
401 outgoingParticles.push_back(*particle);
402 }
403}
404
405
407 G4InuclParticle* target) {
408 if (verboseLevel > 1)
409 G4cout << " >>> G4CollisionOutput::setOnShell" << G4endl;
410
411 const G4double accuracy = 0.00001; // momentum concerves at the level of 10 keV
412
413 on_shell = false;
414
415 G4LorentzVector ini_mom = bullet->getMomentum();
416 G4LorentzVector momt = target->getMomentum();
418 if(verboseLevel > 2){
419 G4cout << " bullet momentum = " << ini_mom.e() <<", "<< ini_mom.x() <<", "<< ini_mom.y()<<", "<< ini_mom.z()<<G4endl;
420 G4cout << " target momentum = " << momt.e()<<", "<< momt.x()<<", "<< momt.y()<<", "<< momt.z()<<G4endl;
421 G4cout << " Fstate momentum = " << out_mom.e()<<", "<< out_mom.x()<<", "<< out_mom.y()<<", "<< out_mom.z()<<G4endl;
422 }
423 // correction for internal conversion
424 G4LorentzVector el4mom(0.,0.,0.,electron_mass_c2/GeV);
425 for(G4int i=0; i < numberOfOutgoingParticles(); ++i) {
426 if (outgoingParticles[i].getDefinition() == G4Electron::Electron())
427 momt += el4mom;
428 }
429
430 ini_mom += momt;
431
432 mom_non_cons = ini_mom - out_mom;
433 G4double pnc = mom_non_cons.rho();
434 G4double enc = mom_non_cons.e();
435
437
438 if(verboseLevel > 2) {
440 G4cout << " momentum non conservation: " << G4endl
441 << " e " << enc << " p " << pnc << G4endl
442 << " remaining exitation " << eex_rest << G4endl;
443 }
444
445 if (std::fabs(enc) <= accuracy && pnc <= accuracy) {
446 on_shell = true;
447 return;
448 }
449
450 // Adjust "last" particle's four-momentum to balance event
451 // ONLY adjust particles with sufficient e or p to remain physical!
452
453 if (verboseLevel > 2) G4cout << " re-balancing four-momenta" << G4endl;
454
457 G4int nfrag = numberOfFragments();
458
459 G4LorentzVector last_mom; // Buffer to reduce memory churn
460
461 if (npart > 0) {
462 for (G4int ip=npart-1; ip>=0; ip--) {
463 if (outgoingParticles[ip].getKineticEnergy()+enc > 0.) {
464 last_mom = outgoingParticles[ip].getMomentum();
465 last_mom += mom_non_cons;
466 outgoingParticles[ip].setMomentum(last_mom);
467 break;
468 }
469 }
470 } else if (nnuc > 0) {
471 for (G4int in=nnuc-1; in>=0; in--) {
472 if (outgoingNuclei[in].getKineticEnergy()+enc > 0.) {
473 last_mom = outgoingNuclei[in].getMomentum();
474 last_mom += mom_non_cons;
475 outgoingNuclei[in].setMomentum(last_mom);
476 break;
477 }
478 }
479 } else if (nfrag > 0) {
480 for (G4int ifr=nfrag-1; ifr>=0; ifr--) {
481 // NOTE: G4Fragment does not use Bertini units!
482 last_mom = recoilFragments[ifr].GetMomentum()/GeV;
483 if ((last_mom.e()-last_mom.m())+enc > 0.) {
484 last_mom += mom_non_cons;
485 recoilFragments[ifr].SetMomentum(last_mom*GeV);
486 break;
487 }
488 }
489 }
490
491 // Recompute momentum non-conservation parameters
492 out_mom = getTotalOutputMomentum();
493 mom_non_cons = ini_mom - out_mom;
494 pnc = mom_non_cons.rho();
495 enc = mom_non_cons.e();
496
497 if(verboseLevel > 2){
499 G4cout << " momentum non conservation after (1): " << G4endl
500 << " e " << enc << " p " << pnc << G4endl;
501 }
502
503 // Can energy be balanced just with nuclear excitation?
504 G4bool need_hard_tuning = true;
505
506 G4double encMeV = mom_non_cons.e() / GeV; // Excitation below is in MeV
507 if (nfrag > 0) {
508 for (G4int i=0; i < nfrag; i++) {
509 G4double eex = recoilFragments[0].GetExcitationEnergy();
510 if (eex > 0.0 && eex + encMeV >= 0.0) {
511 // NOTE: G4Fragment doesn't have function to change excitation energy
512 // ==> theRecoilFragment.SetExcitationEnergy(eex+encMeV);
513
514 G4LorentzVector fragMom = recoilFragments[i].GetMomentum();
515 G4double newMass = fragMom.m() + encMeV; // .m() includes eex
516 fragMom.setVectM(fragMom.vect(), newMass);
517 need_hard_tuning = false;
518 break;
519 }
520 }
521 } else if (nnuc > 0) {
522 for (G4int i=0; i < nnuc; i++) {
523 G4double eex = outgoingNuclei[i].getExitationEnergy();
524 if (eex > 0.0 && eex + encMeV >= 0.0) {
525 outgoingNuclei[i].setExitationEnergy(eex+encMeV);
526 need_hard_tuning = false;
527 break;
528 }
529 }
530 if (need_hard_tuning && encMeV > 0.) {
531 outgoingNuclei[0].setExitationEnergy(encMeV);
532 need_hard_tuning = false;
533 }
534 }
535
536 if (!need_hard_tuning) {
537 on_shell = true;
538 return;
539 }
540
541 // Momentum (hard) tuning required for energy conservation
542 if (verboseLevel > 2) {
543 G4cout << " trying hard (particle-pair) tuning" << G4endl;
544 }
545 /*****
546 // Hard tuning of quasielastic particle against nucleus
547 if (npart == 1) {
548 if (verboseLevel > 2)
549 G4cout << " tuning particle 0 against recoil nucleus" << G4endl;
550
551 G4LorentzVector mom1 = outgoingParticles[0].getMomentum();
552 G4LorentzVector mom2 = recoilFragments[0].GetMomentum()/GeV;
553
554 // Preserve momentum direction of outgoing particle
555 G4ThreeVector pdir = mom1.vect().unit();
556
557 // Two-body kinematics (nucleon against nucleus) in overall CM system
558 G4double ecmsq = ini_mom.m2();
559 G4double a = 0.5 * (ecmsq - mom1.m2() - mom2.m2());
560 G4double pmod = std::sqrt((a*a - mom1.m2()*mom2.m2()) / ecmsq);
561
562 mom1.setVectM(pdir*pmod, mom1.m());
563 mom2.setVectM(-pdir*pmod, mom2.m());
564
565 // Boost CM momenta back into lab frame, assign back to particles
566 mom1.boost(-ini_mom.boostVector());
567 mom2.boost(-ini_mom.boostVector());
568
569 outgoingParticles[0].setMomentum(mom1);
570 recoilFragments[0].SetMomentum(mom2*GeV);
571 } else { // Hard tuning using pair of outgoing particles
572 *****/
573 std::pair<std::pair<G4int, G4int>, G4int> tune_par = selectPairToTune(enc);
574 std::pair<G4int, G4int> tune_particles = tune_par.first;
575 G4int mom_ind = tune_par.second;
576
577 G4bool tuning_possible =
578 (tune_particles.first >= 0 && tune_particles.second >= 0 &&
579 mom_ind >= G4LorentzVector::X);
580
581 if (!tuning_possible) {
582 if (verboseLevel > 2) G4cout << " tuning impossible " << G4endl;
583 return;
584 }
585
586 if(verboseLevel > 2) {
587 G4cout << " p1 " << tune_particles.first << " p2 " << tune_particles.second
588 << " ind " << mom_ind << G4endl;
589 }
590
591 G4LorentzVector mom1 = outgoingParticles[tune_particles.first].getMomentum();
592 G4LorentzVector mom2 = outgoingParticles[tune_particles.second].getMomentum();
593
594 if (tuneSelectedPair(mom1, mom2, mom_ind)) {
595 outgoingParticles[tune_particles.first ].setMomentum(mom1);
596 outgoingParticles[tune_particles.second].setMomentum(mom2);
597 }
598 else return;
599 //***** }
600
601 out_mom = getTotalOutputMomentum();
602 std::sort(outgoingParticles.begin(), outgoingParticles.end(), G4ParticleLargerEkin());
603 mom_non_cons = ini_mom - out_mom;
604 pnc = mom_non_cons.rho();
605 enc = mom_non_cons.e();
606
607 on_shell = (std::fabs(enc) < accuracy || pnc < accuracy);
608
609 if(verboseLevel > 2) {
610 G4cout << " momentum non conservation tuning: " << G4endl
611 << " e " << enc << " p " << pnc
612 << (on_shell?" success":" FAILURE") << G4endl;
613 }
614}
615
616
617// Returns excitation energy in GeV
618
620 eex_rest = 0.;
621 G4int i(0);
622 for (i=0; i < numberOfOutgoingNuclei(); i++)
623 eex_rest += outgoingNuclei[i].getExitationEnergyInGeV();
624 for (i=0; i < numberOfFragments(); i++)
625 eex_rest += recoilFragments[i].GetExcitationEnergy() / GeV;
626}
627
628
629std::pair<std::pair<G4int, G4int>, G4int>
631 if (verboseLevel > 2)
632 G4cout << " >>> G4CollisionOutput::selectPairToTune" << G4endl;
633
634 std::pair<G4int, G4int> tup(-1, -1);
635 G4int i3 = -1;
636 std::pair<std::pair<G4int, G4int>, G4int> tuner(tup, i3); // Set invalid
637
638 if (outgoingParticles.size() < 2) return tuner; // Nothing to do
639
640 G4int ibest1 = -1;
641 G4int ibest2 = -1;
642 G4double pbest = 0.0;
643 G4double pcut = 0.3 * std::sqrt(1.88 * std::fabs(de));
644 G4double p1 = 0.0;
645
646 for (G4int i = 0; i < G4int(outgoingParticles.size())-1; i++) {
647 G4LorentzVector mom1 = outgoingParticles[i].getMomentum();
648
649 for (G4int j = i+1; j < G4int(outgoingParticles.size()); j++) {
650 G4LorentzVector mom2 = outgoingParticles[j].getMomentum();
651
652 for (G4int l = G4LorentzVector::X; l<=G4LorentzVector::Z; l++) {
653 if (mom1[l]*mom2[l]<0.0) {
654 if (std::fabs(mom1[l])>pcut && std::fabs(mom2[l])>pcut) {
655 G4double psum = std::fabs(mom1[l]) + std::fabs(mom2[l]);
656
657 if (psum > pbest) {
658 ibest1 = i;
659 ibest2 = j;
660 i3 = l;
661 p1 = mom1[l];
662 pbest = psum;
663 } // psum > pbest
664 } // mom1 and mom2 > pcut
665 } // mom1 ~ -mom2
666 } // for (G4int l ...
667 } // for (G4int j ...
668 } // for (G4int i ...
669
670 if (i3 < 0) return tuner;
671
672 tuner.second = i3; // Momentum axis for tuning
673
674 // NOTE: Sign of de determines order for special case of p1==0.
675 if (de > 0.0) {
676 tuner.first.first = (p1>0.) ? ibest1 : ibest2;
677 tuner.first.second = (p1>0.) ? ibest2 : ibest1;
678 } else {
679 tuner.first.first = (p1<0.) ? ibest2 : ibest1;
680 tuner.first.second = (p1<0.) ? ibest1 : ibest2;
681 }
682
683 return tuner;
684}
685
686
688 G4LorentzVector& mom2,
689 G4int mom_ind) const {
690 if (verboseLevel > 2)
691 G4cout << " >>> G4CollisionOutput::tuneSelectedPair" << G4endl;
692
693 G4double newE12 = mom1.e() + mom2.e() + mom_non_cons.e();
694 G4double R = 0.5 * (newE12*newE12 + mom2.e()*mom2.e() - mom1.e()*mom1.e()) / newE12;
695 G4double Q = -(mom1[mom_ind] + mom2[mom_ind]) / newE12;
696 G4double UDQ = 1.0 / (Q*Q - 1.0);
697 G4double W = (R * Q + mom2[mom_ind]) * UDQ;
698 G4double V = (mom2.e()*mom2.e() - R*R) * UDQ;
699 G4double DET = W*W + V;
700
701 if (DET < 0.0) {
702 if (verboseLevel > 2) G4cout << " DET < 0 : tuning failed" << G4endl;
703 return false;
704 }
705
706 // Tuning allowed only for non-negative determinant
707 G4double x1 = -(W + std::sqrt(DET));
708 G4double x2 = -(W - std::sqrt(DET));
709
710 // choose the appropriate solution
711 G4bool xset = false;
712 G4double x = 0.0;
713
714 if (mom_non_cons.e() > 0.0) { // x has to be > 0.0
715 if (x1 > 0.0) {
716 if (R + Q * x1 >= 0.0) {
717 x = x1;
718 xset = true;
719 }
720 }
721 if (!xset && x2 > 0.0) {
722 if (R + Q * x2 >= 0.0) {
723 x = x2;
724 xset = true;
725 }
726 }
727 } else {
728 if (x1 < 0.0) {
729 if (R + Q * x1 >= 0.) {
730 x = x1;
731 xset = true;
732 }
733 }
734 if (!xset && x2 < 0.0) {
735 if (R + Q * x2 >= 0.0) {
736 x = x2;
737 xset = true;
738 }
739 }
740 } // if(mom_non_cons.e() > 0.0)
741
742 if (!xset) {
743 if (verboseLevel > 2) G4cout << " no appropriate solution found" << G4endl;
744 return false;
745 }
746
747 mom1[mom_ind] += x; // retune momentums
748 mom2[mom_ind] -= x;
749 return true;
750}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4InuclNuclei >::iterator nucleiIterator
std::vector< G4Fragment >::iterator fragmentIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
static const G4double pos
std::vector< G4ReactionProduct * > G4ReactionProductVector
static constexpr double GeV
Definition: G4SIunits.hh:203
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void setVectM(const Hep3Vector &spatial, double mass)
Hep3Vector vect() const
const G4InuclElementaryParticle & getParticle() const
G4int numberOfOutgoingParticles() const
G4LorentzVector mom_non_cons
G4int getTotalStrangeness() const
std::vector< G4InuclElementaryParticle > outgoingParticles
G4LorentzVector getTotalOutputMomentum() const
std::vector< G4InuclNuclei > outgoingNuclei
G4CollisionOutput & operator=(const G4CollisionOutput &right)
G4int getTotalBaryonNumber() const
std::pair< std::pair< G4int, G4int >, G4int > selectPairToTune(G4double de) const
void removeRecoilFragment(G4int index=-1)
void boostToLabFrame(const G4LorentzConvertor &convertor)
void rotateEvent(const G4LorentzRotation &rotate)
void removeOutgoingNucleus(G4int index)
void printCollisionOutput(std::ostream &os=G4cout) const
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
G4int getTotalCharge() const
const G4Fragment & getRecoilFragment(G4int index=0) const
void setOnShell(G4InuclParticle *bullet, G4InuclParticle *target)
G4int numberOfOutgoingNuclei() const
void add(const G4CollisionOutput &right)
G4bool tuneSelectedPair(G4LorentzVector &mom1, G4LorentzVector &mom2, G4int mom_index) const
void trivialise(G4InuclParticle *bullet, G4InuclParticle *target)
void addOutgoingNuclei(const std::vector< G4InuclNuclei > &nuclea)
G4int numberOfFragments() const
void removeOutgoingParticle(G4int index)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
static const G4Fragment emptyFragment
std::vector< G4Fragment > recoilFragments
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4LorentzVector getMomentum() const
G4bool reflectionNeeded() const
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
G4LorentzVector rotate(const G4LorentzVector &mom) const
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
const G4String & GetParticleName() const
float electron_mass_c2
Definition: hepunit.py:273
static double Q[]