Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
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 ( )

Definition at line 81 of file G4CollisionOutput.cc.

References G4cout, and G4endl.

82  : verboseLevel(0), eex_rest(0), on_shell(false) {
83  if (verboseLevel > 1)
84  G4cout << " >>> G4CollisionOutput::G4CollisionOutput" << G4endl;
85 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Function Documentation

G4bool G4CollisionOutput::acceptable ( ) const
inline

Definition at line 173 of file G4CollisionOutput.hh.

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

173 { return on_shell; };
void G4CollisionOutput::add ( const G4CollisionOutput right)

Definition at line 120 of file G4CollisionOutput.cc.

References addOutgoingNuclei(), and addOutgoingParticles().

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

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

Definition at line 136 of file G4CollisionOutput.cc.

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

136  {
137  outgoingNuclei.insert(outgoingNuclei.end(), nuclea.begin(), nuclea.end());
138 }
void G4CollisionOutput::addOutgoingNucleus ( const G4InuclNuclei nuclei)
inline

Definition at line 85 of file G4CollisionOutput.hh.

Referenced by G4EquilibriumEvaporator::deExcite().

85  {
86  outgoingNuclei.push_back(nuclei);
87  };
void G4CollisionOutput::addOutgoingParticle ( const G4InuclElementaryParticle particle)
inline
void G4CollisionOutput::addOutgoingParticle ( const G4CascadParticle cparticle)

Definition at line 142 of file G4CollisionOutput.cc.

References addOutgoingParticle(), and G4CascadParticle::getParticle().

142  {
143  addOutgoingParticle(cparticle.getParticle());
144 }
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
const G4InuclElementaryParticle & getParticle() const
void G4CollisionOutput::addOutgoingParticles ( const std::vector< G4InuclElementaryParticle > &  particles)

Definition at line 131 of file G4CollisionOutput.cc.

Referenced by add(), G4ElementaryParticleCollider::collide(), G4CascadeCheckBalance::collide(), G4BigBanger::deExcite(), G4CascadeDeexcitation::deExcite(), G4PreCompoundDeexcitation::deExcite(), G4IntraNucleiCascader::finishCascade(), and G4IntraNucleiCascader::setupCascade().

131  {
132  outgoingParticles.insert(outgoingParticles.end(),
133  particles.begin(), particles.end());
134 }
void G4CollisionOutput::addOutgoingParticles ( const std::vector< G4CascadParticle > &  cparticles)

Definition at line 146 of file G4CollisionOutput.cc.

References addOutgoingParticle().

146  {
147  for (unsigned i=0; i<cparticles.size(); i++)
148  addOutgoingParticle(cparticles[i]);
149 }
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
void G4CollisionOutput::addOutgoingParticles ( const G4ReactionProductVector rproducts)

Definition at line 153 of file G4CollisionOutput.cc.

References G4cout, G4endl, G4ParticleDefinition::GetAtomicMass(), G4ParticleDefinition::GetAtomicNumber(), G4ParticleDefinition::GetParticleName(), python.hepunit::GeV, numberOfOutgoingNuclei(), numberOfOutgoingParticles(), G4InuclParticle::PreCompound, and G4InuclElementaryParticle::type().

153  {
154  if (!rproducts) return; // Sanity check, no error if null
155 
156  if (verboseLevel) {
157  G4cout << " >>> G4CollisionOutput::addOutgoingParticles(G4RPVector)"
158  << G4endl;
159  }
160 
161  G4ReactionProductVector::const_iterator j;
162  for (j=rproducts->begin(); j!=rproducts->end(); ++j) {
163  G4ParticleDefinition* pd = (*j)->GetDefinition();
165 
166  // FIXME: Momentum returned by value; extra copying here!
167  G4LorentzVector mom((*j)->GetMomentum(), (*j)->GetTotalEnergy());
168  mom /= GeV; // Convert from GEANT4 to Bertini units
169 
170  if (verboseLevel>1)
171  G4cout << " Processing " << pd->GetParticleName() << " (" << type
172  << "), momentum " << mom << " GeV" << G4endl;
173 
174  // Nucleons and nuclei are jumbled together in the list
175  // NOTE: Resize and set/fill avoid unnecessary temporary copies
176  if (type) {
177  outgoingParticles.resize(numberOfOutgoingParticles()+1);
178  outgoingParticles.back().fill(mom, pd, G4InuclParticle::PreCompound);
179 
180  if (verboseLevel>1) G4cout << outgoingParticles.back() << G4endl;
181  } else {
182  outgoingNuclei.resize(numberOfOutgoingNuclei()+1);
183  outgoingNuclei.back().fill(mom,pd->GetAtomicMass(),pd->GetAtomicNumber(),
185 
186  if (verboseLevel>1) G4cout << outgoingNuclei.back() << G4endl;
187  }
188  }
189 }
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetAtomicNumber() const
G4GLOB_DLL std::ostream G4cout
G4int numberOfOutgoingParticles() const
G4int GetAtomicMass() const
G4int numberOfOutgoingNuclei() const
#define G4endl
Definition: G4ios.hh:61
void G4CollisionOutput::addRecoilFragment ( const G4Fragment aFragment)
inline

Definition at line 98 of file G4CollisionOutput.hh.

Referenced by G4NonEquilibriumEvaporator::deExcite(), G4Fissioner::deExcite(), G4EquilibriumEvaporator::deExcite(), and G4IntraNucleiCascader::finishCascade().

98  {
99  if (aFragment) addRecoilFragment(*aFragment);
100  }
void addRecoilFragment(const G4Fragment *aFragment)
void G4CollisionOutput::addRecoilFragment ( const G4Fragment aFragment)
inline

Definition at line 102 of file G4CollisionOutput.hh.

102  {
103  recoilFragments.push_back(aFragment);
104  }
void G4CollisionOutput::boostToLabFrame ( const G4LorentzConvertor convertor)

Definition at line 318 of file G4CollisionOutput.cc.

References G4cout, G4endl, python.hepunit::GeV, and sort().

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

318  {
319  if (verboseLevel > 1)
320  G4cout << " >>> G4CollisionOutput::boostToLabFrame" << G4endl;
321 
322  particleIterator ipart = outgoingParticles.begin();
323  for(; ipart != outgoingParticles.end(); ipart++) {
324  ipart->setMomentum(boostToLabFrame(ipart->getMomentum(), convertor));
325  }
326 
327  std::sort(outgoingParticles.begin(), outgoingParticles.end(),
329 
330  nucleiIterator inuc = outgoingNuclei.begin();
331  for (; inuc != outgoingNuclei.end(); inuc++) {
332  inuc->setMomentum(boostToLabFrame(inuc->getMomentum(), convertor));
333  }
334 
335  // NOTE: Fragment momentum must be converted to and from Bertini units
336  G4LorentzVector fmom;
337  fragmentIterator ifrag = recoilFragments.begin();
338  for (; ifrag != recoilFragments.end(); ifrag++) {
339  fmom = ifrag->GetMomentum() / GeV;
340  ifrag->SetMomentum(boostToLabFrame(fmom, convertor)*GeV);
341  }
342 }
subroutine sort(A, N)
Definition: dpm25nuc7.f:4670
G4GLOB_DLL std::ostream G4cout
void boostToLabFrame(const G4LorentzConvertor &convertor)
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
std::vector< G4Fragment >::iterator fragmentIterator
#define G4endl
Definition: G4ios.hh:61
G4LorentzVector G4CollisionOutput::boostToLabFrame ( G4LorentzVector  mom,
const G4LorentzConvertor convertor 
) const

Definition at line 346 of file G4CollisionOutput.cc.

References G4LorentzConvertor::backToTheLab(), G4LorentzConvertor::reflectionNeeded(), G4LorentzConvertor::rotate(), CLHEP::HepLorentzVector::setZ(), and CLHEP::HepLorentzVector::z().

347  {
348  if (convertor.reflectionNeeded()) mom.setZ(-mom.z());
349  mom = convertor.rotate(mom);
350  mom = convertor.backToTheLab(mom);
351 
352  return mom;
353 }
G4LorentzVector rotate(const G4LorentzVector &mom) const
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
G4bool reflectionNeeded() const
const std::vector<G4InuclNuclei>& G4CollisionOutput::getOutgoingNuclei ( ) const
inline
std::vector<G4InuclNuclei>& G4CollisionOutput::getOutgoingNuclei ( )
inline

Definition at line 140 of file G4CollisionOutput.hh.

140 { return outgoingNuclei; };
const std::vector<G4InuclElementaryParticle>& G4CollisionOutput::getOutgoingParticles ( ) const
inline
std::vector<G4InuclElementaryParticle>& G4CollisionOutput::getOutgoingParticles ( )
inline

Definition at line 130 of file G4CollisionOutput.hh.

130  {
131  return outgoingParticles;
132  };
const G4Fragment & G4CollisionOutput::getRecoilFragment ( G4int  index = 0) const

Definition at line 112 of file G4CollisionOutput.cc.

References numberOfFragments().

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

112  {
113  return ( (index >= 0 && index < numberOfFragments())
114  ? recoilFragments[index] : emptyFragment);
115 }
G4int numberOfFragments() const
const std::vector<G4Fragment>& G4CollisionOutput::getRecoilFragments ( ) const
inline

Definition at line 146 of file G4CollisionOutput.hh.

146  {
147  return recoilFragments;
148  };
std::vector<G4Fragment>& G4CollisionOutput::getRecoilFragments ( )
inline

Definition at line 150 of file G4CollisionOutput.hh.

150 { return recoilFragments; };
double G4CollisionOutput::getRemainingExitationEnergy ( ) const
inline

Definition at line 172 of file G4CollisionOutput.hh.

172 { return eex_rest; };
G4int G4CollisionOutput::getTotalBaryonNumber ( ) const

Definition at line 267 of file G4CollisionOutput.cc.

References G4InuclParticleNames::baryon(), G4cout, G4endl, numberOfFragments(), numberOfOutgoingNuclei(), and numberOfOutgoingParticles().

Referenced by G4CascadeCheckBalance::collide().

267  {
268  if (verboseLevel > 1)
269  G4cout << " >>> G4CollisionOutput::getTotalBaryonNumber" << G4endl;
270 
271  G4int baryon = 0;
272  G4int i(0);
273  for(i=0; i < numberOfOutgoingParticles(); i++) {
274  baryon += outgoingParticles[i].baryon();
275  }
276  for(i=0; i < numberOfOutgoingNuclei(); i++) {
277  baryon += G4int(outgoingNuclei[i].getA());
278  }
279  for(i=0; i < numberOfFragments(); i++) {
280  baryon += recoilFragments[i].GetA_asInt();
281  }
282 
283  return baryon;
284 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
#define G4endl
Definition: G4ios.hh:61
G4int numberOfFragments() const
G4int G4CollisionOutput::getTotalCharge ( ) const

Definition at line 248 of file G4CollisionOutput.cc.

References G4cout, G4endl, numberOfFragments(), numberOfOutgoingNuclei(), and numberOfOutgoingParticles().

Referenced by G4CascadeCheckBalance::collide().

248  {
249  if (verboseLevel > 1)
250  G4cout << " >>> G4CollisionOutput::getTotalCharge" << G4endl;
251 
252  G4int charge = 0;
253  G4int i(0);
254  for(i=0; i < numberOfOutgoingParticles(); i++) {
255  charge += G4int(outgoingParticles[i].getCharge());
256  }
257  for(i=0; i < numberOfOutgoingNuclei(); i++) {
258  charge += G4int(outgoingNuclei[i].getCharge());
259  }
260  for(i=0; i < numberOfFragments(); i++) {
261  charge += recoilFragments[i].GetZ_asInt();
262  }
263 
264  return charge;
265 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
#define G4endl
Definition: G4ios.hh:61
G4int numberOfFragments() const
G4LorentzVector G4CollisionOutput::getTotalOutputMomentum ( ) const

Definition at line 229 of file G4CollisionOutput.cc.

References G4cout, G4endl, python.hepunit::GeV, numberOfFragments(), numberOfOutgoingNuclei(), and numberOfOutgoingParticles().

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

229  {
230  if (verboseLevel > 1)
231  G4cout << " >>> G4CollisionOutput::getTotalOutputMomentum" << G4endl;
232 
233  G4LorentzVector tot_mom;
234  G4int i(0);
235  for(i=0; i < numberOfOutgoingParticles(); i++) {
236  tot_mom += outgoingParticles[i].getMomentum();
237  }
238  for(i=0; i < numberOfOutgoingNuclei(); i++) {
239  tot_mom += outgoingNuclei[i].getMomentum();
240  }
241  for(i=0; i < numberOfFragments(); i++) {
242  tot_mom += recoilFragments[i].GetMomentum()/GeV; // Need Bertini units!
243  }
244 
245  return tot_mom;
246 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
#define G4endl
Definition: G4ios.hh:61
G4int numberOfFragments() const
G4int G4CollisionOutput::getTotalStrangeness ( ) const

Definition at line 286 of file G4CollisionOutput.cc.

References G4cout, G4endl, and numberOfOutgoingParticles().

Referenced by G4CascadeCheckBalance::collide().

286  {
287  if (verboseLevel > 1)
288  G4cout << " >>> G4CollisionOutput::getTotalStrangeness" << G4endl;
289 
290  G4int strange = 0;
291  G4int i(0);
292  for(i=0; i < numberOfOutgoingParticles(); i++) {
293  strange += outgoingParticles[i].getStrangeness();
294  }
295 
296  return strange;
297 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4int numberOfOutgoingParticles() const
#define G4endl
Definition: G4ios.hh:61
G4int G4CollisionOutput::numberOfFragments ( ) const
inline
G4int G4CollisionOutput::numberOfOutgoingNuclei ( ) const
inline
G4int G4CollisionOutput::numberOfOutgoingParticles ( ) const
inline
G4CollisionOutput & G4CollisionOutput::operator= ( const G4CollisionOutput right)

Definition at line 88 of file G4CollisionOutput.cc.

89 {
90  if (this != &right) {
91  verboseLevel = right.verboseLevel;
92  outgoingParticles = right.outgoingParticles;
93  outgoingNuclei = right.outgoingNuclei;
94  recoilFragments = right.recoilFragments;
95  eex_rest = right.eex_rest;
96  on_shell = right.on_shell;
97  }
98  return *this;
99 }
void G4CollisionOutput::printCollisionOutput ( std::ostream &  os = G4cout) const

Definition at line 300 of file G4CollisionOutput.cc.

References G4endl, numberOfFragments(), numberOfOutgoingNuclei(), and numberOfOutgoingParticles().

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

300  {
301  os << " Output: " << G4endl
302  << " Outgoing Particles: " << numberOfOutgoingParticles() << G4endl;
303 
304  G4int i(0);
305  for(i=0; i < numberOfOutgoingParticles(); i++)
306  os << outgoingParticles[i] << G4endl;
307 
308  os << " Outgoing Nuclei: " << numberOfOutgoingNuclei() << G4endl;
309  for(i=0; i < numberOfOutgoingNuclei(); i++)
310  os << outgoingNuclei[i] << G4endl;
311  for(i=0; i < numberOfFragments(); i++)
312  os << recoilFragments[i] << G4endl;
313 }
int G4int
Definition: G4Types.hh:78
G4int numberOfOutgoingParticles() const
G4int numberOfOutgoingNuclei() const
#define G4endl
Definition: G4ios.hh:61
G4int numberOfFragments() const
void G4CollisionOutput::removeOutgoingNucleus ( G4int  index)

Definition at line 199 of file G4CollisionOutput.cc.

References numberOfOutgoingNuclei().

Referenced by removeOutgoingNucleus().

199  {
200  if (index >= 0 && index < numberOfOutgoingNuclei())
201  outgoingNuclei.erase(outgoingNuclei.begin()+(size_t)index);
202 }
G4int numberOfOutgoingNuclei() const
void G4CollisionOutput::removeOutgoingNucleus ( const G4InuclNuclei nuclei)

Definition at line 212 of file G4CollisionOutput.cc.

References G4InuclParticleNames::nuclei.

212  {
214  std::find(outgoingNuclei.begin(), outgoingNuclei.end(), nuclei);
215  if (pos != outgoingNuclei.end()) outgoingNuclei.erase(pos);
216 }
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
void G4CollisionOutput::removeOutgoingNucleus ( const G4InuclNuclei nuclei)
inline

Definition at line 116 of file G4CollisionOutput.hh.

References removeOutgoingNucleus().

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

Definition at line 194 of file G4CollisionOutput.cc.

References numberOfOutgoingParticles().

Referenced by removeOutgoingParticle().

194  {
195  if (index >= 0 && index < numberOfOutgoingParticles())
196  outgoingParticles.erase(outgoingParticles.begin()+(size_t)index);
197 }
G4int numberOfOutgoingParticles() const
void G4CollisionOutput::removeOutgoingParticle ( const G4InuclElementaryParticle particle)

Definition at line 206 of file G4CollisionOutput.cc.

206  {
208  std::find(outgoingParticles.begin(), outgoingParticles.end(), particle);
209  if (pos != outgoingParticles.end()) outgoingParticles.erase(pos);
210 }
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
void G4CollisionOutput::removeOutgoingParticle ( const G4InuclElementaryParticle particle)
inline

Definition at line 110 of file G4CollisionOutput.hh.

References removeOutgoingParticle().

110  {
111  if (particle) removeOutgoingParticle(*particle);
112  }
void removeOutgoingParticle(G4int index)
void G4CollisionOutput::removeRecoilFragment ( G4int  index = -1)

Definition at line 220 of file G4CollisionOutput.cc.

References numberOfFragments().

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

220  {
221  if (index < 0) recoilFragments.clear();
222  else if (index < numberOfFragments())
223  recoilFragments.erase(recoilFragments.begin()+(size_t)index);
224 }
G4int numberOfFragments() const
void G4CollisionOutput::reset ( )
void G4CollisionOutput::rotateEvent ( const G4LorentzRotation rotate)

Definition at line 357 of file G4CollisionOutput.cc.

References G4cout, G4endl, and rotate().

Referenced by G4CascadeInterface::ApplyYourself().

357  {
358  if (verboseLevel > 1)
359  G4cout << " >>> G4CollisionOutput::rotateEvent" << G4endl;
360 
361  particleIterator ipart = outgoingParticles.begin();
362  for(; ipart != outgoingParticles.end(); ipart++)
363  ipart->setMomentum(ipart->getMomentum()*=rotate);
364 
365  nucleiIterator inuc = outgoingNuclei.begin();
366  for (; inuc != outgoingNuclei.end(); inuc++)
367  inuc->setMomentum(inuc->getMomentum()*=rotate);
368 
369  fragmentIterator ifrag = recoilFragments.begin();
370  for (; ifrag != recoilFragments.end(); ifrag++) {
371  G4LorentzVector mom = ifrag->GetMomentum(); // Need copy
372  ifrag->SetMomentum(mom*=rotate);
373  }
374 }
subroutine rotate
Definition: dpm25nuc2.f:10457
G4GLOB_DLL std::ostream G4cout
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
std::vector< G4Fragment >::iterator fragmentIterator
#define G4endl
Definition: G4ios.hh:61
void G4CollisionOutput::setOnShell ( G4InuclParticle bullet,
G4InuclParticle target 
)

Definition at line 402 of file G4CollisionOutput.cc.

References CLHEP::HepLorentzVector::e(), G4cout, G4endl, G4InuclParticle::getMomentum(), getTotalOutputMomentum(), python.hepunit::GeV, CLHEP::HepLorentzVector::m(), numberOfFragments(), numberOfOutgoingNuclei(), numberOfOutgoingParticles(), printCollisionOutput(), CLHEP::HepLorentzVector::rho(), setRemainingExitationEnergy(), CLHEP::HepLorentzVector::setVectM(), sort(), CLHEP::HepLorentzVector::vect(), test::x, CLHEP::HepLorentzVector::X, CLHEP::HepLorentzVector::x(), CLHEP::HepLorentzVector::y(), and CLHEP::HepLorentzVector::z().

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

403  {
404  if (verboseLevel > 1)
405  G4cout << " >>> G4CollisionOutput::setOnShell" << G4endl;
406 
407  const G4double accuracy = 0.00001; // momentum concerves at the level of 10 keV
408 
409  on_shell = false;
410 
411  G4LorentzVector ini_mom = bullet->getMomentum();
412  G4LorentzVector momt = target->getMomentum();
414  if(verboseLevel > 2){
415  G4cout << " bullet momentum = " << ini_mom.e() <<", "<< ini_mom.x() <<", "<< ini_mom.y()<<", "<< ini_mom.z()<<G4endl;
416  G4cout << " target momentum = " << momt.e()<<", "<< momt.x()<<", "<< momt.y()<<", "<< momt.z()<<G4endl;
417  G4cout << " Fstate momentum = " << out_mom.e()<<", "<< out_mom.x()<<", "<< out_mom.y()<<", "<< out_mom.z()<<G4endl;
418  }
419 
420  ini_mom += momt;
421  G4LorentzVector mom_non_cons = ini_mom - out_mom;
422 
423  G4double pnc = mom_non_cons.rho();
424  G4double enc = mom_non_cons.e();
425 
427 
428  if(verboseLevel > 2) {
430  G4cout << " momentum non conservation: " << G4endl
431  << " e " << enc << " p " << pnc << G4endl
432  << " remaining exitation " << eex_rest << G4endl;
433  }
434 
435  if (std::fabs(enc) <= accuracy && pnc <= accuracy) {
436  on_shell = true;
437  return;
438  }
439 
440  // Adjust "last" particle's four-momentum to balance event
441  // ONLY adjust particles with sufficient e or p to remain physical!
442 
443  if (verboseLevel > 2) G4cout << " re-balancing four-momenta" << G4endl;
444 
446  G4int nnuc = numberOfOutgoingNuclei();
447  G4int nfrag = numberOfFragments();
448 
449  G4LorentzVector last_mom; // Buffer to reduce memory churn
450 
451  if (npart > 0) {
452  for (G4int ip=npart-1; ip>=0; ip--) {
453  if (outgoingParticles[ip].getKineticEnergy()+enc > 0.) {
454  last_mom = outgoingParticles[ip].getMomentum();
455  last_mom += mom_non_cons;
456  outgoingParticles[ip].setMomentum(last_mom);
457  break;
458  }
459  }
460  } else if (nnuc > 0) {
461  for (G4int in=nnuc-1; in>=0; in--) {
462  if (outgoingNuclei[in].getKineticEnergy()+enc > 0.) {
463  last_mom = outgoingNuclei[in].getMomentum();
464  last_mom += mom_non_cons;
465  outgoingNuclei[in].setMomentum(last_mom);
466  break;
467  }
468  }
469  } else if (nfrag > 0) {
470  for (G4int ifr=nfrag-1; ifr>=0; ifr--) {
471  // NOTE: G4Fragment does not use Bertini units!
472  last_mom = recoilFragments[ifr].GetMomentum()/GeV;
473  if ((last_mom.e()-last_mom.m())+enc > 0.) {
474  last_mom += mom_non_cons;
475  recoilFragments[ifr].SetMomentum(last_mom*GeV);
476  break;
477  }
478  }
479  }
480 
481  // Recompute momentum non-conservation parameters
482  out_mom = getTotalOutputMomentum();
483  mom_non_cons = ini_mom - out_mom;
484  pnc = mom_non_cons.rho();
485  enc = mom_non_cons.e();
486 
487  if(verboseLevel > 2){
489  G4cout << " momentum non conservation after (1): " << G4endl
490  << " e " << enc << " p " << pnc << G4endl;
491  }
492 
493  // Can energy be balanced just with nuclear excitation?
494  G4bool need_hard_tuning = true;
495 
496  G4double encMeV = mom_non_cons.e() / GeV; // Excitation below is in MeV
497  if (nfrag > 0) {
498  for (G4int i=0; i < nfrag; i++) {
499  G4double eex = recoilFragments[0].GetExcitationEnergy();
500  if (eex > 0.0 && eex + encMeV >= 0.0) {
501  // NOTE: G4Fragment doesn't have function to change excitation energy
502  // ==> theRecoilFragment.SetExcitationEnergy(eex+encMeV);
503 
504  G4LorentzVector fragMom = recoilFragments[i].GetMomentum();
505  G4double newMass = fragMom.m() + encMeV; // .m() includes eex
506  fragMom.setVectM(fragMom.vect(), newMass);
507  need_hard_tuning = false;
508  break;
509  }
510  }
511  } else if (nnuc > 0) {
512  for (G4int i=0; i < nnuc; i++) {
513  G4double eex = outgoingNuclei[i].getExitationEnergy();
514  if (eex > 0.0 && eex + encMeV >= 0.0) {
515  outgoingNuclei[i].setExitationEnergy(eex+encMeV);
516  need_hard_tuning = false;
517  break;
518  }
519  }
520  if (need_hard_tuning && encMeV > 0.) {
521  outgoingNuclei[0].setExitationEnergy(encMeV);
522  need_hard_tuning = false;
523  }
524  }
525 
526  if (!need_hard_tuning) {
527  on_shell = true;
528  return;
529  }
530 
531  // Momentum (hard) tuning required for energy conservation
532  if (verboseLevel > 2)
533  G4cout << " trying hard (particle-pair) tuning" << G4endl;
534 
535  std::pair<std::pair<G4int, G4int>, G4int> tune_par = selectPairToTune(mom_non_cons.e());
536  std::pair<G4int, G4int> tune_particles = tune_par.first;
537  G4int mom_ind = tune_par.second;
538 
539  if(verboseLevel > 2) {
540  G4cout << " p1 " << tune_particles.first << " p2 " << tune_particles.second
541  << " ind " << mom_ind << G4endl;
542  }
543 
544  G4bool tuning_possible =
545  (tune_particles.first >= 0 && tune_particles.second >= 0 &&
546  mom_ind >= G4LorentzVector::X);
547 
548  if (!tuning_possible) {
549  if (verboseLevel > 2) G4cout << " tuning impossible " << G4endl;
550  return;
551  }
552 
553  G4LorentzVector mom1 = outgoingParticles[tune_particles.first].getMomentum();
554  G4LorentzVector mom2 = outgoingParticles[tune_particles.second].getMomentum();
555  G4double newE12 = mom1.e() + mom2.e() + mom_non_cons.e();
556  G4double R = 0.5 * (newE12 * newE12 + mom2.e() * mom2.e() - mom1.e() * mom1.e()) / newE12;
557  G4double Q = -(mom1[mom_ind] + mom2[mom_ind]) / newE12;
558  G4double UDQ = 1.0 / (Q * Q - 1.0);
559  G4double W = (R * Q + mom2[mom_ind]) * UDQ;
560  G4double V = (mom2.e() * mom2.e() - R * R) * UDQ;
561  G4double DET = W * W + V;
562 
563  if (DET < 0.0) {
564  if (verboseLevel > 2) G4cout << " DET < 0 " << G4endl;
565  return;
566  }
567 
568  // Tuning allowed only for non-negative determinant
569  G4double x1 = -(W + std::sqrt(DET));
570  G4double x2 = -(W - std::sqrt(DET));
571 
572  // choose the appropriate solution
573  G4bool xset = false;
574  G4double x = 0.0;
575 
576  if(mom_non_cons.e() > 0.0) { // x has to be > 0.0
577  if(x1 > 0.0) {
578  if(R + Q * x1 >= 0.0) {
579  x = x1;
580  xset = true;
581  };
582  };
583  if(!xset && x2 > 0.0) {
584  if(R + Q * x2 >= 0.0) {
585  x = x2;
586  xset = true;
587  };
588  };
589  } else {
590  if(x1 < 0.0) {
591  if(R + Q * x1 >= 0.) {
592  x = x1;
593  xset = true;
594  };
595  };
596  if(!xset && x2 < 0.0) {
597  if(R + Q * x2 >= 0.0) {
598  x = x2;
599  xset = true;
600  };
601  };
602  } // if(mom_non_cons.e() > 0.0)
603 
604  if(!xset) {
605  if(verboseLevel > 2)
606  G4cout << " no appropriate solution found " << G4endl;
607  return;
608  } // if(xset)
609 
610  // retune momentums
611  mom1[mom_ind] += x;
612  mom2[mom_ind] -= x;
613  outgoingParticles[tune_particles.first ].setMomentum(mom1);
614  outgoingParticles[tune_particles.second].setMomentum(mom2);
615  out_mom = getTotalOutputMomentum();
616  std::sort(outgoingParticles.begin(), outgoingParticles.end(), G4ParticleLargerEkin());
617  mom_non_cons = ini_mom - out_mom;
618  pnc = mom_non_cons.rho();
619  enc = mom_non_cons.e();
620 
621  on_shell = (std::fabs(enc) < accuracy || pnc < accuracy);
622 
623  if(verboseLevel > 2) {
624  G4cout << " momentum non conservation tuning: " << G4endl
625  << " e " << enc << " p " << pnc
626  << (on_shell?" success":" FAILURE") << G4endl;
627  }
628 }
G4LorentzVector getMomentum() const
void printCollisionOutput(std::ostream &os=G4cout) const
subroutine sort(A, N)
Definition: dpm25nuc7.f:4670
int G4int
Definition: G4Types.hh:78
G4LorentzVector getTotalOutputMomentum() const
void setVectM(const Hep3Vector &spatial, double mass)
Hep3Vector vect() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4int numberOfOutgoingParticles() const
double rho() const
G4int numberOfOutgoingNuclei() const
#define G4endl
Definition: G4ios.hh:61
G4int numberOfFragments() const
double G4double
Definition: G4Types.hh:76
void G4CollisionOutput::setRemainingExitationEnergy ( )

Definition at line 633 of file G4CollisionOutput.cc.

References python.hepunit::GeV, numberOfFragments(), and numberOfOutgoingNuclei().

Referenced by setOnShell().

633  {
634  eex_rest = 0.;
635  G4int i(0);
636  for (i=0; i < numberOfOutgoingNuclei(); i++)
637  eex_rest += outgoingNuclei[i].getExitationEnergyInGeV();
638  for (i=0; i < numberOfFragments(); i++)
639  eex_rest += recoilFragments[i].GetExcitationEnergy() / GeV;
640 }
int G4int
Definition: G4Types.hh:78
G4int numberOfOutgoingNuclei() const
G4int numberOfFragments() const
void G4CollisionOutput::setVerboseLevel ( G4int  verbose)
inline
void G4CollisionOutput::trivialise ( G4InuclParticle bullet,
G4InuclParticle target 
)

Definition at line 377 of file G4CollisionOutput.cc.

References G4cout, G4endl, and reset().

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

378  {
379  if (verboseLevel > 1)
380  G4cout << " >>> G4CollisionOutput::trivialize" << G4endl;
381 
382  reset(); // Discard existing output, replace with bullet/target
383 
384  if (G4InuclNuclei* nuclei_target = dynamic_cast<G4InuclNuclei*>(target)) {
385  outgoingNuclei.push_back(*nuclei_target);
386  } else {
387  G4InuclElementaryParticle* particle =
388  dynamic_cast<G4InuclElementaryParticle*>(target);
389  outgoingParticles.push_back(*particle);
390  }
391 
392  if (G4InuclNuclei* nuclei_bullet = dynamic_cast<G4InuclNuclei*>(bullet)) {
393  outgoingNuclei.push_back(*nuclei_bullet);
394  } else {
395  G4InuclElementaryParticle* particle =
396  dynamic_cast<G4InuclElementaryParticle*>(bullet);
397  outgoingParticles.push_back(*particle);
398  }
399 }
const XML_Char * target
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

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