Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CollisionOutput Class Reference

#include <G4CollisionOutput.hh>

Public Member Functions

 G4CollisionOutput ()
 
G4CollisionOutputoperator= (const G4CollisionOutput &right)
 
void setVerboseLevel (G4int verbose)
 
void reset ()
 
void add (const G4CollisionOutput &right)
 
void addOutgoingParticle (const G4InuclElementaryParticle &particle)
 
void addOutgoingParticles (const std::vector< G4InuclElementaryParticle > &particles)
 
void addOutgoingNucleus (const G4InuclNuclei &nuclei)
 
void addOutgoingNuclei (const std::vector< G4InuclNuclei > &nuclea)
 
void addOutgoingParticle (const G4CascadParticle &cparticle)
 
void addOutgoingParticles (const std::vector< G4CascadParticle > &cparticles)
 
void addOutgoingParticles (const G4ReactionProductVector *rproducts)
 
void addRecoilFragment (const G4Fragment *aFragment)
 
void addRecoilFragment (const G4Fragment &aFragment)
 
void removeOutgoingParticle (G4int index)
 
void removeOutgoingParticle (const G4InuclElementaryParticle &particle)
 
void removeOutgoingParticle (const G4InuclElementaryParticle *particle)
 
void removeOutgoingNucleus (G4int index)
 
void removeOutgoingNucleus (const G4InuclNuclei &nuclei)
 
void removeOutgoingNucleus (const G4InuclNuclei *nuclei)
 
void removeRecoilFragment (G4int index=-1)
 
G4int numberOfOutgoingParticles () const
 
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles () const
 
std::vector< G4InuclElementaryParticle > & getOutgoingParticles ()
 
G4int numberOfOutgoingNuclei () const
 
const std::vector< G4InuclNuclei > & getOutgoingNuclei () const
 
std::vector< G4InuclNuclei > & getOutgoingNuclei ()
 
G4int numberOfFragments () const
 
const G4FragmentgetRecoilFragment (G4int index=0) const
 
const std::vector< G4Fragment > & getRecoilFragments () const
 
std::vector< G4Fragment > & getRecoilFragments ()
 
G4LorentzVector getTotalOutputMomentum () const
 
G4int getTotalCharge () const
 
G4int getTotalBaryonNumber () const
 
G4int getTotalStrangeness () const
 
void printCollisionOutput (std::ostream &os=G4cout) const
 
void boostToLabFrame (const G4LorentzConvertor &convertor)
 
G4LorentzVector boostToLabFrame (G4LorentzVector mom, const G4LorentzConvertor &convertor) const
 
void rotateEvent (const G4LorentzRotation &rotate)
 
void trivialise (G4InuclParticle *bullet, G4InuclParticle *target)
 
void setOnShell (G4InuclParticle *bullet, G4InuclParticle *target)
 
void setRemainingExitationEnergy ()
 
double getRemainingExitationEnergy () const
 
G4bool acceptable () const
 

Detailed Description

Definition at line 66 of file G4CollisionOutput.hh.

Constructor & Destructor Documentation

◆ G4CollisionOutput()

G4CollisionOutput::G4CollisionOutput ( )

Definition at line 84 of file G4CollisionOutput.cc.

85 : verboseLevel(0), eex_rest(0), on_shell(false) {
86 if (verboseLevel > 1)
87 G4cout << " >>> G4CollisionOutput::G4CollisionOutput" << G4endl;
88}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

Member Function Documentation

◆ acceptable()

G4bool G4CollisionOutput::acceptable ( ) const
inline

Definition at line 173 of file G4CollisionOutput.hh.

173{ return on_shell; };

Referenced by G4InuclCollider::collide(), and G4IntraNucleiCascader::finishCascade().

◆ add()

void G4CollisionOutput::add ( const G4CollisionOutput right)

Definition at line 123 of file G4CollisionOutput.cc.

123 {
124 addOutgoingParticles(right.outgoingParticles);
125 addOutgoingNuclei(right.outgoingNuclei);
126 recoilFragments = right.recoilFragments; // Replace, don't combine
127 eex_rest = 0.;
128 on_shell = false;
129}
void addOutgoingNuclei(const std::vector< G4InuclNuclei > &nuclea)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)

Referenced by G4InuclCollider::collide(), G4CascadeCheckBalance::collide(), G4CascadeDeexcitation::deExcite(), G4InuclCollider::deexcite(), G4IntraNucleiCascader::finalize(), and G4InuclCollider::rescatter().

◆ addOutgoingNuclei()

void G4CollisionOutput::addOutgoingNuclei ( const std::vector< G4InuclNuclei > &  nuclea)

Definition at line 139 of file G4CollisionOutput.cc.

139 {
140 outgoingNuclei.insert(outgoingNuclei.end(), nuclea.begin(), nuclea.end());
141}

Referenced by add(), and G4CascadeCheckBalance::collide().

◆ addOutgoingNucleus()

void G4CollisionOutput::addOutgoingNucleus ( const G4InuclNuclei nuclei)
inline

Definition at line 85 of file G4CollisionOutput.hh.

85 {
86 outgoingNuclei.push_back(nuclei);
87 };

Referenced by G4EquilibriumEvaporator::deExcite().

◆ addOutgoingParticle() [1/2]

void G4CollisionOutput::addOutgoingParticle ( const G4CascadParticle cparticle)

Definition at line 145 of file G4CollisionOutput.cc.

145 {
146 addOutgoingParticle(cparticle.getParticle());
147}
const G4InuclElementaryParticle & getParticle() const
void addOutgoingParticle(const G4InuclElementaryParticle &particle)

◆ addOutgoingParticle() [2/2]

◆ addOutgoingParticles() [1/3]

void G4CollisionOutput::addOutgoingParticles ( const G4ReactionProductVector rproducts)

Definition at line 156 of file G4CollisionOutput.cc.

156 {
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) {
180 outgoingParticles.resize(numberOfOutgoingParticles()+1);
181 outgoingParticles.back().fill(mom, pd, G4InuclParticle::PreCompound);
182
183 if (verboseLevel>1) G4cout << outgoingParticles.back() << G4endl;
184 } else {
185 outgoingNuclei.resize(numberOfOutgoingNuclei()+1);
186 outgoingNuclei.back().fill(mom,pd->GetAtomicMass(),pd->GetAtomicNumber(),
188
189 if (verboseLevel>1) G4cout << outgoingNuclei.back() << G4endl;
190 }
191 }
192}
int G4int
Definition: G4Types.hh:85
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
const G4String & GetParticleName() const

◆ addOutgoingParticles() [2/3]

void G4CollisionOutput::addOutgoingParticles ( const std::vector< G4CascadParticle > &  cparticles)

Definition at line 149 of file G4CollisionOutput.cc.

149 {
150 for (unsigned i=0; i<cparticles.size(); i++)
151 addOutgoingParticle(cparticles[i]);
152}

◆ addOutgoingParticles() [3/3]

void G4CollisionOutput::addOutgoingParticles ( const std::vector< G4InuclElementaryParticle > &  particles)

◆ addRecoilFragment() [1/2]

void G4CollisionOutput::addRecoilFragment ( const G4Fragment aFragment)
inline

Definition at line 102 of file G4CollisionOutput.hh.

102 {
103 recoilFragments.push_back(aFragment);
104 }

◆ addRecoilFragment() [2/2]

void G4CollisionOutput::addRecoilFragment ( const G4Fragment aFragment)
inline

◆ boostToLabFrame() [1/2]

void G4CollisionOutput::boostToLabFrame ( const G4LorentzConvertor convertor)

Definition at line 322 of file G4CollisionOutput.cc.

322 {
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}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4Fragment >::iterator fragmentIterator
void boostToLabFrame(const G4LorentzConvertor &convertor)

Referenced by boostToLabFrame(), G4InuclCollider::collide(), and G4EquilibriumEvaporator::deExcite().

◆ boostToLabFrame() [2/2]

G4LorentzVector G4CollisionOutput::boostToLabFrame ( G4LorentzVector  mom,
const G4LorentzConvertor convertor 
) const

Definition at line 350 of file G4CollisionOutput.cc.

351 {
352 if (convertor.reflectionNeeded()) mom.setZ(-mom.z());
353 mom = convertor.rotate(mom);
354 mom = convertor.backToTheLab(mom);
355
356 return mom;
357}
G4bool reflectionNeeded() const
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
G4LorentzVector rotate(const G4LorentzVector &mom) const

◆ getOutgoingNuclei() [1/2]

std::vector< G4InuclNuclei > & G4CollisionOutput::getOutgoingNuclei ( )
inline

Definition at line 140 of file G4CollisionOutput.hh.

140{ return outgoingNuclei; };

◆ getOutgoingNuclei() [2/2]

◆ getOutgoingParticles() [1/2]

std::vector< G4InuclElementaryParticle > & G4CollisionOutput::getOutgoingParticles ( )
inline

Definition at line 130 of file G4CollisionOutput.hh.

130 {
131 return outgoingParticles;
132 };

◆ getOutgoingParticles() [2/2]

◆ getRecoilFragment()

const G4Fragment & G4CollisionOutput::getRecoilFragment ( G4int  index = 0) const

Definition at line 115 of file G4CollisionOutput.cc.

115 {
116 return ( (index >= 0 && index < numberOfFragments())
117 ? recoilFragments[index] : emptyFragment);
118}
G4int numberOfFragments() const

Referenced by G4InuclCollider::collide(), G4CascadeDeexcitation::deExcite(), G4EquilibriumEvaporator::deExcite(), G4NonEquilibriumEvaporator::deExcite(), and G4InuclCollider::rescatter().

◆ getRecoilFragments() [1/2]

std::vector< G4Fragment > & G4CollisionOutput::getRecoilFragments ( )
inline

Definition at line 150 of file G4CollisionOutput.hh.

150{ return recoilFragments; };

◆ getRecoilFragments() [2/2]

const std::vector< G4Fragment > & G4CollisionOutput::getRecoilFragments ( ) const
inline

Definition at line 146 of file G4CollisionOutput.hh.

146 {
147 return recoilFragments;
148 };

◆ getRemainingExitationEnergy()

double G4CollisionOutput::getRemainingExitationEnergy ( ) const
inline

Definition at line 172 of file G4CollisionOutput.hh.

172{ return eex_rest; };

◆ getTotalBaryonNumber()

G4int G4CollisionOutput::getTotalBaryonNumber ( ) const

Definition at line 271 of file G4CollisionOutput.cc.

271 {
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}

Referenced by G4CascadeCheckBalance::collide().

◆ getTotalCharge()

G4int G4CollisionOutput::getTotalCharge ( ) const

Definition at line 251 of file G4CollisionOutput.cc.

251 {
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}

Referenced by G4CascadeCheckBalance::collide().

◆ getTotalOutputMomentum()

G4LorentzVector G4CollisionOutput::getTotalOutputMomentum ( ) const

Definition at line 232 of file G4CollisionOutput.cc.

232 {
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}

Referenced by G4CascadeCheckBalance::collide(), and setOnShell().

◆ getTotalStrangeness()

G4int G4CollisionOutput::getTotalStrangeness ( ) const

Definition at line 290 of file G4CollisionOutput.cc.

290 {
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}

Referenced by G4CascadeCheckBalance::collide().

◆ numberOfFragments()

G4int G4CollisionOutput::numberOfFragments ( ) const
inline

◆ numberOfOutgoingNuclei()

◆ numberOfOutgoingParticles()

◆ operator=()

G4CollisionOutput & G4CollisionOutput::operator= ( const G4CollisionOutput right)

Definition at line 91 of file G4CollisionOutput.cc.

92{
93 if (this != &right) {
94 verboseLevel = right.verboseLevel;
95 outgoingParticles = right.outgoingParticles;
96 outgoingNuclei = right.outgoingNuclei;
97 recoilFragments = right.recoilFragments;
98 eex_rest = right.eex_rest;
99 on_shell = right.on_shell;
100 }
101 return *this;
102}

◆ printCollisionOutput()

void G4CollisionOutput::printCollisionOutput ( std::ostream &  os = G4cout) const

Definition at line 304 of file G4CollisionOutput.cc.

304 {
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}

Referenced by G4CascadeInterface::ApplyYourself(), G4CascadeDeexcitation::deExcite(), G4EvaporationInuclCollider::deExcite(), G4CascadeCoalescence::FindClusters(), G4IntraNucleiCascader::finishCascade(), G4NucleiModel::generateParticleFate(), G4CascadeInterface::Propagate(), setOnShell(), G4CascadeInterface::throwNonConservationFailure(), and G4CascadeColliderBase::validateOutput().

◆ removeOutgoingNucleus() [1/3]

void G4CollisionOutput::removeOutgoingNucleus ( const G4InuclNuclei nuclei)

Definition at line 215 of file G4CollisionOutput.cc.

215 {
217 std::find(outgoingNuclei.begin(), outgoingNuclei.end(), nuclei);
218 if (pos != outgoingNuclei.end()) outgoingNuclei.erase(pos);
219}

◆ removeOutgoingNucleus() [2/3]

void G4CollisionOutput::removeOutgoingNucleus ( const G4InuclNuclei nuclei)
inline

Definition at line 116 of file G4CollisionOutput.hh.

116 {
117 if (nuclei) removeOutgoingNucleus(*nuclei);
118 }
void removeOutgoingNucleus(G4int index)

◆ removeOutgoingNucleus() [3/3]

void G4CollisionOutput::removeOutgoingNucleus ( G4int  index)

Definition at line 202 of file G4CollisionOutput.cc.

202 {
203 if (index >= 0 && index < numberOfOutgoingNuclei())
204 outgoingNuclei.erase(outgoingNuclei.begin()+(size_t)index);
205}

Referenced by removeOutgoingNucleus().

◆ removeOutgoingParticle() [1/3]

void G4CollisionOutput::removeOutgoingParticle ( const G4InuclElementaryParticle particle)

Definition at line 209 of file G4CollisionOutput.cc.

209 {
211 std::find(outgoingParticles.begin(), outgoingParticles.end(), particle);
212 if (pos != outgoingParticles.end()) outgoingParticles.erase(pos);
213}

◆ removeOutgoingParticle() [2/3]

void G4CollisionOutput::removeOutgoingParticle ( const G4InuclElementaryParticle particle)
inline

Definition at line 110 of file G4CollisionOutput.hh.

110 {
111 if (particle) removeOutgoingParticle(*particle);
112 }
void removeOutgoingParticle(G4int index)

◆ removeOutgoingParticle() [3/3]

void G4CollisionOutput::removeOutgoingParticle ( G4int  index)

Definition at line 197 of file G4CollisionOutput.cc.

197 {
198 if (index >= 0 && index < numberOfOutgoingParticles())
199 outgoingParticles.erase(outgoingParticles.begin()+(size_t)index);
200}

Referenced by removeOutgoingParticle().

◆ removeRecoilFragment()

void G4CollisionOutput::removeRecoilFragment ( G4int  index = -1)

Definition at line 223 of file G4CollisionOutput.cc.

223 {
224 if (index < 0) recoilFragments.clear();
225 else if (index < numberOfFragments())
226 recoilFragments.erase(recoilFragments.begin()+(size_t)index);
227}

Referenced by G4InuclCollider::collide(), and G4InuclCollider::rescatter().

◆ reset()

◆ rotateEvent()

void G4CollisionOutput::rotateEvent ( const G4LorentzRotation rotate)

Definition at line 361 of file G4CollisionOutput.cc.

361 {
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}

◆ setOnShell()

void G4CollisionOutput::setOnShell ( G4InuclParticle bullet,
G4InuclParticle target 
)

Definition at line 406 of file G4CollisionOutput.cc.

407 {
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}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
void setVectM(const Hep3Vector &spatial, double mass)
Hep3Vector vect() const
G4LorentzVector getTotalOutputMomentum() const
void printCollisionOutput(std::ostream &os=G4cout) const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4LorentzVector getMomentum() const

Referenced by G4InuclCollider::collide(), and G4IntraNucleiCascader::finishCascade().

◆ setRemainingExitationEnergy()

void G4CollisionOutput::setRemainingExitationEnergy ( )

Definition at line 619 of file G4CollisionOutput.cc.

619 {
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}

Referenced by setOnShell().

◆ setVerboseLevel()

void G4CollisionOutput::setVerboseLevel ( G4int  verbose)
inline

◆ trivialise()

void G4CollisionOutput::trivialise ( G4InuclParticle bullet,
G4InuclParticle target 
)

Definition at line 381 of file G4CollisionOutput.cc.

382 {
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}

Referenced by G4InuclCollider::collide(), G4LightTargetCollider::collide(), and G4IntraNucleiCascader::finalize().


The documentation for this class was generated from the following files: