Geant4-11
Public Member Functions | Protected Member Functions | Private Attributes
G4HadPhaseSpaceKopylov Class Reference

#include <G4HadPhaseSpaceKopylov.hh>

Inheritance diagram for G4HadPhaseSpaceKopylov:
G4VHadPhaseSpaceAlgorithm G4VHadDecayAlgorithm

Public Member Functions

 G4HadPhaseSpaceKopylov (G4int verbose=0)
 
void Generate (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
const G4StringGetName () const
 
G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int verbose)
 
virtual ~G4HadPhaseSpaceKopylov ()
 

Protected Member Functions

G4double BetaKopylov (G4int K) const
 
virtual void GenerateMultiBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
virtual void GenerateTwoBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
virtual G4bool IsDecayAllowed (G4double initialMass, const std::vector< G4double > &masses) const
 
void PrintVector (const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
 
G4double TwoBodyMomentum (G4double M0, G4double M1, G4double M2) const
 
G4double UniformPhi () const
 
G4double UniformTheta () const
 
G4ThreeVector UniformVector (G4double mag=1.) const
 

Private Attributes

G4String name
 
G4int verboseLevel
 

Detailed Description

Definition at line 37 of file G4HadPhaseSpaceKopylov.hh.

Constructor & Destructor Documentation

◆ G4HadPhaseSpaceKopylov()

G4HadPhaseSpaceKopylov::G4HadPhaseSpaceKopylov ( G4int  verbose = 0)
inline

Definition at line 39 of file G4HadPhaseSpaceKopylov.hh.

40 : G4VHadPhaseSpaceAlgorithm("G4HadPhaseSpaceKopylov",verbose) {;}
G4VHadPhaseSpaceAlgorithm(const G4String &algName, G4int verbose=0)

◆ ~G4HadPhaseSpaceKopylov()

virtual G4HadPhaseSpaceKopylov::~G4HadPhaseSpaceKopylov ( )
inlinevirtual

Definition at line 41 of file G4HadPhaseSpaceKopylov.hh.

41{;}

Member Function Documentation

◆ BetaKopylov()

G4double G4HadPhaseSpaceKopylov::BetaKopylov ( G4int  K) const
protected

Definition at line 92 of file G4HadPhaseSpaceKopylov.cc.

92 {
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 ) {
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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
#define G4UniformRand()
Definition: Randomize.hh:52
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:166

References G4endl, G4Exception(), G4UniformRand, G4Pow::GetInstance(), JustWarning, and G4Pow::powN().

Referenced by GenerateMultiBody().

◆ Generate()

void G4VHadDecayAlgorithm::Generate ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
inherited

Definition at line 48 of file G4VHadDecayAlgorithm.cc.

50 {
51 if (verboseLevel) G4cout << GetName() << "::Generate" << G4endl;
52
53 // Initialization and sanity check
54 finalState.clear();
55 if (!IsDecayAllowed(initialMass, masses)) return;
56
57 // Allow different procedures for two-body or N-body distributions
58 if (masses.size() == 2U)
59 GenerateTwoBody(initialMass, masses, finalState);
60 else
61 GenerateMultiBody(initialMass, masses, finalState);
62}
G4GLOB_DLL std::ostream G4cout
virtual void GenerateTwoBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0
const G4String & GetName() const
virtual G4bool IsDecayAllowed(G4double initialMass, const std::vector< G4double > &masses) const
virtual void GenerateMultiBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0

References G4cout, G4endl, G4VHadDecayAlgorithm::GenerateMultiBody(), G4VHadDecayAlgorithm::GenerateTwoBody(), G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::IsDecayAllowed(), and G4VHadDecayAlgorithm::verboseLevel.

Referenced by G4HadDecayGenerator::Generate().

◆ GenerateMultiBody()

void G4HadPhaseSpaceKopylov::GenerateMultiBody ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protectedvirtual

Implements G4VHadDecayAlgorithm.

Definition at line 43 of file G4HadPhaseSpaceKopylov.cc.

46 {
47 if (GetVerboseLevel()) G4cout << GetName() << "::GenerateMultiBody" << G4endl;
48
49 finalState.clear();
50
51 size_t N = 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 (size_t 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}
G4double BetaKopylov(G4int K) const
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const
G4ThreeVector UniformVector(G4double mag=1.) const

References BetaKopylov(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), CLHEP::HepLorentzVector::setVectM(), G4VHadDecayAlgorithm::TwoBodyMomentum(), and G4VHadPhaseSpaceAlgorithm::UniformVector().

◆ GenerateTwoBody()

void G4VHadPhaseSpaceAlgorithm::GenerateTwoBody ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protectedvirtualinherited

Implements G4VHadDecayAlgorithm.

Definition at line 50 of file G4VHadPhaseSpaceAlgorithm.cc.

53 {
54 if (GetVerboseLevel()>1)
55 G4cout << " >>> G4HadDecayGenerator::FillTwoBody" << G4endl;
56
57 // Initialization and sanity check
58 finalState.clear();
59 if (masses.size() != 2U) return; // Should not have been called
60
61 // Momentum of final state (energy balance has already been checked)
62 G4double p = TwoBodyMomentum(initialMass,masses[0],masses[1]);
63 if (GetVerboseLevel()>2) G4cout << " finalState momentum = " << p << G4endl;
64
65 finalState.resize(2); // Allows filling by index
66 finalState[0].setVectM(UniformVector(p), masses[0]);
67 finalState[1].setVectM(-finalState[0].vect(), masses[1]);
68}

References G4cout, G4endl, G4VHadDecayAlgorithm::GetVerboseLevel(), G4VHadDecayAlgorithm::TwoBodyMomentum(), and G4VHadPhaseSpaceAlgorithm::UniformVector().

◆ GetName()

const G4String & G4VHadDecayAlgorithm::GetName ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4VHadDecayAlgorithm::GetVerboseLevel ( ) const
inlineinherited

◆ IsDecayAllowed()

G4bool G4VHadDecayAlgorithm::IsDecayAllowed ( G4double  initialMass,
const std::vector< G4double > &  masses 
) const
protectedvirtualinherited

Definition at line 67 of file G4VHadDecayAlgorithm.cc.

69 {
70 G4bool okay =
71 (initialMass > 0. && masses.size() >= 2 &&
72 initialMass >= std::accumulate(masses.begin(),masses.end(),0.));
73
74 if (verboseLevel) {
75 G4cout << GetName() << "::IsDecayAllowed? initialMass " << initialMass
76 << " " << masses.size() << " masses sum "
77 << std::accumulate(masses.begin(),masses.end(),0.) << G4endl;
78
79 if (verboseLevel>1) PrintVector(masses," ",G4cout);
80
81 G4cout << " Returning " << okay << G4endl;
82 }
83
84 return okay;
85}
bool G4bool
Definition: G4Types.hh:86
void PrintVector(const std::vector< G4double > &v, const G4String &name, std::ostream &os) const

References G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::PrintVector(), and G4VHadDecayAlgorithm::verboseLevel.

Referenced by G4VHadDecayAlgorithm::Generate().

◆ PrintVector()

void G4VHadDecayAlgorithm::PrintVector ( const std::vector< G4double > &  v,
const G4String name,
std::ostream &  os 
) const
protectedinherited

Definition at line 121 of file G4VHadDecayAlgorithm.cc.

123 {
124 os << " " << vname << "(" << v.size() << ") ";
125 std::copy(v.begin(), v.end(), std::ostream_iterator<G4double>(os, " "));
126 os << std::endl;
127}
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

References field_utils::copy().

Referenced by G4HadPhaseSpaceGenbod::FillEnergySteps(), G4HadPhaseSpaceGenbod::FillRandomBuffer(), G4HadPhaseSpaceNBodyAsai::GenerateMultiBody(), G4HadPhaseSpaceGenbod::Initialize(), and G4VHadDecayAlgorithm::IsDecayAllowed().

◆ SetVerboseLevel()

virtual void G4VHadDecayAlgorithm::SetVerboseLevel ( G4int  verbose)
inlinevirtualinherited

◆ TwoBodyMomentum()

G4double G4VHadDecayAlgorithm::TwoBodyMomentum ( G4double  M0,
G4double  M1,
G4double  M2 
) const
protectedinherited

Definition at line 90 of file G4VHadDecayAlgorithm.cc.

91 {
92 G4double PSQ = (M0+M1+M2)*(M0+M1-M2)*(M0-M1+M2)*(M0-M1-M2);
93 if (PSQ < 0.) {
94 G4cout << GetName() << ": problem of decay of M(GeV) " << M0/GeV
95 << " to M1(GeV) " << M1/GeV << " and M2(GeV) " << M2/GeV
96 << " PSQ(MeV) " << PSQ/MeV << " < 0" << G4endl;
97 // exception only if the problem is numerically significant
98 if (PSQ < -CLHEP::eV) {
99 throw G4HadronicException(__FILE__, __LINE__,"Error in decay kinematics");
100 }
101
102 PSQ = 0.;
103 }
104
105 return std::sqrt(PSQ)/(2.*M0);
106}
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double eV

References CLHEP::eV, G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), GeV, and MeV.

Referenced by G4HadPhaseSpaceGenbod::ComputeWeightScale(), G4HadPhaseSpaceGenbod::FillEnergySteps(), G4CascadeFinalStateAlgorithm::FillUsingKopylov(), GenerateMultiBody(), G4HadPhaseSpaceNBodyAsai::GenerateMultiBody(), G4CascadeFinalStateAlgorithm::GenerateTwoBody(), and G4VHadPhaseSpaceAlgorithm::GenerateTwoBody().

◆ UniformPhi()

G4double G4VHadDecayAlgorithm::UniformPhi ( ) const
protectedinherited

◆ UniformTheta()

G4double G4VHadDecayAlgorithm::UniformTheta ( ) const
protectedinherited

Definition at line 110 of file G4VHadDecayAlgorithm.cc.

110 {
111 return std::acos(2.0*G4UniformRand() - 1.0);
112}

References G4UniformRand.

Referenced by G4CascadeFinalStateAlgorithm::FillUsingKopylov(), and G4VHadPhaseSpaceAlgorithm::UniformVector().

◆ UniformVector()

G4ThreeVector G4VHadPhaseSpaceAlgorithm::UniformVector ( G4double  mag = 1.) const
protectedinherited

Definition at line 73 of file G4VHadPhaseSpaceAlgorithm.cc.

73 {
74 // FIXME: Should this be made a static thread-local buffer?
77 return v;
78}
void setRThetaPhi(double r, double theta, double phi)
G4double UniformTheta() const

References CLHEP::Hep3Vector::setRThetaPhi(), G4VHadDecayAlgorithm::UniformPhi(), and G4VHadDecayAlgorithm::UniformTheta().

Referenced by GenerateMultiBody(), G4HadPhaseSpaceNBodyAsai::GenerateMultiBody(), and G4VHadPhaseSpaceAlgorithm::GenerateTwoBody().

Field Documentation

◆ name

G4String G4VHadDecayAlgorithm::name
privateinherited

◆ verboseLevel

G4int G4VHadDecayAlgorithm::verboseLevel
privateinherited

The documentation for this class was generated from the following files: