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, which provides multiple algorithms 28 // for sampling. Momentum vectors are generated in the center-of-mass 29 // frame of the decay, and returned in a user-supplied buffer. A sampling 30 // algorithm is specified via constructor argument. 31 // 32 // Author: Michael Kelsey (SLAC) <kelsey@slac.stanford.edu> 33 34 #include "G4HadDecayGenerator.hh" 35 #include "G4VHadDecayAlgorithm.hh" 36 #include "G4HadPhaseSpaceKopylov.hh" 37 #include "G4HadPhaseSpaceGenbod.hh" 38 #include "G4HadPhaseSpaceNBodyAsai.hh" 39 #include "G4HadronicException.hh" 40 #include "G4LorentzVector.hh" 41 #include "G4ParticleDefinition.hh" 42 #include "G4PhysicalConstants.hh" 43 #include "G4SystemOfUnits.hh" 44 #include "G4ThreeVector.hh" 45 #include "Randomize.hh" 46 #include <vector> 47 #include <algorithm> 48 #include <numeric> 49 #include <iterator> 50 #include <iostream> 51 52 53 // Constructors and destructor 54 55 G4HadDecayGenerator::G4HadDecayGenerator(Algorithm alg, G4int verbose) 56 : verboseLevel(verbose), theAlgorithm(0) { 57 switch (alg) { 58 case Kopylov: theAlgorithm = new G4HadPhaseSpaceKopylov(verboseLevel); break; 59 case GENBOD: theAlgorithm = new G4HadPhaseSpaceGenbod(verboseLevel); break; 60 case NBody: theAlgorithm = new G4HadPhaseSpaceNBodyAsai(verboseLevel); break; 61 case NONE: theAlgorithm = 0; break; // User may explicitly set no algorithm 62 default: ReportInvalidAlgorithm(alg); 63 } 64 65 if (verboseLevel) { 66 G4cout << " >>> G4HadDecayGenerator"; 67 if (theAlgorithm) G4cout << " using " << theAlgorithm->GetName(); 68 G4cout << G4endl; 69 } 70 } 71 72 G4HadDecayGenerator::G4HadDecayGenerator(G4VHadDecayAlgorithm* alg, 73 G4int verbose) 74 : verboseLevel(verbose), theAlgorithm(alg) { 75 if (verboseLevel) { 76 G4cout << " >>> G4HadDecayGenerator"; 77 if (theAlgorithm) G4cout << " using " << theAlgorithm->GetName(); 78 G4cout << G4endl; 79 } 80 } 81 82 G4HadDecayGenerator::~G4HadDecayGenerator() { 83 delete theAlgorithm; 84 theAlgorithm = 0; 85 } 86 87 88 // Sanity checks -- throws exception if no algorithm chosen 89 90 void G4HadDecayGenerator::ReportInvalidAlgorithm(Algorithm alg) const { 91 if (verboseLevel) 92 G4cerr << "G4HadDecayGenerator: bad algorithm code " << alg << G4endl; 93 94 throw G4HadronicException(__FILE__, __LINE__, "Invalid algorithm code"); 95 } 96 97 void G4HadDecayGenerator::ReportMissingAlgorithm() const { 98 if (verboseLevel) 99 G4cerr << "G4HadDecayGenerator: no algorithm specified" << G4endl; 100 101 throw G4HadronicException(__FILE__, __LINE__, "Null algorithm pointer"); 102 } 103 104 105 // Enable (or disable if 0) diagnostic messages 106 void G4HadDecayGenerator::SetVerboseLevel(G4int verbose) { 107 verboseLevel = verbose; 108 if (theAlgorithm) theAlgorithm->SetVerboseLevel(verbose); 109 } 110 111 const G4String& G4HadDecayGenerator::GetAlgorithmName() const { 112 static const G4String& none = "NONE"; 113 return (theAlgorithm ? theAlgorithm->GetName() : none); 114 } 115 116 117 // Initial state (rest mass) and list of final masses 118 119 G4bool 120 G4HadDecayGenerator::Generate(G4double initialMass, 121 const std::vector<G4double>& masses, 122 std::vector<G4LorentzVector>& finalState) { 123 if (verboseLevel) 124 G4cout << " >>> G4HadDecayGenerator::Generate (mass)" << G4endl; 125 126 if (!theAlgorithm) ReportMissingAlgorithm(); 127 128 if (masses.size() == 1U) 129 return GenerateOneBody(initialMass, masses, finalState); 130 131 theAlgorithm->Generate(initialMass, masses, finalState); 132 return !finalState.empty(); // Generator failure returns empty state 133 } 134 135 // Initial state particle and list of final masses 136 137 G4bool 138 G4HadDecayGenerator::Generate(const G4ParticleDefinition* initialPD, 139 const std::vector<G4double>& masses, 140 std::vector<G4LorentzVector>& finalState) { 141 if (verboseLevel) 142 G4cout << " >>> G4HadDecayGenerator::Generate (particle)" << G4endl; 143 144 return (initialPD && Generate(initialPD->GetPDGMass(), masses, finalState)); 145 } 146 147 // Final state particles will be boosted to initial-state frame 148 149 G4bool 150 G4HadDecayGenerator::Generate(const G4LorentzVector& initialState, 151 const std::vector<G4double>& masses, 152 std::vector<G4LorentzVector>& finalState) { 153 if (verboseLevel) 154 G4cout << " >>> G4HadDecayGenerator::Generate (frame)" << G4endl; 155 156 G4bool good = Generate(initialState.m(), masses, finalState); 157 if (good) { 158 G4ThreeVector bv = initialState.boostVector(); 159 for (size_t i=0; i<finalState.size(); i++) { 160 finalState[i].boost(bv); 161 } 162 } 163 164 return good; 165 } 166 167 168 // Handle special case of "one body decay" (used for kaon mixing) 169 170 G4bool G4HadDecayGenerator:: 171 GenerateOneBody(G4double initialMass, 172 const std::vector<G4double>& masses, 173 std::vector<G4LorentzVector>& finalState) const { 174 if (verboseLevel>1) 175 G4cout << " >>> G4HadDecayGenerator::GenerateOneBody" << G4endl; 176 177 // Initialization and sanity checks 178 finalState.clear(); 179 180 if (masses.size() != 1U) return false; // Should not have been called 181 if (std::fabs(initialMass-masses[0]) > eV) return false; 182 183 if (verboseLevel>2) G4cout << " finalState mass = " << masses[0] << G4endl; 184 185 finalState.push_back(G4LorentzVector(0.,0.,0.,masses[0])); 186 return true; 187 } 188