Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4InuclParticle.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4InuclParticle.cc 71719 2013-06-21 00:01:54Z mkelsey $
27 //
28 // 20100409 M. Kelsey -- Drop unused string argument from ctors.
29 // 20110721 M. Kelsey -- Add model ID as optional ctor argument (so subclasses
30 // don't have to call SetModel()).
31 // 20110922 M. Kelsey -- Add stream argument to printParticle() => print()
32 // 20130620 M. Kelsey -- Address Coverity #37504, check self in op=()
33 // 20130806 M. Kelsey -- Per A. Dotti, move empty G4DynPart to file scope to
34 // address thread-collision problem.
35 
36 #include <cmath>
37 
38 #include "G4InuclParticle.hh"
39 #include "G4ios.hh"
40 #include "G4SystemOfUnits.hh"
41 
42 // Internal constructor only usable by subclasses
44  const G4LorentzVector& mom,
46  : modelId(model) {
47  setDefinition(pd);
48  setMomentum(mom);
49 }
50 
51 
52 // Assignment operator for use with std::sort()
54  if (this != &right) {
55  pDP = right.pDP;
56  modelId = right.modelId;
57  }
58 
59  return *this;
60 }
61 
62 
63 // Set particle definition allowing for null pointer to erase DynPart content
64 
65 namespace {
66  static const G4DynamicParticle empty; // To zero out everything
67 }
68 
70  if (pd) pDP.SetDefinition(pd);
71  else pDP = empty;
72 }
73 
74 
75 // WARNING! Bertini code doesn't do four-vectors; repair mass before use!
77  G4double mass = getMass();
78  if (std::fabs(mass-mom.m()) <= 1e-5)
79  pDP.Set4Momentum(mom*GeV/MeV); // From Bertini to G4 units
80  else
81  pDP.SetMomentum(mom.vect()*GeV/MeV); // Don't change current mass!
82 }
83 
84 
85 // Proper stream output (just calls print())
86 
87 std::ostream& operator<<(std::ostream& os, const G4InuclParticle& part) {
88  part.print(os);
89  return os;
90 }
91 
92 void G4InuclParticle::print(std::ostream& os) const {
94  os << " px " << mom.px() << " py " << mom.py() << " pz " << mom.pz()
95  << " pmod " << mom.rho() << " E " << mom.e()
96  << " creator model " << modelId;
97 }
98 
void SetMomentum(const G4ThreeVector &momentum)
G4LorentzVector getMomentum() const
Hep3Vector vect() const
virtual void print(std::ostream &os) const
double py() const
double px() const
double rho() const
const XML_Char XML_Content * model
void setDefinition(G4ParticleDefinition *pd)
void Set4Momentum(const G4LorentzVector &momentum)
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
G4InuclParticle & operator=(const G4InuclParticle &right)
double pz() const
void setMomentum(const G4LorentzVector &mom)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double getMass() const