Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ReactionKinematics.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 // CERN Geneva Switzerland
29 //
30 // ------------ G4ReactionDynamics::TwoBody ------
31 // new method TwoBodyScattering
32 // by Christian V"olcker (CERN-Munich), September 1997
33 // E-mail: Christian.Volcker@cern.ch
34 // ************************************************************
35 //-----------------------------------------------------------------------------
36 
37 #include "G4ReactionKinematics.hh"
38 #include "G4PhysicalConstants.hh"
39 
40 // ************************************************************
42  const G4DynamicParticle* pIn1, const G4DynamicParticle* pIn2,
43  G4DynamicParticle* pOut1, G4DynamicParticle* pOut2)
44 // ************************************************************
45 {
46 // initial particles:
47 
48 // - total invariant mass
49  G4LorentzVector sumIn(pIn1->Get4Momentum()+pIn2->Get4Momentum());
50  G4double invariantMass=sumIn.mag();
51 
52 // - beta of center-of-mass system
53  G4ThreeVector betaCMS=sumIn.boostVector();
54 
55 // final particles:
56 
57 // - get final particle masses
58  G4double massOut1=pOut1->GetMass();
59  G4double massOut2=pOut2->GetMass();
60 
61 // - calculate breakup momentum:
62  G4double breakupMomentum=BreakupMomentum(invariantMass, massOut1, massOut2);
63 
64 // - isotropic decay angle
65  G4double costheta = 2.0*G4UniformRand() - 1.0;
66  G4double sintheta = std::sqrt(1.0 - costheta*costheta);
67  G4double phi = 2.0*pi*G4UniformRand();
68 
69 // - setup LorentzVectors
70  G4double pz=costheta*breakupMomentum;
71  G4double px=sintheta*std::cos(phi)*breakupMomentum;
72  G4double py=sintheta*std::sin(phi)*breakupMomentum;
73 
74  G4double breakupMomentumSquared=breakupMomentum*breakupMomentum;
75  G4double energy1=std::sqrt(breakupMomentumSquared+massOut1*massOut1);
76  G4double energy2=std::sqrt(breakupMomentumSquared+massOut2*massOut2);
77 
78  G4LorentzVector lorentz1(px, py, pz, energy1);
79  G4LorentzVector lorentz2(px, py, pz, energy2);
80 
81 // - back into lab system
82 
83  lorentz1.boost(betaCMS);
84  lorentz2.boost(betaCMS);
85 
86 // fill in new particles:
87 
88  pOut1->Set4Momentum(lorentz1);
89  pOut2->Set4Momentum(lorentz2);
90 
91  return;
92 }
#define G4UniformRand()
Definition: Randomize.hh:87
G4double GetMass() const
HepLorentzVector & boost(double, double, double)
G4double invariantMass(const G4double E, const ThreeVector &p)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double BreakupMomentum(G4double totalMass, G4double m1, G4double m2)
void TwoBodyScattering(const G4DynamicParticle *pIn1, const G4DynamicParticle *pIn2, G4DynamicParticle *pOut1, G4DynamicParticle *pOut2)
double G4double
Definition: G4Types.hh:76