source: trunk/source/processes/hadronic/stopping/src/G4MuonMinusCaptureAtRest.cc @ 1050

Last change on this file since 1050 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

File size: 13.1 KB
Line 
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,v 1.53 2008/10/02 20:57:52 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29//   G4MuonMinusCaptureAtRest physics process
30//   Larry Felawka (TRIUMF) and Art Olin (TRIUMF)
31//   April 1998
32//---------------------------------------------------------------------
33//
34// Modifications:
35// 18/08/2000  V.Ivanchenko Update description
36// 12/12/2003  H.P.Wellisch Completly rewrite mu-nuclear part
37// 17/05/2006  V.Ivanchenko Cleanup
38// 15/11/2006  V.Ivanchenko Review and rewrite all kinematics
39// 24/01/2007  V.Ivanchenko Force to work with integer Z and A
40//
41//-----------------------------------------------------------------------------
42
43#include "G4MuonMinusCaptureAtRest.hh"
44#include "G4DynamicParticle.hh"
45#include "Randomize.hh"
46#include "G4He3.hh"
47#include "G4NeutrinoMu.hh"
48#include "G4Fragment.hh"
49#include "G4ReactionProductVector.hh"
50#include "G4Proton.hh"
51#include "G4PionPlus.hh"
52#include "G4GHEKinematicsVector.hh"
53#include "G4Fancy3DNucleus.hh"
54#include "G4ExcitationHandler.hh"
55#include "G4HadronicProcessStore.hh"
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59G4MuonMinusCaptureAtRest::G4MuonMinusCaptureAtRest(const G4String& processName,
60                                                   G4ProcessType   aType ) :
61  G4VRestProcess (processName, aType), nCascade(0), targetZ(0), targetA(0), 
62  isInitialised(false)
63{
64  SetProcessSubType(fHadronAtRest);
65  Cascade    = new G4GHEKinematicsVector [17];
66  pSelector  = new G4StopElementSelector();
67  pEMCascade = new G4MuMinusCaptureCascade();
68  theN       = new G4Fancy3DNucleus();
69  theHandler = new G4ExcitationHandler();
70  G4HadronicProcessStore::Instance()->RegisterExtraProcess(this);
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
75G4MuonMinusCaptureAtRest::~G4MuonMinusCaptureAtRest()
76{
77  delete [] Cascade;
78  delete pSelector;
79  delete pEMCascade;
80  delete theN;
81  delete theHandler;
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
86G4bool G4MuonMinusCaptureAtRest::IsApplicable(const G4ParticleDefinition& p)
87{
88  return ( &p == G4MuonMinus::MuonMinus() );
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
93void G4MuonMinusCaptureAtRest::PreparePhysicsTable(const G4ParticleDefinition& p) 
94{
95  G4HadronicProcessStore::Instance()->RegisterParticleForExtraProcess(this, &p);
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99
100void G4MuonMinusCaptureAtRest::BuildPhysicsTable(const G4ParticleDefinition& p) 
101{
102  G4HadronicProcessStore::Instance()->PrintInfo(&p);
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
107G4VParticleChange* G4MuonMinusCaptureAtRest::AtRestDoIt(const G4Track& track, 
108                                                        const G4Step&)
109{
110  //
111  // Handles MuonMinuss at rest; a MuonMinus can either create secondaries or
112  // do nothing (in which case it should be sent back to decay-handling
113  // section
114  //
115  aParticleChange.Initialize(track);
116
117  // select element and get Z,A.
118  G4Element* aEle = pSelector->GetElement(track.GetMaterial());
119  targetZ = aEle->GetZ();
120  targetA = G4double(G4int(aEle->GetN()+0.5)); 
121  G4int ni = 0;
122
123  G4IsotopeVector* isv = aEle->GetIsotopeVector();
124  if(isv) ni = isv->size();
125
126  if(ni == 1) {
127    targetA = G4double(aEle->GetIsotope(0)->GetN());
128  } else if(ni > 1) {
129    G4double* ab = aEle->GetRelativeAbundanceVector();
130    G4double y = G4UniformRand();
131    G4int j = -1;
132    ni--;
133    do {
134      j++;
135      y -= ab[j];
136    } while (y > 0.0 && j < ni);
137    targetA = G4double(aEle->GetIsotope(j)->GetN());
138  }
139 
140  // Do the electromagnetic cascade of the muon in the nuclear field.
141  nCascade   = 0;
142  targetMass = G4NucleiProperties::GetNuclearMass(targetA, targetZ);
143  nCascade   = pEMCascade->DoCascade(targetZ, targetMass, Cascade);
144
145  // Decide on Decay or Capture, and doit.
146  G4double lambdac  = pSelector->GetMuonCaptureRate(targetZ, targetA);
147  G4double lambdad  = pSelector->GetMuonDecayRate(targetZ, targetA);
148  G4double lambda   = lambdac + lambdad;
149
150  // ===  Throw for capture  time.
151
152  G4double tDelay = -std::log(G4UniformRand()) / lambda;
153 
154  G4ReactionProductVector * captureResult = 0;
155  G4int nEmSecondaries = nCascade;
156  G4int nSecondaries = nCascade;
157  /*
158  G4cout << "lambda= " << lambda << " lambdac= " << lambdac
159         << " nem= " << nEmSecondaries << G4endl;
160  */
161  if( G4UniformRand()*lambda > lambdac) 
162    pEMCascade->DoBoundMuonMinusDecay(targetZ, &nEmSecondaries, Cascade);
163  else 
164    captureResult = DoMuCapture();
165 
166  // fill the final state
167  if(captureResult) nSecondaries += captureResult->size();
168  else nSecondaries = nEmSecondaries;
169  //G4cout << " nsec= " << nSecondaries << " nem= " << nEmSecondaries << G4endl;
170
171  aParticleChange.SetNumberOfSecondaries( nSecondaries );
172
173  G4double globalTime = track.GetGlobalTime();
174  G4ThreeVector position = track.GetPosition();
175  // Store nuclear cascade
176  if(captureResult) {
177    G4int n = captureResult->size();
178    for ( G4int isec = 0; isec < n; isec++ ) {
179      G4ReactionProduct* aParticle = captureResult->operator[](isec);
180      G4DynamicParticle * aNewParticle = new G4DynamicParticle();
181      aNewParticle->SetDefinition( aParticle->GetDefinition() );
182      G4LorentzVector itV(aParticle->GetTotalEnergy(), aParticle->GetMomentum());
183      aNewParticle->SetMomentum(itV.vect());
184      G4double localtime = globalTime + tDelay + aParticle->GetTOF();
185      G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position);
186      aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
187      aParticleChange.AddSecondary( aNewTrack );
188      delete aParticle;
189    }
190    delete captureResult;
191  }
192 
193  // Store electromagnetic cascade
194
195  if(nEmSecondaries > 0) {
196
197    for ( G4int isec = 0; isec < nEmSecondaries; isec++ ) {
198      G4ParticleDefinition* pd = Cascade[isec].GetParticleDef();
199      G4double localtime = globalTime;
200      if(isec >= nCascade) localtime += tDelay;
201      if(pd) {
202        G4DynamicParticle* aNewParticle = new G4DynamicParticle;
203        aNewParticle->SetDefinition( pd );
204        aNewParticle->SetMomentum( Cascade[isec].GetMomentum() );
205
206        G4Track* aNewTrack = new G4Track( aNewParticle, localtime, position );
207        aNewTrack->SetTouchableHandle(track.GetTouchableHandle());
208        aParticleChange.AddSecondary( aNewTrack );
209      }
210    }
211  }
212
213  aParticleChange.ProposeLocalEnergyDeposit(0.0);
214  aParticleChange.ProposeTrackStatus(fStopAndKill); 
215
216  return &aParticleChange;
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
220
221G4ReactionProductVector* G4MuonMinusCaptureAtRest::DoMuCapture()
222{
223  G4double mumass = G4MuonMinus::MuonMinus()->GetPDGMass();
224  G4double muBindingEnergy = pEMCascade->GetKShellEnergy(targetZ);
225  /*
226  G4cout << "G4MuonMinusCaptureAtRest::DoMuCapture called Emu= "
227         << muBindingEnergy << G4endl;
228  */
229  // Energy on K-shell
230  G4double muEnergy = mumass + muBindingEnergy;
231  G4double muMom = std::sqrt(muBindingEnergy*(muBindingEnergy + 2.0*mumass));
232  G4double availableEnergy = targetMass + mumass - muBindingEnergy;
233  G4LorentzVector momInitial(0.0,0.0,0.0,availableEnergy);
234  G4LorentzVector momResidual;
235
236  G4ThreeVector vmu = muMom*pEMCascade->GetRandomVec();
237  G4LorentzVector aMuMom(vmu, muEnergy);
238
239  G4double residualMass = 
240    G4NucleiProperties::GetNuclearMass(targetA, targetZ - 1.0);
241
242  G4ReactionProductVector* aPreResult = 0;
243  G4ReactionProduct* aNu = new G4ReactionProduct();
244  aNu->SetDefinition( G4NeutrinoMu::NeutrinoMu() );
245
246  G4int iz = G4int(targetZ);
247  G4int ia = G4int(targetA);
248
249  // proton as a target
250  if(iz <= 2) {
251
252    if(ia > 1) {
253      if(iz == 1 && ia == 2) { 
254        availableEnergy -= neutron_mass_c2;
255      } else if(iz == 1 && ia == 3) {
256        availableEnergy -= 2.0*neutron_mass_c2;
257      } else if(iz == 2) {
258        G4ParticleDefinition* pd = 0;
259        if(ia == 3) pd = G4Deuteron::Deuteron();
260        if(ia == 4) pd = G4Triton::Triton();
261        else 
262          pd = G4ParticleTable::GetParticleTable()->FindIon(1,ia-1,0,1);
263
264        //      G4cout << "Extra " << pd->GetParticleName() << G4endl;
265        availableEnergy -= pd->GetPDGMass();
266      }
267    }
268    //
269    //  Computation in assumption of CM collision of mu and nucleaon
270    // 
271    G4double Enu  = 0.5*(availableEnergy - 
272                         neutron_mass_c2*neutron_mass_c2/availableEnergy);
273
274    // make the nu, and transform to lab;
275    G4ThreeVector nu3Mom = Enu*pEMCascade->GetRandomVec();
276
277    G4ReactionProduct* aN = new G4ReactionProduct();
278    aN->SetDefinition( G4Neutron::Neutron() );
279    aN->SetTotalEnergy( availableEnergy - Enu );
280    aN->SetMomentum( -nu3Mom );
281
282    aNu->SetTotalEnergy( Enu );
283    aNu->SetMomentum( nu3Mom );
284    aPreResult = new G4ReactionProductVector();
285
286    aPreResult->push_back(aN ); 
287    aPreResult->push_back(aNu);
288
289    if(verboseLevel > 1)
290      G4cout << "DoMuCapture on H or He" 
291             <<" EkinN(MeV)= " << (availableEnergy - Enu - neutron_mass_c2)/MeV
292             <<" Enu(MeV)= "<<aNu->GetTotalEnergy()/MeV
293             <<" n= " << aPreResult->size()
294             <<G4endl;
295
296    return aPreResult;
297  }
298
299  // pick random proton inside nucleus
300  G4double eEx;
301  do {
302    theN->Init(targetA, targetZ); 
303    G4LorentzVector thePMom;
304    G4Nucleon * aNucleon = 0;
305    G4int theProtonCounter = G4int( targetZ * G4UniformRand() );
306    G4int counter = 0;
307    theN->StartLoop();
308
309    while( (aNucleon=theN->GetNextNucleon()) ) {
310
311      if( aNucleon->GetDefinition() == G4Proton::Proton() ) {
312        counter++;
313        if(counter == theProtonCounter) {
314          thePMom  = aNucleon->GetMomentum();
315          break;
316        }
317      }
318    }
319
320    // Get the nu momentum in the CMS
321    G4LorentzVector theCMS = thePMom + aMuMom;
322    G4ThreeVector bst = theCMS.boostVector();
323
324    G4double Ecms = theCMS.mag();
325    G4double Enu  = 0.5*(Ecms - neutron_mass_c2*neutron_mass_c2/Ecms);
326    eEx = 0.0;
327
328    if(Enu > 0.0) {
329      // make the nu, and transform to lab;
330      G4ThreeVector nu3Mom = Enu*pEMCascade->GetRandomVec();
331      G4LorentzVector nuMom(nu3Mom, Enu);
332
333      // nu in lab.
334      nuMom.boost(bst);
335      aNu->SetTotalEnergy( nuMom.e() );
336      aNu->SetMomentum( nuMom.vect() );
337   
338      // make residual
339      momResidual = momInitial - nuMom;
340
341      // Call pre-compound on the rest.
342      eEx = momResidual.mag();
343      if(verboseLevel > 1)
344        G4cout << "G4MuonMinusCaptureAtRest::DoMuCapture: " 
345               << " Eex(MeV)= " << (eEx-residualMass)/MeV
346               << " Enu(MeV)= "<<aNu->GetTotalEnergy()/MeV
347               <<G4endl;
348    }
349  } while(eEx <= residualMass);
350
351//  G4cout << "muonCapture : " << eEx << " " << residualMass
352//         << " A,Z= " << targetA << ", "<< targetZ
353//       << "  " << G4int(targetA) << ", " << G4int(targetZ) << G4endl;
354
355  //
356  // Start Deexcitation
357  //
358  G4ThreeVector fromBreit = momResidual.boostVector();
359  G4LorentzVector fscm(0.0,0.0,0.0, eEx);
360  G4Fragment anInitialState;
361  anInitialState.SetA(targetA);
362  anInitialState.SetZ(G4double(iz - 1));
363  anInitialState.SetNumberOfParticles(2);
364  anInitialState.SetNumberOfCharged(0);
365  anInitialState.SetNumberOfHoles(1);
366  anInitialState.SetMomentum(fscm);
367  aPreResult = theHandler->BreakItUp(anInitialState);
368
369  G4ReactionProductVector::iterator ires;
370  G4double eBal = availableEnergy - aNu->GetTotalEnergy();
371  for(ires=aPreResult->begin(); ires!=aPreResult->end(); ires++) {
372    G4LorentzVector itV((*ires)->GetTotalEnergy(), (*ires)->GetMomentum());
373    itV.boost(fromBreit);
374    (*ires)->SetTotalEnergy(itV.t());
375    (*ires)->SetMomentum(itV.vect());
376    eBal -= itV.t();
377  }
378  //
379  // fill neutrino into result
380  //
381  aPreResult->push_back(aNu);
382 
383  if(verboseLevel > 1)
384    G4cout << "DoMuCapture:  Nsec= " 
385           << aPreResult->size() << " Ebalance(MeV)= " << eBal/MeV
386           <<" E0(MeV)= " <<availableEnergy/MeV
387           <<" Mres(GeV)= " <<residualMass/GeV
388           <<G4endl;
389
390  return aPreResult;
391} 
392
Note: See TracBrowser for help on using the repository browser.