Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QCaptureAtRest.hh
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// ---------------- G4QCaptureAtRest header ----------------
29// by Mikhail Kossov, December 2003.
30// Header of G4QCaptureAtRest class of the CHIPS Simulation Branch in GEANT4
31// -------------------------------------------------------------------------------
32// At present (May 2009) only pi-, K- and antiNucleon capture are tested, which
33// are the most crucial for the in matter simulation. The hyperon capture (Sigma-,
34// Xi-, Omega-, antiSigma+) is implemented, but not tested and it is not clear how
35// frequently this kind of interaction takes place in the simulation of the hadronic
36// showers. The antiNeutron Capture At Rest is implemented by this G4QCaptureAtRest
37// class, but it is not clear how the anti-neutrons are stopped in Geant4 tracking.
38// It can be stopped only by interactions with electrons, as the annihilation cross
39// section is huge and any interaction with nucleus results in annihilation. The
40// mu- & tau- Capture At Rest (mu-,nu) & (mu-,nu) are weak processes, which must
41// be simulated together with the reversed Betha decay (e-,nu). While mu- capture is
42// similar to the pi- capture from the nuclear fragmentation point of view (the energy
43// scale is shrinked because m_mu < m_pi and a part of the energy is lost because of
44// the neutrino radiation), the time scale of the mu- capture process is not exact,
45// but it is clear, that it is well delayed. By this reason the mu- capture can be
46// excluded from the G4QCaptureAtRest and can be implemented in the "LongLivingDecay"
47// branch of simulation, which includes excited states of nuclei and short living
48// isotopes. On the "Fast Simulation" Level all radioactive isotopes, long living
49// nuclear excitations, mu-atoms etc, which can be important for the background
50// signals, must be collected in the continuous database and simulated separately.
51// CHIPS is SU(3) event generator, so it does not include reactions with the heavy
52// (c,b,t) quarks involved such as antiDs-, which can be simulated only by SU(6)
53// QUIPS (QUark Invariant Phase Space) model. - May 2009, M.Kossov.-
54// -------------------------------------------------------------------------------
55// All algorithms are similar: the captured particle is absorbed by a nuclear cluster
56// with the subsequent Quark Exchange nuclear fragmentation. The Anti-Proton (antiSigma+)
57// Capture algorithm is more complicated: the anti-baryon annihilates with the quasyfree
58// nucleons on the nuclear periphery. The peripheral interaction results in a number
59// of mesons. A part of them misses the nucleus and comes directly to the output,
60// while others create Multy Quasmon Excitation in the nucleus with the subsequent
61// Quark Excange Fragmentation of the nucleus. At present the two step mechanism of
62// the antiProton-Nucleus interaction is hardwired in the G4QEnvironment class, but
63// with time the first step of the interaction can be moved to this G4QCaptureAtRest
64// class, to make the G4QEnvirement class simpler and better defined. This is
65// necessary because the G4QEnvironment class is going to loos the previlage of
66// the CHIPS Head Class (as previously the G4Quasmon class lost it) and G4QCollision
67// class is going to be the CHIPS Head Class, where a few Nuclear Environments can
68// exist (e.g. the Nuclear Environment of the Projectile Nucleus and the Nuclear
69// Environment of the Target Nucleus). By the way, the antiProton-H1 interaction At
70// Rest (CHIPSI) can be still simulated with only the G4Quasmon class, as this
71// reaction does not have any nuclear environment.- May 2009, Mikhail Kossov.-
72// --------------------------------------------------------------------------------
73// ****************************************************************************************
74// This Header is a part of the CHIPS physics package (author: M. Kosov)
75// ****************************************************************************************
76// Short Description: This is a universal process for nuclear capture
77// (including annihilation) of all negative particles (cold neutrons, negative
78// hadrons, negative leptons: mu- & tau-). It can be used for the cold neutron
79// capture, but somebody should decide what is the probability (defined
80// by the capture cross-section and atomic material properties) to switch
81// the cold neutron to the at-rest neutron. - M.K. 2009.
82// ----------------------------------------------------------------------
83
84#ifndef G4QCaptureAtRest_hh
85#define G4QCaptureAtRest_hh
86
87// GEANT4 Headers
88#include "globals.hh"
89#include "G4ios.hh"
90#include "G4VRestProcess.hh"
91#include "G4ParticleTypes.hh"
92#include "G4VParticleChange.hh"
94#include "G4DynamicParticle.hh"
95#include "Randomize.hh"
96#include "G4ThreeVector.hh"
97#include "G4LorentzVector.hh"
98#include "G4RandomDirection.hh"
99
100// CHIPS Headers
101#include "G4QEnvironment.hh"
102#include "G4QIsotope.hh"
103#include "G4QPDGToG4Particle.hh"
104
106{
107private:
108
109 // Hide assignment operator as private
110 G4QCaptureAtRest& operator=(const G4QCaptureAtRest &right);
111
112 // Copy constructor
114
115public:
116
117 // Constructor
118 G4QCaptureAtRest(const G4String& processName ="CHIPSNuclearCaptureAtRest");
119
120 // Destructor
121 virtual ~G4QCaptureAtRest();
122
123 virtual G4bool IsApplicable(const G4ParticleDefinition& particle);
124
125 G4VParticleChange* AtRestDoIt(const G4Track& aTrack, const G4Step& aStep);
126
128
130
131 // Static functions
132 static void SetManual();
133 static void SetStandard();
134 static void SetParameters(G4double temper=180., G4double ssin2g=.1, G4double etaetap=.3,
135 G4double fN=0., G4double fD=0., G4double cP=1., G4double mR=1.,
136 G4int npCHIPSWorld=234, G4double solAn=.5, G4bool efFlag=false,
137 G4double piTh=141.4,G4double mpi2=20000.,G4double dinum=1880.);
138
139protected:
140
141 // zero mean lifetime
143 G4double RandomizeDecayElectron(G4int Z); // Randomize energy of decay electron (in MeV)
144private:
145
146 G4bool RandomizeMuDecayOrCapture(G4int Z, G4int N); // true=MuCapture, false=MuDecay
147 void CalculateEnergyDepositionOfMuCapture(G4int Z); // (2p->1s, MeV) @@ Now N-independent
148 G4bool RandomizeTauDecayOrCapture(G4int Z, G4int N);// true=TauCapture, false=TauDecay
149 void CalculateEnergyDepositionOfTauCapture(G4int Z);// (2p->1s, MeV) @@N-independ,Improve
150
151// BODY
152private:
153 // Static Parameters
154 static G4bool manualFlag; // If false then standard parameters are used
155 static G4int nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
156 // -> Parameters of the G4Quasmon class:
157 static G4double Temperature; // Quasmon Temperature
158 static G4double SSin2Gluons; // Percent of ssbar sea in a constituen gluon
159 static G4double EtaEtaprime; // Part of eta-prime in all etas
160 // -> Parameters of the G4QNucleus class:
161 static G4double freeNuc; // probability of the quasi-free baryon on surface
162 static G4double freeDib; // probability of the quasi-free dibaryon on surface
163 static G4double clustProb; // clusterization probability in dense region
164 static G4double mediRatio; // relative vacuum hadronization probability
165 // -> Parameters of the G4QEnvironment class:
166 static G4bool EnergyFlux; // Flag for Energy Flux use instead of Multy Quasmon
167 static G4double SolidAngle; // Part of Solid Angle to capture secondaries(@@A-dep)
168 static G4double PiPrThresh; // Pion Production Threshold for gammas
169 static G4double M2ShiftVir; // Shift for M2=-Q2=m_pi^2 of the virtual gamma
170 static G4double DiNuclMass; // Double Nucleon Mass for virtual normalization
171 //
172 // Working parameters
173 G4LorentzVector EnMomConservation; // Residual of Energy/Momentum Cons.
174 G4int nOfNeutrons; // #of neutrons in the target nucleus
175 // Modifires for the reaction
176 G4double Time; // Time shift of the capture reaction
177 G4double EnergyDeposition; // Energy deposited in the reaction
178
179};
180#endif
G4ForceCondition
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4int GetNumberOfNeutronsInTarget()
G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
virtual ~G4QCaptureAtRest()
static void SetStandard()
static void SetParameters(G4double temper=180., G4double ssin2g=.1, G4double etaetap=.3, G4double fN=0., G4double fD=0., G4double cP=1., G4double mR=1., G4int npCHIPSWorld=234, G4double solAn=.5, G4bool efFlag=false, G4double piTh=141.4, G4double mpi2=20000., G4double dinum=1880.)
G4LorentzVector GetEnegryMomentumConservation()
G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *)
virtual G4bool IsApplicable(const G4ParticleDefinition &particle)
static void SetManual()
G4double RandomizeDecayElectron(G4int Z)
Definition: G4Step.hh:78