Geant4 Cross Reference |
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 "phase space" generators. Subclasses 28 // implement a specific algorithm, such as Kopylov, GENBOD, or Makoto's 29 // NBody. Subclasses are used by G4HadPhaseSpaceGenerator. 30 // 31 // Author: Michael Kelsey (SLAC) <kelsey@slac.stanford.edu> 32 33 #include "G4VHadDecayAlgorithm.hh" 34 #include "G4HadronicException.hh" 35 #include "G4PhysicalConstants.hh" 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 // Initial state (rest mass) and list of final masses 47 48 void G4VHadDecayAlgorithm::Generate(G4double initialMass, 49 const std::vector<G4double>& masses, 50 std::vector<G4LorentzVector>& finalState) { 51 if (verboseLevel) G4cout << GetName() << "::Generate" << G4endl; 52 53 // Initialization and sanity check 54 finalState.clear(); 55 if (!IsDecayAllowed(initialMass, masses)) return; 56 57 // Allow different procedures for two-body or N-body distributions 58 if (masses.size() == 2U) 59 GenerateTwoBody(initialMass, masses, finalState); 60 else 61 GenerateMultiBody(initialMass, masses, finalState); 62 } 63 64 65 // Base class does very simple validation of configuration 66 67 G4bool G4VHadDecayAlgorithm:: 68 IsDecayAllowed(G4double initialMass, 69 const std::vector<G4double>& masses) const { 70 G4bool okay = 71 (initialMass > 0. && masses.size() >= 2 && 72 initialMass >= std::accumulate(masses.begin(),masses.end(),0.)); 73 74 if (verboseLevel) { 75 G4cout << GetName() << "::IsDecayAllowed? initialMass " << initialMass 76 << " " << masses.size() << " masses sum " 77 << std::accumulate(masses.begin(),masses.end(),0.) << G4endl; 78 79 if (verboseLevel>1) PrintVector(masses," ",G4cout); 80 81 G4cout << " Returning " << okay << G4endl; 82 } 83 84 return okay; 85 } 86 87 88 // Momentum function (c.f. PDK() function from CERNLIB W515) 89 90 G4double G4VHadDecayAlgorithm::TwoBodyMomentum(G4double M0, G4double M1, 91 G4double M2) const { 92 G4double PSQ = (M0+M1+M2)*(M0+M1-M2)*(M0-M1+M2)*(M0-M1-M2); 93 if (PSQ < 0.) { 94 G4cout << GetName() << ": problem of decay of M(GeV) " << M0/GeV 95 << " to M1(GeV) " << M1/GeV << " and M2(GeV) " << M2/GeV 96 << " PSQ(MeV) " << PSQ/MeV << " < 0" << G4endl; 97 // exception only if the problem is numerically significant 98 if (PSQ < -CLHEP::eV) { 99 throw G4HadronicException(__FILE__, __LINE__,"Error in decay kinematics"); 100 } 101 102 PSQ = 0.; 103 } 104 105 return std::sqrt(PSQ)/(2.*M0); 106 } 107 108 // Convenience functions for uniform angular distributions 109 110 G4double G4VHadDecayAlgorithm::UniformTheta() const { 111 return std::acos(2.0*G4UniformRand() - 1.0); 112 } 113 114 G4double G4VHadDecayAlgorithm::UniformPhi() const { 115 return twopi*G4UniformRand(); 116 } 117 118 119 // Dump contents of vector to output 120 121 void G4VHadDecayAlgorithm:: 122 PrintVector(const std::vector<G4double>& v, 123 const G4String& vname, std::ostream& os) const { 124 os << " " << vname << "(" << v.size() << ") "; 125 std::copy(v.begin(), v.end(), std::ostream_iterator<G4double>(os, " ")); 126 os << std::endl; 127 } 128