Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VHadPhaseSpaceAlgorithm.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// Abstract base class for multibody uniform phase space generators.
28// Subclasses implement a specific algorithm, such as Kopylov, GENBOD,
29// or Makoto's NBody. Subclasses are used by G4HadDecayGenerator.
30//
31// Author: Michael Kelsey (SLAC) <kelsey@slac.stanford.edu>
32
36#include "G4SystemOfUnits.hh"
37#include "G4ThreeVector.hh"
38#include "Randomize.hh"
39#include <algorithm>
40#include <iostream>
41#include <iterator>
42#include <numeric>
43#include <vector>
44
45
46
47
48// Two body decay with uniform angular distribution
49
51GenerateTwoBody(G4double initialMass,
52 const std::vector<G4double>& masses,
53 std::vector<G4LorentzVector>& finalState) {
54 if (GetVerboseLevel()>1)
55 G4cout << " >>> G4HadDecayGenerator::FillTwoBody" << G4endl;
56
57 // Initialization and sanity check
58 finalState.clear();
59 if (masses.size() != 2U) return; // Should not have been called
60
61 // Momentum of final state (energy balance has already been checked)
62 G4double p = TwoBodyMomentum(initialMass,masses[0],masses[1]);
63 if (GetVerboseLevel()>2) G4cout << " finalState momentum = " << p << G4endl;
64
65 finalState.resize(2); // Allows filling by index
66 finalState[0].setVectM(UniformVector(p), masses[0]);
67 finalState[1].setVectM(-finalState[0].vect(), masses[1]);
68}
69
70
71// Samples a random vector with given magnitude
72
74 // FIXME: Should this be made a static thread-local buffer?
77 return v;
78}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void setRThetaPhi(double r, double theta, double phi)
G4double UniformTheta() const
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const
virtual void GenerateTwoBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
G4ThreeVector UniformVector(G4double mag=1.) const