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 // Multibody "phase space" generator using Kopylov's algorithm 28 // 29 // Author: Michael Kelsey (SLAC) <kelsey@slac.stanford.edu> 30 31 #include "G4HadPhaseSpaceKopylov.hh" 32 #include "G4LorentzVector.hh" 33 #include "G4Pow.hh" 34 #include "Randomize.hh" 35 #include <vector> 36 #include <algorithm> 37 #include <numeric> 38 #include <cmath> 39 40 41 // Generator 42 43 void G4HadPhaseSpaceKopylov:: 44 GenerateMultiBody(G4double initialMass, 45 const std::vector<G4double>& masses, 46 std::vector<G4LorentzVector>& finalState) { 47 if (GetVerboseLevel()) G4cout << GetName() << "::GenerateMultiBody" << G4endl; 48 49 finalState.clear(); 50 51 G4int N = (G4int)masses.size(); 52 finalState.resize(N); 53 54 G4double mtot = std::accumulate(masses.begin(), masses.end(), 0.0); 55 G4double mu = mtot; 56 G4double PFragMagCM = 0.0; 57 G4double Mass = initialMass; 58 G4double T = Mass-mtot; 59 G4LorentzVector PFragCM(0.0,0.0,0.0,0.0); 60 G4LorentzVector PRestCM(0.0,0.0,0.0,0.0); 61 G4LorentzVector PRestLab(0.0,0.0,0.0,Mass); 62 63 for (G4int k=N-1; k>0; --k) { 64 mu -= masses[k]; 65 T *= (k>1) ? BetaKopylov(k) : 0.; 66 67 G4double RestMass = mu + T; 68 69 PFragMagCM = TwoBodyMomentum(Mass,masses[k],RestMass); 70 71 // Create a unit vector with a random direction isotropically distributed 72 G4ThreeVector RandVector = UniformVector(PFragMagCM); 73 74 PFragCM.setVectM(RandVector,masses[k]); 75 PRestCM.setVectM(-RandVector,RestMass); 76 77 G4ThreeVector BoostV = PRestLab.boostVector(); 78 79 PFragCM.boost(BoostV); 80 PRestCM.boost(BoostV); 81 PRestLab = PRestCM; 82 Mass = RestMass; 83 finalState[k] = PFragCM; 84 } 85 86 finalState[0] = PRestLab; 87 } 88 89 90 // Generate scale factor for final state particle 91 92 G4double G4HadPhaseSpaceKopylov::BetaKopylov(G4int K) const { 93 G4Pow* g4pow = G4Pow::GetInstance(); 94 95 G4int N = 3*K - 5; 96 G4double xN = G4double(N); 97 G4double Fmax = std::sqrt(g4pow->powN(xN/(xN+1.),N)/(xN+1.)); 98 99 G4double F, chi; 100 const G4int maxNumberOfLoops = 10000; 101 G4int loopCounter = 0; 102 do { 103 chi = G4UniformRand(); 104 F = std::sqrt(g4pow->powN(chi,N)*(1.-chi)); 105 } while ( ( Fmax*G4UniformRand() > F ) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 02.11.2015, A.Ribon */ 106 if ( loopCounter >= maxNumberOfLoops ) { 107 G4ExceptionDescription ed; 108 ed << " Failed sampling after maxNumberOfLoops attempts : forced exit" << G4endl; 109 G4Exception( " G4HadPhaseSpaceKopylov::BetaKopylov ", "HAD_KOPYLOV_001", JustWarning, ed ); 110 } 111 112 return chi; 113 } 114