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

#include <G4HadPhaseSpaceNBodyAsai.hh>

Inheritance diagram for G4HadPhaseSpaceNBodyAsai:
G4VHadPhaseSpaceAlgorithm G4VHadDecayAlgorithm

Public Member Functions

 G4HadPhaseSpaceNBodyAsai (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 ~G4HadPhaseSpaceNBodyAsai ()
 

Protected Member Functions

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 G4HadPhaseSpaceNBodyAsai.hh.

Constructor & Destructor Documentation

◆ G4HadPhaseSpaceNBodyAsai()

G4HadPhaseSpaceNBodyAsai::G4HadPhaseSpaceNBodyAsai ( G4int  verbose = 0)
inline

Definition at line 39 of file G4HadPhaseSpaceNBodyAsai.hh.

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

◆ ~G4HadPhaseSpaceNBodyAsai()

virtual G4HadPhaseSpaceNBodyAsai::~G4HadPhaseSpaceNBodyAsai ( )
inlinevirtual

Definition at line 41 of file G4HadPhaseSpaceNBodyAsai.hh.

41{;}

Member Function Documentation

◆ 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}
#define G4endl
Definition: G4ios.hh:57
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 G4HadPhaseSpaceNBodyAsai::GenerateMultiBody ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protectedvirtual

Implements G4VHadDecayAlgorithm.

Definition at line 50 of file G4HadPhaseSpaceNBodyAsai.cc.

53 {
54 if (GetVerboseLevel()) G4cout << GetName() << "::GenerateMultiBody" << G4endl;
55
56 finalState.clear();
57
58 //daughters' mass
59 G4int numberOfDaughters = masses.size();
60 G4double sumofmasses =
61 std::accumulate(masses.begin(), masses.end(), 0.);
62
63 //Calculate daughter momentum
64 std::vector<G4double> daughtermomentum(numberOfDaughters);
65 std::vector<G4double> sm(numberOfDaughters);
66 G4double tmas;
67 G4double weight = 1.0;
68 G4int numberOfTry = 0;
69 G4int i;
70
71 std::vector<G4double> rd(numberOfDaughters);
72 do {
73 //Generate random number in descending order
74 rd[0] = 1.0;
75 std::generate(rd.begin()+1, rd.end(), uniformRand);
76 std::sort(rd.begin(), rd.end(), std::greater<G4double>());
77
78 if (GetVerboseLevel()>1) PrintVector(rd,"rd",G4cout);
79
80 //calcurate virtual mass
81 tmas = initialMass - sumofmasses;
82 G4double temp = sumofmasses;
83 for(i =0; i < numberOfDaughters; i++) {
84 sm[i] = rd[i]*tmas + temp;
85 temp -= masses[i];
86 if (GetVerboseLevel()>1) {
87 G4cout << i << " random number:" << rd[i]
88 << " virtual mass:" << sm[i]/GeV << " GeV/c2" <<G4endl;
89 }
90 }
91
92 //Calculate daughter momentum
93 weight = 1.0;
94 i = numberOfDaughters-1;
95 daughtermomentum[i] = TwoBodyMomentum(sm[i-1],masses[i-1],sm[i]);
96 if (GetVerboseLevel()>1) {
97 G4cout << " daughter " << i << ": momentum "
98 << daughtermomentum[i]/GeV << " GeV/c" <<G4endl;
99 }
100 for(i =numberOfDaughters-2; i>=0; i--) {
101 // calculate
102 daughtermomentum[i] = TwoBodyMomentum(sm[i],masses[i],sm[i+1]);
103 if(daughtermomentum[i] < 0.0) {
104 // !!! illegal momentum !!!
105 if (GetVerboseLevel()>0) {
106 G4cout << "G4HadPhaseSpaceNBodyAsai::Generate "
107 << " can not calculate daughter momentum "
108 << "\n initialMass " << initialMass/GeV << " GeV/c2"
109 << "\n daughter " << i << ": mass "
110 << masses[i]/GeV << " GeV/c2; momentum "
111 << daughtermomentum[i]/GeV << " GeV/c" << G4endl;
112 }
113 return; // Error detection
114 }
115
116 // calculate weight of this events
117 weight *= daughtermomentum[i]/sm[i];
118 if (GetVerboseLevel()>1) {
119 G4cout << " daughter " << i << ": momentum "
120 << daughtermomentum[i]/GeV << " GeV/c" <<G4endl;
121 }
122 }
123 if (GetVerboseLevel()>1) {
124 G4cout << " weight: " << weight <<G4endl;
125 }
126
127 // exit if number of Try exceeds 100
128 if (numberOfTry++ > 100) {
129 if (GetVerboseLevel()>0) {
130 G4cout << "G4HadPhaseSpaceNBodyAsai::Generate "
131 << " can not determine Decay Kinematics " << G4endl;
132 }
133 return; // Error detection
134 }
135 } while (weight > G4UniformRand()); /* Loop checking, 02.11.2015, A.Ribon */
136
137 if (GetVerboseLevel()>1) {
138 G4cout << "Start calculation of daughters momentum vector "<<G4endl;
139 }
140
142
143 finalState.resize(numberOfDaughters);
144
145 i = numberOfDaughters-2;
146
147 G4ThreeVector direction = UniformVector(daughtermomentum[i]);
148
149 finalState[i].setVectM(direction, masses[i]);
150 finalState[i+1].setVectM(-direction, masses[i+1]);
151
152 for (i = numberOfDaughters-3; i >= 0; i--) {
153 direction = UniformVector();
154
155 //create daughter particle
156 finalState[i].setVectM(-daughtermomentum[i]*direction, masses[i]);
157
158 // boost already created particles
159 beta = daughtermomentum[i];
160 beta /= std::sqrt(beta*beta + sm[i+1]*sm[i+1]);
161 for (G4int j = i+1; j<numberOfDaughters; j++) {
162 finalState[j].boost(beta*direction);
163 }
164 }
165}
static constexpr double GeV
Definition: G4SIunits.hh:203
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
void PrintVector(const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const
G4ThreeVector UniformVector(G4double mag=1.) const
void generate(const G4double sqrtS, ParticleList &particles)
Generate an event in the CM system.

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, G4cout, G4endl, G4UniformRand, G4INCL::PhaseSpaceGenerator::generate(), G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), GeV, G4VHadDecayAlgorithm::PrintVector(), G4InuclParticleNames::sm, G4VHadDecayAlgorithm::TwoBodyMomentum(), anonymous_namespace{G4HadPhaseSpaceNBodyAsai.cc}::uniformRand(), 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

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(), 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 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(), G4HadPhaseSpaceKopylov::GenerateMultiBody(), 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 G4HadPhaseSpaceKopylov::GenerateMultiBody(), 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: