Geant4  10.00.p02
G4MuonMinusCaptureAtRest.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4MuonMinusCaptureAtRest.cc 81888 2014-06-06 13:17:25Z gcosmo $
27 //
28 // G4MuonMinusCaptureAtRest physics process
29 // Larry Felawka (TRIUMF) and Art Olin (TRIUMF)
30 // April 1998
31 //---------------------------------------------------------------------
32 //
33 // Modifications:
34 // 18/08/2000 V.Ivanchenko Update description
35 // 12/12/2003 H.P.Wellisch Completly rewrite mu-nuclear part
36 // 17/05/2006 V.Ivanchenko Cleanup
37 // 15/11/2006 V.Ivanchenko Review and rewrite all kinematics
38 // 24/01/2007 V.Ivanchenko Force to work with integer Z and A
39 // 23/01/2009 V.Ivanchenko Add deregistration
40 //
41 //-----------------------------------------------------------------------------
42 
44 #include "G4DynamicParticle.hh"
45 #include "Randomize.hh"
46 #include "G4PhysicalConstants.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "G4He3.hh"
49 #include "G4NeutrinoMu.hh"
50 #include "G4Fragment.hh"
52 #include "G4Proton.hh"
53 #include "G4PionPlus.hh"
54 #include "G4GHEKinematicsVector.hh"
55 #include "G4Fancy3DNucleus.hh"
56 #include "G4ExcitationHandler.hh"
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60 
62  G4ProcessType aType ) :
63  G4VRestProcess (processName, aType), nCascade(0), targetZ(0), targetA(0),
64  isInitialised(false)
65 {
67  Cascade = new G4GHEKinematicsVector [17];
70  theN = new G4Fancy3DNucleus();
73  targetMass = 0.0;
74 }
75 
76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77 
79 {
81  delete [] Cascade;
82  delete pSelector;
83  delete pEMCascade;
84  delete theN;
85  delete theHandler;
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
91 {
92  return ( &p == G4MuonMinus::MuonMinus() );
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96 
98 {
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
105 {
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110 
112  const G4Step&)
113 {
114  //
115  // Handles MuonMinuss at rest; a MuonMinus can either create secondaries or
116  // do nothing (in which case it should be sent back to decay-handling
117  // section
118  //
120 
121  // select element and get Z,A.
122  G4Element* aEle = pSelector->GetElement(track.GetMaterial());
123  targetZ = aEle->GetZ();
124  targetA = G4double(G4int(aEle->GetN()+0.5));
125  G4int ni = 0;
126 
127  G4IsotopeVector* isv = aEle->GetIsotopeVector();
128  if(isv) ni = isv->size();
129 
130  if(ni == 1) {
131  targetA = G4double(aEle->GetIsotope(0)->GetN());
132  } else if(ni > 1) {
134  G4double y = G4UniformRand();
135  G4int j = -1;
136  ni--;
137  do {
138  j++;
139  y -= ab[j];
140  } while (y > 0.0 && j < ni);
141  targetA = G4double(aEle->GetIsotope(j)->GetN());
142  }
143 
144  // Do the electromagnetic cascade of the muon in the nuclear field.
145  nCascade = 0;
148 
149  // Decide on Decay or Capture, and doit.
152  G4double lambda = lambdac + lambdad;
153 
154  // === Throw for capture time.
155 
156  G4double tDelay = -std::log(G4UniformRand()) / lambda;
157 
158  G4ReactionProductVector * captureResult = 0;
159  G4int nEmSecondaries = nCascade;
160  G4int nSecondaries = nCascade;
161  /*
162  G4cout << "lambda= " << lambda << " lambdac= " << lambdac
163  << " nem= " << nEmSecondaries << G4endl;
164  */
165  if( G4UniformRand()*lambda > lambdac)
166  pEMCascade->DoBoundMuonMinusDecay(targetZ, &nEmSecondaries, Cascade);
167  else
168  captureResult = DoMuCapture();
169 
170  // fill the final state
171  if(captureResult) nSecondaries += captureResult->size();
172  else nSecondaries = nEmSecondaries;
173  //G4cout << " nsec= " << nSecondaries << " nem= " << nEmSecondaries << G4endl;
174 
175  aParticleChange.SetNumberOfSecondaries( nSecondaries );
176 
177  G4double globalTime = track.GetGlobalTime();
179  // Store nuclear cascade
180  if(captureResult) {
181  G4int n = captureResult->size();
182  for ( G4int isec = 0; isec < n; isec++ ) {
183  G4ReactionProduct* aParticle = captureResult->operator[](isec);
184  G4DynamicParticle * aNewParticle = new G4DynamicParticle();
185  aNewParticle->SetDefinition( aParticle->GetDefinition() );
186  G4LorentzVector itV(aParticle->GetTotalEnergy(), aParticle->GetMomentum());
187  aNewParticle->SetMomentum(itV.vect());
188  G4double localtime = globalTime + tDelay + aParticle->GetTOF();
189  G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position);
190  aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
191  aParticleChange.AddSecondary( aNewTrack );
192  delete aParticle;
193  }
194  delete captureResult;
195  }
196 
197  // Store electromagnetic cascade
198 
199  if(nEmSecondaries > 0) {
200 
201  for ( G4int isec = 0; isec < nEmSecondaries; isec++ ) {
203  G4double localtime = globalTime;
204  if(isec >= nCascade) localtime += tDelay;
205  if(pd) {
206  G4DynamicParticle* aNewParticle = new G4DynamicParticle;
207  aNewParticle->SetDefinition( pd );
208  aNewParticle->SetMomentum( Cascade[isec].GetMomentum() );
209 
210  G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position );
211  aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
212  aParticleChange.AddSecondary( aNewTrack );
213  }
214  }
215  }
216 
219 
220  return &aParticleChange;
221 }
222 
223 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224 
226 {
228  G4double muBindingEnergy = pEMCascade->GetKShellEnergy(targetZ);
229  /*
230  G4cout << "G4MuonMinusCaptureAtRest::DoMuCapture called Emu= "
231  << muBindingEnergy << G4endl;
232  */
233  // Energy on K-shell
234  G4double muEnergy = mumass + muBindingEnergy;
235  G4double muMom = std::sqrt(muBindingEnergy*(muBindingEnergy + 2.0*mumass));
236  G4double availableEnergy = targetMass + mumass - muBindingEnergy;
237  G4LorentzVector momInitial(0.0,0.0,0.0,availableEnergy);
238  G4LorentzVector momResidual;
239 
240  G4ThreeVector vmu = muMom*pEMCascade->GetRandomVec();
241  G4LorentzVector aMuMom(vmu, muEnergy);
242 
243  G4double residualMass =
245 
246  G4ReactionProductVector* aPreResult = 0;
249 
250  G4int iz = G4int(targetZ);
251  G4int ia = G4int(targetA);
252 
253  // p, d, t, 3He or alpha as target
254  if(iz <= 2) {
255 
256  if(ia > 1) {
257  if(iz == 1 && ia == 2) {
258  availableEnergy -= neutron_mass_c2;
259  } else if(iz == 1 && ia == 3) {
260  availableEnergy -= 2.0*neutron_mass_c2;
261  } else if(iz == 2) {
262  G4ParticleDefinition* pd = 0;
263  if (ia == 3) {
264  pd = G4Deuteron::Deuteron();
265  } else if(ia == 4) {
266  pd = G4Triton::Triton();
267  } else {
269  }
270 
271  // G4cout << "Extra " << pd->GetParticleName() << G4endl;
272  availableEnergy -= pd->GetPDGMass();
273  }
274  }
275  //
276  // Computation in assumption of CM collision of mu and nucleaon
277  //
278  G4double Enu = 0.5*(availableEnergy -
279  neutron_mass_c2*neutron_mass_c2/availableEnergy);
280 
281  // make the nu, and transform to lab;
282  G4ThreeVector nu3Mom = Enu*pEMCascade->GetRandomVec();
283 
286  aN->SetTotalEnergy( availableEnergy - Enu );
287  aN->SetMomentum( -nu3Mom );
288 
289  aNu->SetTotalEnergy( Enu );
290  aNu->SetMomentum( nu3Mom );
291  aPreResult = new G4ReactionProductVector();
292 
293  aPreResult->push_back(aN );
294  aPreResult->push_back(aNu);
295 
296  if(verboseLevel > 1)
297  G4cout << "DoMuCapture on H or He"
298  <<" EkinN(MeV)= " << (availableEnergy - Enu - neutron_mass_c2)/MeV
299  <<" Enu(MeV)= "<<aNu->GetTotalEnergy()/MeV
300  <<" n= " << aPreResult->size()
301  <<G4endl;
302 
303  return aPreResult;
304  }
305 
306  // pick random proton inside nucleus
307  G4double eEx;
308  do {
309  theN->Init(ia, iz);
310  G4LorentzVector thePMom;
311  G4Nucleon * aNucleon = 0;
312  G4int theProtonCounter = G4int( targetZ * G4UniformRand() );
313  G4int counter = 0;
314  theN->StartLoop();
315 
316  while( (aNucleon=theN->GetNextNucleon()) ) {
317 
318  if( aNucleon->GetDefinition() == G4Proton::Proton() ) {
319  counter++;
320  if(counter == theProtonCounter) {
321  thePMom = aNucleon->GetMomentum();
322  break;
323  }
324  }
325  }
326 
327  // Get the nu momentum in the CMS
328  G4LorentzVector theCMS = thePMom + aMuMom;
329  G4ThreeVector bst = theCMS.boostVector();
330 
331  G4double Ecms = theCMS.mag();
332  G4double Enu = 0.5*(Ecms - neutron_mass_c2*neutron_mass_c2/Ecms);
333  eEx = 0.0;
334 
335  if(Enu > 0.0) {
336  // make the nu, and transform to lab;
337  G4ThreeVector nu3Mom = Enu*pEMCascade->GetRandomVec();
338  G4LorentzVector nuMom(nu3Mom, Enu);
339 
340  // nu in lab.
341  nuMom.boost(bst);
342  aNu->SetTotalEnergy( nuMom.e() );
343  aNu->SetMomentum( nuMom.vect() );
344 
345  // make residual
346  momResidual = momInitial - nuMom;
347 
348  // Call pre-compound on the rest.
349  eEx = momResidual.mag();
350  if(verboseLevel > 1)
351  G4cout << "G4MuonMinusCaptureAtRest::DoMuCapture: "
352  << " Eex(MeV)= " << (eEx-residualMass)/MeV
353  << " Enu(MeV)= "<<aNu->GetTotalEnergy()/MeV
354  <<G4endl;
355  }
356  } while(eEx <= residualMass);
357 
358 // G4cout << "muonCapture : " << eEx << " " << residualMass
359 // << " A,Z= " << targetA << ", "<< targetZ
360 // << " " << G4int(targetA) << ", " << G4int(targetZ) << G4endl;
361 
362  //
363  // Start Deexcitation
364  //
365  G4ThreeVector fromBreit = momResidual.boostVector();
366  G4LorentzVector fscm(0.0,0.0,0.0, eEx);
367  G4Fragment anInitialState;
368  anInitialState.SetA(targetA);
369  anInitialState.SetZ(G4double(iz - 1));
370  anInitialState.SetNumberOfParticles(2);
371  anInitialState.SetNumberOfCharged(0);
372  anInitialState.SetNumberOfHoles(1);
373  anInitialState.SetMomentum(fscm);
374  aPreResult = theHandler->BreakItUp(anInitialState);
375 
376  G4ReactionProductVector::iterator ires;
377  G4double eBal = availableEnergy - aNu->GetTotalEnergy();
378  for(ires=aPreResult->begin(); ires!=aPreResult->end(); ires++) {
379  G4LorentzVector itV((*ires)->GetTotalEnergy(), (*ires)->GetMomentum());
380  itV.boost(fromBreit);
381  (*ires)->SetTotalEnergy(itV.t());
382  (*ires)->SetMomentum(itV.vect());
383  eBal -= itV.t();
384  }
385  //
386  // fill neutrino into result
387  //
388  aPreResult->push_back(aNu);
389 
390  if(verboseLevel > 1)
391  G4cout << "DoMuCapture: Nsec= "
392  << aPreResult->size() << " Ebalance(MeV)= " << eBal/MeV
393  <<" E0(MeV)= " <<availableEnergy/MeV
394  <<" Mres(GeV)= " <<residualMass/GeV
395  <<G4endl;
396 
397  return aPreResult;
398 }
399 
void DeRegisterExtraProcess(G4VProcess *)
G4MuonMinusCaptureAtRest(const G4String &processName="muMinusCaptureAtRest", G4ProcessType aType=fHadronic)
static const double MeV
Definition: G4SIunits.hh:193
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::vector< G4Isotope * > G4IsotopeVector
void SetMomentum(const G4ThreeVector &momentum)
const G4LorentzVector & GetMomentum() const
Definition: G4Nucleon.hh:71
G4int verboseLevel
Definition: G4VProcess.hh:368
CLHEP::Hep3Vector G4ThreeVector
static G4HadronicProcessStore * Instance()
G4double GetN() const
Definition: G4Element.hh:134
G4ParticleDefinition * GetParticleDef()
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4MuMinusCaptureCascade * pEMCascade
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
const G4ThreeVector & GetPosition() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetMuonCaptureRate(G4double Z, G4double A)
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:355
G4double GetMuonDecayRate(G4double Z, G4double A)
void SetTouchableHandle(const G4TouchableHandle &apValue)
int G4int
Definition: G4Types.hh:78
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState) const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleDefinition * GetDefinition() const
void SetA(G4double value)
Definition: G4Fragment.hh:314
std::vector< G4ReactionProduct * > G4ReactionProductVector
virtual G4bool IsApplicable(const G4ParticleDefinition &)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
G4IonTable * GetIonTable() const
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4StopElementSelector * pSelector
static G4NeutrinoMu * NeutrinoMu()
Definition: G4NeutrinoMu.cc:85
bool G4bool
Definition: G4Types.hh:79
G4int DoCascade(const G4double Z, const G4double A, G4GHEKinematicsVector *Cascade)
G4double iz
Definition: TRTMaterials.hh:39
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:276
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void SetTotalEnergy(const G4double en)
void Init(G4int theA, G4int theZ)
static G4Triton * Triton()
Definition: G4Triton.cc:95
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double GeV
Definition: G4SIunits.hh:196
Definition: G4Step.hh:76
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4int n
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:364
G4double GetGlobalTime() const
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
const G4TouchableHandle & GetTouchableHandle() const
G4Material * GetMaterial() const
void RegisterExtraProcess(G4VProcess *)
virtual void Initialize(const G4Track &)
G4double GetTotalEnergy() const
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
static const G4double ab
G4double GetKShellEnergy(G4double Z)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
void SetNumberOfSecondaries(G4int totSecondaries)
int position
Definition: filter.cc:7
G4double GetTOF() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4Element * GetElement(const G4Material *aMaterial)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
void DoBoundMuonMinusDecay(G4double Z, G4int *nCascade, G4GHEKinematicsVector *Cascade)
void AddSecondary(G4Track *aSecondary)
G4ThreeVector GetMomentum() const
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
#define G4endl
Definition: G4ios.hh:61
virtual G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:85
void SetZ(G4double value)
Definition: G4Fragment.hh:308
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:369
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4GHEKinematicsVector * Cascade
G4Nucleon * GetNextNucleon()
G4ReactionProductVector * DoMuCapture()
void PrintInfo(const G4ParticleDefinition *)
CLHEP::HepLorentzVector G4LorentzVector
G4ProcessType