Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuMinusCapturePrecompound.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$
27//
28//-----------------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32// File name: G4MuMinusCapturePrecompound
33//
34// Author: V.Ivanchenko ([email protected])
35//
36// Creation date: 22 April 2012 on base of G4MuMinusCaptureCascade
37//
38//
39//-----------------------------------------------------------------------------
40//
41// Modifications:
42//
43//-----------------------------------------------------------------------------
44
46#include "Randomize.hh"
47#include "G4RandomDirection.hh"
49#include "G4SystemOfUnits.hh"
50#include "G4MuonMinus.hh"
51#include "G4NeutrinoMu.hh"
52#include "G4Neutron.hh"
53#include "G4Proton.hh"
54#include "G4Triton.hh"
55#include "G4LorentzVector.hh"
57#include "G4NucleiProperties.hh"
59#include "G4PreCompoundModel.hh"
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
66 : G4HadronicInteraction("muMinusNuclearCapture")
67{
68 fMuMass = G4MuonMinus::MuonMinus()->GetPDGMass();
69 fProton = G4Proton::Proton();
70 fNeutron = G4Neutron::Neutron();
71 fThreshold = 10*MeV;
72 fPreCompound = ptr;
73 if(!ptr) {
76 ptr = static_cast<G4VPreCompoundModel*>(p);
77 fPreCompound = ptr;
78 if(!ptr) { fPreCompound = new G4PreCompoundModel(); }
79 }
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
85{
86 result.Clear();
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90
93 G4Nucleus& targetNucleus)
94{
95 result.Clear();
97 fTime = projectile.GetGlobalTime();
98 G4double time0 = fTime;
99
100 G4double muBindingEnergy = projectile.GetBoundEnergy();
101
102 G4int Z = targetNucleus.GetZ_asInt();
103 G4int A = targetNucleus.GetA_asInt();
105
106 /*
107 G4cout << "G4MuMinusCapturePrecompound::ApplyYourself: Emu= "
108 << muBindingEnergy << G4endl;
109 */
110 // Energy on K-shell
111 G4double muEnergy = fMuMass + muBindingEnergy;
112 G4double muMom = std::sqrt(muBindingEnergy*(muBindingEnergy + 2.0*fMuMass));
113 G4double availableEnergy = massA + fMuMass - muBindingEnergy;
114 G4double residualMass = G4NucleiProperties::GetNuclearMass(A, Z - 1);
115
116 G4ThreeVector vmu = muMom*G4RandomDirection();
117 G4LorentzVector aMuMom(vmu, muEnergy);
118
119 // p or 3He as a target
120 // two body reaction mu- + A(Z,A) -> nuMu + A(Z-1,A)
121 if((1 == Z && 1 == A) || (2 == Z && 3 == A)) {
122
123 G4ParticleDefinition* pd = 0;
124 if(1 == Z) { pd = fNeutron; }
125 else { pd = G4Triton::Triton(); }
126
127 //
128 // Computation in assumption of CM reaction
129 //
130 G4double e = 0.5*(availableEnergy -
131 residualMass*residualMass/availableEnergy);
132
134 AddNewParticle(G4NeutrinoMu::NeutrinoMu(), nudir, e);
135 nudir *= -1.0;
136 AddNewParticle(pd, nudir, availableEnergy - e - residualMass);
137
138
139 } else {
140 // sample mu- + p -> nuMu + n reaction in CM of muonic atom
141
142 // muon
143//
144// NOTE by K.Genser and J.Yarba:
145// The code below isn't working because emu always turns smaller than fMuMass
146// For this reason the sqrt is producing a NaN
147//
148// G4double emu = (availableEnergy*availableEnergy - massA*massA
149// + fMuMass*fMuMass)/(2*availableEnergy);
150// G4ThreeVector mudir = G4RandomDirection();
151// G4LorentzVector momMuon(std::sqrt(emu*emu - fMuMass*fMuMass)*mudir, emu);
152
153 // nucleus
154 G4LorentzVector momInitial(0.0,0.0,0.0,availableEnergy);
155 G4LorentzVector momResidual, momNu;
156
157 // pick random proton inside nucleus
158 G4double eEx;
159 fNucleus.Init(A, Z);
160 const std::vector<G4Nucleon>& nucleons= fNucleus.GetNucleons();
162
163 G4int nneutrons = 1;
164 G4int reentryCount = 0;
165
166 do {
167 ++reentryCount;
168 G4int index = 0;
169 do {
170 index=G4int(A*G4UniformRand());
171 pDef = nucleons[index].GetDefinition();
172 } while(pDef != fProton);
173 G4LorentzVector momP = nucleons[index].Get4Momentum();
174
175 // Get CMS kinematics
176 G4LorentzVector theCMS = momP + aMuMom;
177 G4ThreeVector bst = theCMS.boostVector();
178
179 G4double Ecms = theCMS.mag();
180 G4double Enu = 0.5*(Ecms - neutron_mass_c2*neutron_mass_c2/Ecms);
181 eEx = 0.0;
182
183 if(Enu > 0.0) {
184 // make the nu, and transform to lab;
185 momNu.set(Enu*G4RandomDirection(), Enu);
186
187 // nu in lab.
188 momNu.boost(bst);
189 momResidual = momInitial - momNu;
190 eEx = momResidual.mag() - residualMass;
191
192 // release neutron
193
194 if(eEx > 0.0) {
195 G4double eth = residualMass - massA + fThreshold + 2*neutron_mass_c2;
196 if(Ecms - Enu > eth) {
197 theCMS -= momNu;
198 G4double ekin = theCMS.e() - eth;
199 G4ThreeVector dir = theCMS.vect().unit();
200 AddNewParticle(fNeutron, dir, ekin);
201 momResidual -=
202 result.GetSecondary(0)->GetParticle()->Get4Momentum();
203 --Z;
204 --A;
205 residualMass = G4NucleiProperties::GetNuclearMass(A, Z);
206 nneutrons = 0;
207 }
208 }
209 }
210 if(Enu <= 0.0 && eEx <= 0.0 && reentryCount > 100) {
212 ed << "Call for " << GetModelName() << G4endl;
213 ed << "Target Z= " << Z
214 << " A= " << A << G4endl;
215 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
216 G4Exception("G4MuMinusCapturePrecompound::AtRestDoIt", "had006",
217 FatalException, ed);
218 }
219 } while(eEx <= 0.0);
220
221 G4ThreeVector dir = momNu.vect().unit();
222 AddNewParticle(G4NeutrinoMu::NeutrinoMu(), dir, momNu.e());
223
224 G4Fragment initialState(A, Z, momResidual);
225 initialState.SetNumberOfExcitedParticle(nneutrons,0);
226 initialState.SetNumberOfHoles(1,1);
227
228 // decay time for pre-compound/de-excitation starts from zero
229 G4ReactionProductVector* rpv = fPreCompound->DeExcite(initialState);
230 size_t n = rpv->size();
231 for(size_t i=0; i<n; ++i) {
232 G4ReactionProduct* rp = (*rpv)[i];
233
234 // reaction time
235 fTime = time0 + rp->GetTOF();
236 G4ThreeVector direction = rp->GetMomentum().unit();
237 AddNewParticle(rp->GetDefinition(), direction, rp->GetKineticEnergy());
238 delete rp;
239 }
240 delete rpv;
241 }
242 if(verboseLevel > 1)
243 G4cout << "G4MuMinusCapturePrecompound::ApplyYourself: Nsec= "
244 << result.GetNumberOfSecondaries()
245 <<" E0(MeV)= " <<availableEnergy/MeV
246 <<" Mres(GeV)= " <<residualMass/GeV
247 <<G4endl;
248
249 return &result;
250}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
253
254void G4MuMinusCapturePrecompound::ModelDescription(std::ostream& outFile) const
255{
256 outFile << "Sampling of mu- capture by atomic nucleus from K-shell"
257 << " mesoatom orbit.\n"
258 << "Primary reaction mu- + p -> n + neutrino, neutron providing\n"
259 << " initial excitation of the target nucleus and PreCompound"
260 << " model samples final state\n";
261}
262
263//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ FatalException
@ stopAndKill
G4ThreeVector G4RandomDirection()
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
void set(double x, double y, double z, double t)
G4LorentzVector Get4Momentum() const
const std::vector< G4Nucleon > & GetNucleons()
void Init(G4int theA, G4int theZ)
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:335
void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
Definition: G4Fragment.hh:316
void SetStatusChange(G4HadFinalStateStatus aS)
G4int GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
G4double GetBoundEnergy() const
G4double GetGlobalTime() const
G4DynamicParticle * GetParticle()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
const G4String & GetModelName() const
G4MuMinusCapturePrecompound(G4VPreCompoundModel *ptr=0)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void ModelDescription(std::ostream &outFile) const
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4NeutrinoMu * NeutrinoMu()
Definition: G4NeutrinoMu.cc:85
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTOF() const
G4ParticleDefinition * GetDefinition() const
static G4Triton * Triton()
Definition: G4Triton.cc:95
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76