Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4StopElementSelector Class Reference

#include <G4StopElementSelector.hh>

Public Member Functions

 G4StopElementSelector ()
 
 ~G4StopElementSelector ()
 
G4ElementGetElement (const G4Material *aMaterial)
 
G4double GetMuonCaptureRate (G4double Z, G4double A)
 
G4double GetMuonDecayRate (G4double Z, G4double A)
 

Detailed Description

Definition at line 54 of file G4StopElementSelector.hh.

Constructor & Destructor Documentation

G4StopElementSelector::G4StopElementSelector ( )

Definition at line 53 of file G4StopElementSelector.cc.

54 { }
G4StopElementSelector::~G4StopElementSelector ( )

Definition at line 59 of file G4StopElementSelector.cc.

60 { }

Member Function Documentation

G4Element * G4StopElementSelector::GetElement ( const G4Material aMaterial)

Definition at line 65 of file G4StopElementSelector.cc.

References G4UniformRand, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), and G4Material::GetNumberOfElements().

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

66 {
67  // Fermi-Teller Z-low of mu- capture and exceptions
68  // for halogens and oxigen.
69  // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
70  G4int i;
71  G4double Z;
72  const G4int numberOfElements = aMaterial->GetNumberOfElements();
73  const G4ElementVector* theElementVector = aMaterial->GetElementVector();
74 
75  if(1 == numberOfElements) return (*theElementVector)[0];
76 
77  const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
78 
79  G4double sum = 0.0;
80  for ( i=0; i < numberOfElements; i++ ) {
81 
82  Z = (*theElementVector)[i]->GetZ();
83 
84  // Halogens
85  if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
86  sum += 0.66 * Z * theAtomicNumberDensity[i] ;
87 
88  // Oxigen
89  } else if( 8.0 == Z ) {
90  sum += 0.56 * Z * theAtomicNumberDensity[i] ;
91 
92  // Others
93  } else {
94  sum += Z * theAtomicNumberDensity[i] ;
95  }
96  }
97 
98  G4double random = G4UniformRand() * sum;
99  sum = 0.0 ;
100  i = -1;
101 
102  // Selection of element
103  do {
104  i++;
105  Z = (*theElementVector)[i]->GetZ();
106 
107  // Galogens
108  if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
109  sum += 0.66 * Z * theAtomicNumberDensity[i] ;
110 
111  // Oxigen
112  } else if( 8.0 == Z ) {
113  sum += 0.56 * Z * theAtomicNumberDensity[i] ;
114 
115  // Others
116  } else {
117  sum += Z * theAtomicNumberDensity[i] ;
118  }
119  } while ( (sum < random) && (i < numberOfElements - 1) );
120 
121  return (*theElementVector)[i];
122 }
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
G4double G4StopElementSelector::GetMuonCaptureRate ( G4double  Z,
G4double  A 
)

Definition at line 126 of file G4StopElementSelector.cc.

References G4MuonMinusBoundDecay::GetMuonCaptureRate().

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

127 {
129 }
int G4int
Definition: G4Types.hh:78
static G4double GetMuonCaptureRate(G4int Z, G4int A)
G4double G4StopElementSelector::GetMuonDecayRate ( G4double  Z,
G4double  A 
)

Definition at line 133 of file G4StopElementSelector.cc.

References G4MuonMinusBoundDecay::GetMuonDecayRate().

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

134 {
136 }
int G4int
Definition: G4Types.hh:78
static G4double GetMuonDecayRate(G4int Z)

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