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.

00054 { }

G4StopElementSelector::~G4StopElementSelector (  ) 

Definition at line 59 of file G4StopElementSelector.cc.

00060 { }


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 G4ProtonAntiProtonAtRestChips::AtRestDoIt(), G4PionMinusNuclearAtRestChips::AtRestDoIt(), and G4MuonMinusCaptureAtRest::AtRestDoIt().

00066 {
00067   // Fermi-Teller Z-low of mu- capture and exceptions 
00068   // for halogens and oxigen.
00069   // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
00070   G4int i;
00071   G4double Z;
00072   const G4int numberOfElements = aMaterial->GetNumberOfElements();
00073   const G4ElementVector* theElementVector = aMaterial->GetElementVector();
00074 
00075   if(1 == numberOfElements) return (*theElementVector)[0];
00076     
00077   const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
00078 
00079   G4double sum = 0.0;
00080   for ( i=0; i < numberOfElements; i++ ) {
00081 
00082     Z = (*theElementVector)[i]->GetZ();
00083 
00084       // Halogens
00085     if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
00086       sum += 0.66 * Z * theAtomicNumberDensity[i] ; 
00087 
00088       // Oxigen
00089     } else if( 8.0 == Z ) {
00090       sum += 0.56 * Z * theAtomicNumberDensity[i] ; 
00091 
00092       // Others
00093     } else {
00094       sum +=        Z * theAtomicNumberDensity[i] ; 
00095     }
00096   }
00097 
00098   G4double random = G4UniformRand() * sum;
00099   sum = 0.0 ;
00100   i   = -1;
00101 
00102   // Selection of element
00103   do {
00104     i++;
00105     Z = (*theElementVector)[i]->GetZ();
00106 
00107       // Galogens
00108     if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
00109       sum += 0.66 * Z * theAtomicNumberDensity[i] ; 
00110 
00111       // Oxigen
00112     } else if( 8.0 == Z ) {
00113       sum += 0.56 * Z * theAtomicNumberDensity[i] ; 
00114 
00115       // Others
00116     } else {
00117       sum +=        Z * theAtomicNumberDensity[i] ; 
00118     }
00119   } while ( (sum < random) && (i < numberOfElements - 1) );
00120 
00121   return (*theElementVector)[i];
00122 }

G4double G4StopElementSelector::GetMuonCaptureRate ( G4double  Z,
G4double  A 
)

Definition at line 126 of file G4StopElementSelector.cc.

References G4MuonMinusBoundDecay::GetMuonCaptureRate().

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

00127 {
00128   return G4MuonMinusBoundDecay::GetMuonCaptureRate(G4int(Z),G4int(A));
00129 }

G4double G4StopElementSelector::GetMuonDecayRate ( G4double  Z,
G4double  A 
)

Definition at line 133 of file G4StopElementSelector.cc.

References G4MuonMinusBoundDecay::GetMuonDecayRate().

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt().

00134 {
00135   return G4MuonMinusBoundDecay::GetMuonDecayRate(G4int(Z));
00136 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:27 2013 for Geant4 by  doxygen 1.4.7