Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VElasticCollision.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
28#include "globals.hh"
30#include "G4KineticTrack.hh"
32#include "G4Proton.hh"
33#include "G4Neutron.hh"
34#include "G4XNNElastic.hh"
36#include "G4ThreeVector.hh"
37#include "G4LorentzVector.hh"
38#include "G4LorentzRotation.hh"
40#include "G4AngularDistributionNP.hh" // np scattering
41#include "G4AngularDistributionPP.hh" // nn and pp scattering
42#include <typeinfo>
43
45{
46}
47
48
50{ }
51
52
54 const G4KineticTrack& trk2) const
55{
56 const G4VAngularDistribution* angDistribution;
57
58 angDistribution = GetAngularDistribution();
59
60
61 G4LorentzVector pCM=trk1.Get4Momentum() + trk2.Get4Momentum();
62
63 G4LorentzRotation toLabFrame(pCM.boostVector());
64 G4LorentzVector Ptmp=toLabFrame.inverse() * trk1.Get4Momentum(); //trk1 in CMS
66 toZ.rotateZ(-Ptmp.phi());
67 toZ.rotateY(-Ptmp.theta());
68 toLabFrame *= toZ.inverse();
69
70 G4double S = pCM.mag2();
71 G4double m10 = trk1.GetDefinition()->GetPDGMass();
72 G4double m20 = trk2.GetDefinition()->GetPDGMass();
73 if(S-(m10+m20)*(m10+m20) < 0) return new G4KineticTrackVector;
74
75 G4double m_1 = trk1.GetActualMass();
76 G4double m_2 = trk2.GetActualMass();
77
78 // Angles of outgoing particles
79 G4double cosTheta = angDistribution->CosTheta(S,m_1,m_2);
80
81 if ( (trk1.GetDefinition() == G4Proton::Proton() || trk1.GetDefinition() == G4Neutron::Neutron() )
82 &&(trk2.GetDefinition() == G4Proton::Proton() || trk2.GetDefinition() == G4Neutron::Neutron() ) )
83 {
84 if ( trk1.GetDefinition() == trk2.GetDefinition() )
85 {
86 if ( trk1.GetDefinition() == G4Proton::Proton() )
87 {
88// G4cout << "scatterangle pp " << cosTheta
89// << " " << typeid(*angDistribution).name() << G4endl;
90 } else {
91// G4cout << "scatterangle nn " << cosTheta
92// << " " << typeid(*angDistribution).name() << G4endl;
93 }
94 } else {
95// G4cout << "scatterangle pn " << cosTheta
96// << " " << typeid(*angDistribution).name() << G4endl;
97 }
98 } else {
99// G4cout << "scatterangle other " << cosTheta
100// << " " << typeid(*angDistribution).name() << G4endl;
101 }
102
103 G4double phi = angDistribution->Phi();
104 G4double Theta = std::acos(cosTheta);
105
106 // Unit vector of three-momentum
107 G4ThreeVector pFinal1(std::sin(Theta)*std::cos(phi), std::sin(Theta)*std::sin(phi), cosTheta);
108 // Three momentum in cm system
109 G4double pInCM = std::sqrt((S-(m10+m20)*(m10+m20))*(S-(m10-m20)*(m10-m20))/(4.*S));
110 pFinal1 = pFinal1 * pInCM;
111 G4ThreeVector pFinal2 = -pFinal1;
112
113 G4double eFinal1 = std::sqrt(pFinal1.mag2() + m10*m10);
114 G4double eFinal2 = std::sqrt(pFinal2.mag2() + m20*m20);
115
116 G4LorentzVector p4Final1(pFinal1, eFinal1);
117 G4LorentzVector p4Final2(pFinal2, eFinal2);
118
119 // Lorentz transformation
120 p4Final1 *= toLabFrame;
121 p4Final2 *= toLabFrame;
122
123 // Final tracks are copies of incoming ones, with modified 4-momenta
124 G4KineticTrack* final1 = new G4KineticTrack(trk1);
125 final1->Set4Momentum(p4Final1);
126 G4KineticTrack* final2 = new G4KineticTrack(trk2);
127 final2->Set4Momentum(p4Final2);
128
130 finalTracks->push_back(final1);
131 finalTracks->push_back(final2);
132
133 return finalTracks;
134}
G4double S(G4double temp)
double G4double
Definition: G4Types.hh:83
double mag2() const
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
double theta() const
Hep3Vector boostVector() const
void Set4Momentum(const G4LorentzVector &a4Momentum)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetActualMass() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92
virtual G4double Phi() const
virtual G4double CosTheta(G4double s, G4double m1, G4double m2) const =0
virtual const G4VAngularDistribution * GetAngularDistribution() const =0
virtual G4KineticTrackVector * FinalState(const G4KineticTrack &trk1, const G4KineticTrack &trk2) const