Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VHadPhaseSpaceAlgorithm.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 uniform phase space generators.
29 // Subclasses implement a specific algorithm, such as Kopylov, GENBOD,
30 // or Makoto's NBody. Subclasses are used by G4HadDecayGenerator.
31 //
32 // Author: Michael Kelsey (SLAC) <kelsey@slac.stanford.edu>
33 
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 
48 
49 // Two body decay with uniform angular distribution
50 
53  const std::vector<G4double>& masses,
54  std::vector<G4LorentzVector>& finalState) {
55  if (GetVerboseLevel()>1)
56  G4cout << " >>> G4HadDecayGenerator::FillTwoBody" << G4endl;
57 
58  // Initialization and sanity check
59  finalState.clear();
60  if (masses.size() != 2U) return; // Should not have been called
61 
62  // Momentum of final state (energy balance has already been checked)
63  G4double p = TwoBodyMomentum(initialMass,masses[0],masses[1]);
64  if (GetVerboseLevel()>2) G4cout << " finalState momentum = " << p << G4endl;
65 
66  finalState.resize(2); // Allows filling by index
67  finalState[0].setVectM(UniformVector(p), masses[0]);
68  finalState[1].setVectM(-finalState[0].vect(), masses[1]);
69 }
70 
71 
72 // Samples a random vector with given magnitude
73 
75  // FIXME: Should this be made a static thread-local buffer?
76  G4ThreeVector v;
78  return v;
79 }
G4double UniformPhi() const
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
G4int GetVerboseLevel() const
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const
G4ThreeVector UniformVector(G4double mag=1.) const
double G4double
Definition: G4Types.hh:76
virtual void GenerateTwoBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
void setRThetaPhi(double r, double theta, double phi)
G4GLOB_DLL std::ostream G4cout
G4double UniformTheta() const