Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4VElasticCollision Class Reference

#include <G4VElasticCollision.hh>

Inheritance diagram for G4VElasticCollision:
G4VCollision G4CollisionMesonBaryonElastic G4CollisionNNElastic G4CollisionnpElastic

Public Member Functions

 G4VElasticCollision ()
 
virtual ~G4VElasticCollision ()
 
G4bool operator== (const G4VElasticCollision &right) const
 
G4bool operator!= (const G4VElasticCollision &right) const
 
virtual G4KineticTrackVectorFinalState (const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
 
- Public Member Functions inherited from G4VCollision
 G4VCollision ()
 
void establish_G4MT_TLS_G4VCollision ()
 
 G4VCollision (void *s1, void *s2, void *s3, void *s4, void *s5, void *s6, void *s7)
 
virtual ~G4VCollision ()
 
G4bool operator== (const G4VCollision &right) const
 
G4bool operator!= (const G4VCollision &right) const
 
virtual G4double CrossSection (const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
 
virtual G4bool IsInCharge (const G4KineticTrack &trk1, const G4KineticTrack &trk2) const =0
 
virtual G4String GetName () const =0
 
virtual void Print () const
 
virtual void Print (const G4KineticTrack &trk1, const G4KineticTrack &trk2) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCollision
G4int GetNumberOfPartons (G4ParticleDefinition *aP) const
 
virtual const G4CollisionVectorGetComponents () const
 
virtual const
G4VCrossSectionSource
GetCrossSectionSource () const =0
 
virtual const
G4VAngularDistribution
GetAngularDistribution () const =0
 
virtual const std::vector
< G4String > & 
GetListOfColliders (G4int whichOne) const =0
 

Detailed Description

Definition at line 40 of file G4VElasticCollision.hh.

Constructor & Destructor Documentation

G4VElasticCollision::G4VElasticCollision ( )

Definition at line 45 of file G4VElasticCollision.cc.

46 {
47 }
G4VElasticCollision::~G4VElasticCollision ( )
virtual

Definition at line 50 of file G4VElasticCollision.cc.

51 { }

Member Function Documentation

G4KineticTrackVector * G4VElasticCollision::FinalState ( const G4KineticTrack trk1,
const G4KineticTrack trk2 
) const
virtual

Implements G4VCollision.

Definition at line 54 of file G4VElasticCollision.cc.

References CLHEP::HepLorentzVector::boostVector(), G4VAngularDistribution::CosTheta(), G4KineticTrack::Get4Momentum(), G4KineticTrack::GetActualMass(), G4VCollision::GetAngularDistribution(), G4KineticTrack::GetDefinition(), G4ParticleDefinition::GetPDGMass(), CLHEP::HepLorentzRotation::inverse(), CLHEP::Hep3Vector::mag2(), CLHEP::HepLorentzVector::mag2(), G4Neutron::Neutron(), G4VAngularDistribution::Phi(), G4Proton::Proton(), CLHEP::HepLorentzRotation::rotateY(), CLHEP::HepLorentzRotation::rotateZ(), and G4KineticTrack::Set4Momentum().

56 {
57  const G4VAngularDistribution* angDistribution;
58 
59  angDistribution = GetAngularDistribution();
60 
61 
62  G4LorentzVector pCM=trk1.Get4Momentum() + trk2.Get4Momentum();
63 
64  G4LorentzRotation toLabFrame(pCM.boostVector());
65  G4LorentzVector Ptmp=toLabFrame.inverse() * trk1.Get4Momentum(); //trk1 in CMS
67  toZ.rotateZ(-Ptmp.phi());
68  toZ.rotateY(-Ptmp.theta());
69  toLabFrame *= toZ.inverse();
70 
71  G4double S = pCM.mag2();
72  G4double m10 = trk1.GetDefinition()->GetPDGMass();
73  G4double m20 = trk2.GetDefinition()->GetPDGMass();
74  if(S-(m10+m20)*(m10+m20) < 0) return new G4KineticTrackVector;
75 
76  G4double m_1 = trk1.GetActualMass();
77  G4double m_2 = trk2.GetActualMass();
78 
79  // Angles of outgoing particles
80  G4double cosTheta = angDistribution->CosTheta(S,m_1,m_2);
81 
82  if ( (trk1.GetDefinition() == G4Proton::Proton() || trk1.GetDefinition() == G4Neutron::Neutron() )
83  &&(trk2.GetDefinition() == G4Proton::Proton() || trk2.GetDefinition() == G4Neutron::Neutron() ) )
84  {
85  if ( trk1.GetDefinition() == trk2.GetDefinition() )
86  {
87  if ( trk1.GetDefinition() == G4Proton::Proton() )
88  {
89 // G4cout << "scatterangle pp " << cosTheta
90 // << " " << typeid(*angDistribution).name() << G4endl;
91  } else {
92 // G4cout << "scatterangle nn " << cosTheta
93 // << " " << typeid(*angDistribution).name() << G4endl;
94  }
95  } else {
96 // G4cout << "scatterangle pn " << cosTheta
97 // << " " << typeid(*angDistribution).name() << G4endl;
98  }
99  } else {
100 // G4cout << "scatterangle other " << cosTheta
101 // << " " << typeid(*angDistribution).name() << G4endl;
102  }
103 
104  G4double phi = angDistribution->Phi();
105  G4double Theta = std::acos(cosTheta);
106 
107  // Unit vector of three-momentum
108  G4ThreeVector pFinal1(std::sin(Theta)*std::cos(phi), std::sin(Theta)*std::sin(phi), cosTheta);
109  // Three momentum in cm system
110  G4double pInCM = std::sqrt((S-(m10+m20)*(m10+m20))*(S-(m10-m20)*(m10-m20))/(4.*S));
111  pFinal1 = pFinal1 * pInCM;
112  G4ThreeVector pFinal2 = -pFinal1;
113 
114  G4double eFinal1 = std::sqrt(pFinal1.mag2() + m10*m10);
115  G4double eFinal2 = std::sqrt(pFinal2.mag2() + m20*m20);
116 
117  G4LorentzVector p4Final1(pFinal1, eFinal1);
118  G4LorentzVector p4Final2(pFinal2, eFinal2);
119 
120  // Lorentz transformation
121  p4Final1 *= toLabFrame;
122  p4Final2 *= toLabFrame;
123 
124  // Final tracks are copies of incoming ones, with modified 4-momenta
125  G4KineticTrack* final1 = new G4KineticTrack(trk1);
126  final1->Set4Momentum(p4Final1);
127  G4KineticTrack* final2 = new G4KineticTrack(trk2);
128  final2->Set4Momentum(p4Final2);
129 
130  G4KineticTrackVector* finalTracks = new G4KineticTrackVector;
131  finalTracks->push_back(final1);
132  finalTracks->push_back(final2);
133 
134  return finalTracks;
135 }
Hep3Vector boostVector() const
G4double GetActualMass() const
virtual const G4VAngularDistribution * GetAngularDistribution() const =0
virtual G4double CosTheta(G4double s, G4double m1, G4double m2) const =0
HepLorentzRotation & rotateY(double delta)
G4ParticleDefinition * GetDefinition() const
virtual G4double Phi() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void Set4Momentum(const G4LorentzVector &a4Momentum)
HepLorentzRotation & rotateZ(double delta)
G4double GetPDGMass() const
double mag2() const
double mag2() const
HepLorentzRotation inverse() const
double G4double
Definition: G4Types.hh:76
const G4LorentzVector & Get4Momentum() const
G4bool G4VElasticCollision::operator!= ( const G4VElasticCollision right) const
G4bool G4VElasticCollision::operator== ( const G4VElasticCollision right) const

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