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