G4MuonMinusCaptureAtRest Class Reference

#include <G4MuonMinusCaptureAtRest.hh>

Inheritance diagram for G4MuonMinusCaptureAtRest:

G4VRestProcess G4VProcess

Public Member Functions

 G4MuonMinusCaptureAtRest (const G4String &processName="muMinusCaptureAtRest", G4ProcessType aType=fHadronic)
virtual ~G4MuonMinusCaptureAtRest ()
virtual G4bool IsApplicable (const G4ParticleDefinition &)
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
G4double GetMeanLifeTime (const G4Track &, G4ForceCondition *)

Detailed Description

Definition at line 64 of file G4MuonMinusCaptureAtRest.hh.


Constructor & Destructor Documentation

G4MuonMinusCaptureAtRest::G4MuonMinusCaptureAtRest ( const G4String processName = "muMinusCaptureAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 61 of file G4MuonMinusCaptureAtRest.cc.

References fHadronAtRest, G4HadronicProcessStore::Instance(), G4HadronicProcessStore::RegisterExtraProcess(), and G4VProcess::SetProcessSubType().

00062                                                                            :
00063   G4VRestProcess (processName, aType), nCascade(0), targetZ(0), targetA(0), 
00064   isInitialised(false)
00065 {
00066   SetProcessSubType(fHadronAtRest);
00067   Cascade    = new G4GHEKinematicsVector [17];
00068   pSelector  = new G4StopElementSelector();
00069   pEMCascade = new G4MuMinusCaptureCascade();
00070   theN       = new G4Fancy3DNucleus();
00071   theHandler = new G4ExcitationHandler();
00072   G4HadronicProcessStore::Instance()->RegisterExtraProcess(this);
00073 }

G4MuonMinusCaptureAtRest::~G4MuonMinusCaptureAtRest (  )  [virtual]

Definition at line 77 of file G4MuonMinusCaptureAtRest.cc.

References G4HadronicProcessStore::DeRegisterExtraProcess(), and G4HadronicProcessStore::Instance().

00078 {
00079   G4HadronicProcessStore::Instance()->DeRegisterExtraProcess(this);
00080   delete [] Cascade;
00081   delete pSelector;
00082   delete pEMCascade;
00083   delete theN;
00084   delete theHandler;
00085 }


Member Function Documentation

G4VParticleChange * G4MuonMinusCaptureAtRest::AtRestDoIt ( const G4Track ,
const G4Step  
) [virtual]

Reimplemented from G4VRestProcess.

Definition at line 110 of file G4MuonMinusCaptureAtRest.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, G4MuMinusCaptureCascade::DoBoundMuonMinusDecay(), G4MuMinusCaptureCascade::DoCascade(), fStopAndKill, G4UniformRand, G4StopElementSelector::GetElement(), G4Track::GetGlobalTime(), G4Element::GetIsotope(), G4Element::GetIsotopeVector(), G4Track::GetMaterial(), G4StopElementSelector::GetMuonCaptureRate(), G4StopElementSelector::GetMuonDecayRate(), G4Isotope::GetN(), G4Element::GetN(), G4NucleiProperties::GetNuclearMass(), G4Track::GetPosition(), G4Element::GetRelativeAbundanceVector(), G4Track::GetTouchableHandle(), G4Element::GetZ(), G4ParticleChange::Initialize(), G4InuclParticleNames::lambda, CLHEP::detail::n, position, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::SetNumberOfSecondaries(), and G4Track::SetTouchableHandle().

00112 {
00113   //
00114   // Handles MuonMinuss at rest; a MuonMinus can either create secondaries or
00115   // do nothing (in which case it should be sent back to decay-handling
00116   // section
00117   //
00118   aParticleChange.Initialize(track);
00119 
00120   // select element and get Z,A.
00121   G4Element* aEle = pSelector->GetElement(track.GetMaterial());
00122   targetZ = aEle->GetZ();
00123   targetA = G4double(G4int(aEle->GetN()+0.5)); 
00124   G4int ni = 0;
00125 
00126   G4IsotopeVector* isv = aEle->GetIsotopeVector();
00127   if(isv) ni = isv->size();
00128 
00129   if(ni == 1) {
00130     targetA = G4double(aEle->GetIsotope(0)->GetN());
00131   } else if(ni > 1) {
00132     G4double* ab = aEle->GetRelativeAbundanceVector();
00133     G4double y = G4UniformRand();
00134     G4int j = -1;
00135     ni--;
00136     do {
00137       j++;
00138       y -= ab[j];
00139     } while (y > 0.0 && j < ni);
00140     targetA = G4double(aEle->GetIsotope(j)->GetN());
00141   }
00142   
00143   // Do the electromagnetic cascade of the muon in the nuclear field.
00144   nCascade   = 0;
00145   targetMass = G4NucleiProperties::GetNuclearMass(targetA, targetZ);
00146   nCascade   = pEMCascade->DoCascade(targetZ, targetMass, Cascade);
00147 
00148   // Decide on Decay or Capture, and doit.
00149   G4double lambdac  = pSelector->GetMuonCaptureRate(targetZ, targetA);
00150   G4double lambdad  = pSelector->GetMuonDecayRate(targetZ, targetA);
00151   G4double lambda   = lambdac + lambdad;
00152 
00153   // ===  Throw for capture  time.
00154 
00155   G4double tDelay = -std::log(G4UniformRand()) / lambda;
00156   
00157   G4ReactionProductVector * captureResult = 0;
00158   G4int nEmSecondaries = nCascade;
00159   G4int nSecondaries = nCascade;
00160   /*
00161   G4cout << "lambda= " << lambda << " lambdac= " << lambdac 
00162          << " nem= " << nEmSecondaries << G4endl;
00163   */
00164   if( G4UniformRand()*lambda > lambdac) 
00165     pEMCascade->DoBoundMuonMinusDecay(targetZ, &nEmSecondaries, Cascade);
00166   else 
00167     captureResult = DoMuCapture();
00168   
00169   // fill the final state
00170   if(captureResult) nSecondaries += captureResult->size();
00171   else nSecondaries = nEmSecondaries;
00172   //G4cout << " nsec= " << nSecondaries << " nem= " << nEmSecondaries << G4endl;
00173 
00174   aParticleChange.SetNumberOfSecondaries( nSecondaries );
00175 
00176   G4double globalTime = track.GetGlobalTime();
00177   G4ThreeVector position = track.GetPosition();
00178   // Store nuclear cascade
00179   if(captureResult) {
00180     G4int n = captureResult->size();
00181     for ( G4int isec = 0; isec < n; isec++ ) {
00182       G4ReactionProduct* aParticle = captureResult->operator[](isec);
00183       G4DynamicParticle * aNewParticle = new G4DynamicParticle();
00184       aNewParticle->SetDefinition( aParticle->GetDefinition() );
00185       G4LorentzVector itV(aParticle->GetTotalEnergy(), aParticle->GetMomentum());
00186       aNewParticle->SetMomentum(itV.vect());
00187       G4double localtime = globalTime + tDelay + aParticle->GetTOF();
00188       G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position);
00189       aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
00190       aParticleChange.AddSecondary( aNewTrack );
00191       delete aParticle;
00192     }
00193     delete captureResult;
00194   }
00195   
00196   // Store electromagnetic cascade
00197 
00198   if(nEmSecondaries > 0) {
00199 
00200     for ( G4int isec = 0; isec < nEmSecondaries; isec++ ) {
00201       G4ParticleDefinition* pd = Cascade[isec].GetParticleDef();
00202       G4double localtime = globalTime;
00203       if(isec >= nCascade) localtime += tDelay;
00204       if(pd) {
00205         G4DynamicParticle* aNewParticle = new G4DynamicParticle;
00206         aNewParticle->SetDefinition( pd );
00207         aNewParticle->SetMomentum( Cascade[isec].GetMomentum() );
00208 
00209         G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position );
00210         aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
00211         aParticleChange.AddSecondary( aNewTrack );
00212       }
00213     }
00214   }
00215 
00216   aParticleChange.ProposeLocalEnergyDeposit(0.0);
00217   aParticleChange.ProposeTrackStatus(fStopAndKill); 
00218 
00219   return &aParticleChange;
00220 }

void G4MuonMinusCaptureAtRest::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 103 of file G4MuonMinusCaptureAtRest.cc.

References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::PrintInfo().

00104 {
00105   G4HadronicProcessStore::Instance()->PrintInfo(&p);
00106 }

G4double G4MuonMinusCaptureAtRest::GetMeanLifeTime ( const G4Track ,
G4ForceCondition  
) [inline, virtual]

Implements G4VRestProcess.

Definition at line 82 of file G4MuonMinusCaptureAtRest.hh.

00083   {return 0;};

G4bool G4MuonMinusCaptureAtRest::IsApplicable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 89 of file G4MuonMinusCaptureAtRest.cc.

References G4MuonMinus::MuonMinus().

00090 {
00091   return ( &p == G4MuonMinus::MuonMinus() );
00092 }

void G4MuonMinusCaptureAtRest::PreparePhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 96 of file G4MuonMinusCaptureAtRest.cc.

References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::RegisterParticleForExtraProcess().


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