Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VHadDecayAlgorithm Class Referenceabstract

#include <G4VHadDecayAlgorithm.hh>

Inheritance diagram for G4VHadDecayAlgorithm:

Public Member Functions

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

Protected Member Functions

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
 
virtual G4bool IsDecayAllowed (G4double initialMass, const std::vector< G4double > &masses) const
 
G4double TwoBodyMomentum (G4double M0, G4double M1, G4double M2) const
 
G4double UniformTheta () const
 
G4double UniformPhi () const
 
void PrintVector (const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
 

Detailed Description

Definition at line 43 of file G4VHadDecayAlgorithm.hh.

Constructor & Destructor Documentation

G4VHadDecayAlgorithm::G4VHadDecayAlgorithm ( const G4String algName,
G4int  verbose = 0 
)
inline

Definition at line 45 of file G4VHadDecayAlgorithm.hh.

46  : name(algName), verboseLevel(verbose) {;}
const XML_Char * name
Definition: expat.h:151
virtual G4VHadDecayAlgorithm::~G4VHadDecayAlgorithm ( )
inlinevirtual

Definition at line 47 of file G4VHadDecayAlgorithm.hh.

47 {;}

Member Function Documentation

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

Definition at line 49 of file G4VHadDecayAlgorithm.cc.

51  {
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 }
virtual void GenerateTwoBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0
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
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

virtual void G4VHadDecayAlgorithm::GenerateMultiBody ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protectedpure virtual

Implemented in G4CascadeFinalStateAlgorithm, G4HadPhaseSpaceKopylov, G4HadPhaseSpaceNBodyAsai, and G4HadPhaseSpaceGenbod.

Here is the caller graph for this function:

virtual void G4VHadDecayAlgorithm::GenerateTwoBody ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protectedpure virtual

Implemented in G4CascadeFinalStateAlgorithm, and G4VHadPhaseSpaceAlgorithm.

Here is the caller graph for this function:

const G4String& G4VHadDecayAlgorithm::GetName ( ) const
inline

Definition at line 57 of file G4VHadDecayAlgorithm.hh.

57 { return name; }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

G4int G4VHadDecayAlgorithm::GetVerboseLevel ( ) const
inline

Definition at line 56 of file G4VHadDecayAlgorithm.hh.

56 { return verboseLevel; }

Here is the caller graph for this function:

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

Definition at line 69 of file G4VHadDecayAlgorithm.cc.

70  {
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 }
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void PrintVector(const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
const G4String & GetName() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 123 of file G4VHadDecayAlgorithm.cc.

124  {
125  os << " " << vname << "(" << v.size() << ") ";
126  std::copy(v.begin(), v.end(), std::ostream_iterator<G4double>(os, " "));
127  os << std::endl;
128 }

Here is the caller graph for this function:

virtual void G4VHadDecayAlgorithm::SetVerboseLevel ( G4int  verbose)
inlinevirtual

Reimplemented in G4CascadeFinalStateAlgorithm.

Definition at line 55 of file G4VHadDecayAlgorithm.hh.

55 { verboseLevel = verbose; }

Here is the caller graph for this function:

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

Definition at line 91 of file G4VHadDecayAlgorithm.cc.

92  {
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 }
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
const G4String & GetName() const
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VHadDecayAlgorithm::UniformPhi ( ) const
protected

Definition at line 115 of file G4VHadDecayAlgorithm.cc.

115  {
116  return twopi*G4UniformRand();
117 }
static constexpr double twopi
Definition: G4SIunits.hh:76
#define G4UniformRand()
Definition: Randomize.hh:97

Here is the caller graph for this function:

G4double G4VHadDecayAlgorithm::UniformTheta ( ) const
protected

Definition at line 111 of file G4VHadDecayAlgorithm.cc.

111  {
112  return std::acos(2.0*G4UniformRand() - 1.0);
113 }
#define G4UniformRand()
Definition: Randomize.hh:97

Here is the caller graph for this function:


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