Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VHadDecayAlgorithm.cc
이 파일의 문서화 페이지로 가기
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 // Abstract base class for multibody "phase space" generators. Subclasses
29 // implement a specific algorithm, such as Kopylov, GENBOD, or Makoto's
30 // NBody. Subclasses are used by G4HadPhaseSpaceGenerator.
31 //
32 // Author: Michael Kelsey (SLAC) <kelsey@slac.stanford.edu>
33 
34 #include "G4VHadDecayAlgorithm.hh"
35 #include "G4HadronicException.hh"
36 #include "G4PhysicalConstants.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4ThreeVector.hh"
39 #include "Randomize.hh"
40 #include <algorithm>
41 #include <iostream>
42 #include <iterator>
43 #include <numeric>
44 #include <vector>
45 
46 
47 // Initial state (rest mass) and list of final masses
48 
50  const std::vector<G4double>& masses,
51  std::vector<G4LorentzVector>& finalState) {
52  if (verboseLevel) G4cout << GetName() << "::Generate" << G4endl;
53 
54  // Initialization and sanity check
55  finalState.clear();
56  if (!IsDecayAllowed(initialMass, masses)) return;
57 
58  // Allow different procedures for two-body or N-body distributions
59  if (masses.size() == 2U)
60  GenerateTwoBody(initialMass, masses, finalState);
61  else
62  GenerateMultiBody(initialMass, masses, finalState);
63 }
64 
65 
66 // Base class does very simple validation of configuration
67 
70  const std::vector<G4double>& masses) const {
71  G4bool okay =
72  (initialMass > 0. && masses.size() >= 2 &&
73  initialMass >= std::accumulate(masses.begin(),masses.end(),0.));
74 
75  if (verboseLevel) {
76  G4cout << GetName() << "::IsDecayAllowed? initialMass " << initialMass
77  << " " << masses.size() << " masses sum "
78  << std::accumulate(masses.begin(),masses.end(),0.) << G4endl;
79 
80  if (verboseLevel>1) PrintVector(masses," ",G4cout);
81 
82  G4cout << " Returning " << okay << G4endl;
83  }
84 
85  return okay;
86 }
87 
88 
89 // Momentum function (c.f. PDK() function from CERNLIB W515)
90 
92  G4double M2) const {
93  G4double PSQ = (M0+M1+M2)*(M0+M1-M2)*(M0-M1+M2)*(M0-M1-M2);
94  if (PSQ < 0.) {
95  G4cout << GetName() << ": problem of decay of M(GeV) " << M0/GeV
96  << " to M1(GeV) " << M1/GeV << " and M2(GeV) " << M2/GeV
97  << " PSQ(MeV) " << PSQ/MeV << " < 0" << G4endl;
98  // exception only if the problem is numerically significant
99  if (PSQ < -CLHEP::eV) {
100  throw G4HadronicException(__FILE__, __LINE__,"Error in decay kinematics");
101  }
102 
103  PSQ = 0.;
104  }
105 
106  return std::sqrt(PSQ)/(2.*M0);
107 }
108 
109 // Convenience functions for uniform angular distributions
110 
112  return std::acos(2.0*G4UniformRand() - 1.0);
113 }
114 
116  return twopi*G4UniformRand();
117 }
118 
119 
120 // Dump contents of vector to output
121 
123 PrintVector(const std::vector<G4double>& v,
124  const G4String& vname, std::ostream& os) const {
125  os << " " << vname << "(" << v.size() << ") ";
126  std::copy(v.begin(), v.end(), std::ostream_iterator<G4double>(os, " "));
127  os << std::endl;
128 }
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double UniformPhi() const
#define G4endl
Definition: G4ios.hh:61
virtual G4bool IsDecayAllowed(G4double initialMass, const std::vector< G4double > &masses) const
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const
virtual void GenerateTwoBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0
virtual void GenerateMultiBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0
void PrintVector(const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double eV
G4GLOB_DLL std::ostream G4cout
void Generate(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
static constexpr double GeV
Definition: G4SIunits.hh:217
const G4String & GetName() const
G4double UniformTheta() const