Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AnnihiToMuPair Class Reference

#include <G4AnnihiToMuPair.hh>

+ Inheritance diagram for G4AnnihiToMuPair:

Public Member Functions

 G4AnnihiToMuPair (const G4String &processName="AnnihiToMuPair", G4ProcessType type=fElectromagnetic)
 
 ~G4AnnihiToMuPair () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void PrintInfoDefinition ()
 
void SetCrossSecFactor (G4double fac)
 
G4double GetCrossSecFactor ()
 
G4double CrossSectionPerVolume (G4double positronEnergy, const G4Material *)
 
G4double ComputeCrossSectionPerElectron (const G4double energy)
 
G4double ComputeCrossSectionPerAtom (const G4double energy, const G4double Z)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
G4AnnihiToMuPairoperator= (const G4AnnihiToMuPair &right)=delete
 
 G4AnnihiToMuPair (const G4AnnihiToMuPair &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 57 of file G4AnnihiToMuPair.hh.

Constructor & Destructor Documentation

◆ G4AnnihiToMuPair() [1/2]

G4AnnihiToMuPair::G4AnnihiToMuPair ( const G4String processName = "AnnihiToMuPair",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 62 of file G4AnnihiToMuPair.cc.

63 :G4VDiscreteProcess (processName, type)
64{
65 //e+ Energy threshold
66 if(processName == "AnnihiToTauPair") {
68 part1 = G4TauPlus::TauPlus();
69 part2 = G4TauMinus::TauMinus();
70 fInfo = "e+e->tau+tau-";
71 } else {
73 part1 = G4MuonPlus::MuonPlus();
74 part2 = G4MuonMinus::MuonMinus();
75 }
76 fMass = part1->GetPDGMass();
77 fLowEnergyLimit =
78 2.*fMass*fMass/CLHEP::electron_mass_c2 - CLHEP::electron_mass_c2;
79
80 //model is ok up to 1000 TeV due to neglected Z-interference
81 fHighEnergyLimit = 1000.*TeV;
82
83 fCurrentSigma = 0.0;
84 fCrossSecFactor = 1.;
86 fManager->Register(this);
87}
@ fAnnihilationToTauTau
@ fAnnihilationToMuMu
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:134
static G4TauPlus * TauPlus()
Definition: G4TauPlus.cc:133
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4AnnihiToMuPair()

G4AnnihiToMuPair::~G4AnnihiToMuPair ( )
override

Definition at line 91 of file G4AnnihiToMuPair.cc.

92{
93 fManager->DeRegister(this);
94}
void DeRegister(G4VEnergyLossProcess *p)

◆ G4AnnihiToMuPair() [2/2]

G4AnnihiToMuPair::G4AnnihiToMuPair ( const G4AnnihiToMuPair )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4AnnihiToMuPair::BuildPhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 105 of file G4AnnihiToMuPair.cc.

106{
108}

◆ ComputeCrossSectionPerAtom()

G4double G4AnnihiToMuPair::ComputeCrossSectionPerAtom ( const G4double  energy,
const G4double  Z 
)

Definition at line 152 of file G4AnnihiToMuPair.cc.

154{
155 return ComputeCrossSectionPerElectron(energy)*Z;
156}
const G4int Z[17]
G4double ComputeCrossSectionPerElectron(const G4double energy)

◆ ComputeCrossSectionPerElectron()

G4double G4AnnihiToMuPair::ComputeCrossSectionPerElectron ( const G4double  energy)

Definition at line 122 of file G4AnnihiToMuPair.cc.

125{
126 G4double rmuon = CLHEP::elm_coupling/fMass; //classical particle radius
127 G4double sig0 = CLHEP::pi*rmuon*rmuon/3.; //constant in crossSection
128 const G4double pial = CLHEP::pi*CLHEP::fine_structure_const; // pi * alphaQED
129
130 if (e <= fLowEnergyLimit) return 0.0;
131
132 const G4double xi = fLowEnergyLimit/e;
133 const G4double piaxi = pial * std::sqrt(xi);
134 G4double sigma = sig0 * xi * (1. + xi*0.5);
135 //G4cout << "### xi= " << xi << " piaxi=" << piaxi << G4endl;
136
137 // argument of the exponent below 0.1 or above 10
138 // Sigma per electron * number of electrons per atom
139 if(xi <= 1.0 - 100*piaxi*piaxi) {
140 sigma *= std::sqrt(1.0 - xi);
141 } else if( xi >= 1.0 - 0.01*piaxi*piaxi) {
142 sigma *= piaxi;
143 } else {
144 sigma *= piaxi/(1. - G4Exp( -piaxi/std::sqrt(1-xi) ));
145 }
146 //G4cout << "### sigma= " << sigma << G4endl;
147 return sigma;
148}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
double G4double
Definition: G4Types.hh:83

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ CrossSectionPerVolume()

G4double G4AnnihiToMuPair::CrossSectionPerVolume ( G4double  positronEnergy,
const G4Material aMaterial 
)

Definition at line 160 of file G4AnnihiToMuPair.cc.

162{
164}
G4double GetTotNbOfElectPerVolume() const
Definition: G4Material.hh:207

Referenced by GetMeanFreePath(), and PostStepDoIt().

◆ GetCrossSecFactor()

G4double G4AnnihiToMuPair::GetCrossSecFactor ( )
inline

Definition at line 81 of file G4AnnihiToMuPair.hh.

81{return fCrossSecFactor;};

◆ GetMeanFreePath()

G4double G4AnnihiToMuPair::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition  
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 168 of file G4AnnihiToMuPair.cc.

171{
172 const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
173 G4double energy = aDynamicPositron->GetTotalEnergy();
174 const G4Material* aMaterial = aTrack.GetMaterial();
175
176 // cross section before step
177 fCurrentSigma = CrossSectionPerVolume(energy, aMaterial);
178
179 // increase the CrossSection by CrossSecFactor (default 1)
180 return (fCurrentSigma > 0.0) ? 1.0/(fCurrentSigma*fCrossSecFactor) : DBL_MAX;
181}
G4double CrossSectionPerVolume(G4double positronEnergy, const G4Material *)
G4double GetTotalEnergy() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double energy(const ThreeVector &p, const G4double m)
#define DBL_MAX
Definition: templates.hh:62

◆ IsApplicable()

G4bool G4AnnihiToMuPair::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 98 of file G4AnnihiToMuPair.cc.

99{
100 return ( &particle == G4Positron::Positron() );
101}
static G4Positron * Positron()
Definition: G4Positron.cc:93

◆ operator=()

G4AnnihiToMuPair & G4AnnihiToMuPair::operator= ( const G4AnnihiToMuPair right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4AnnihiToMuPair::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 185 of file G4AnnihiToMuPair.cc.

190{
192
193 // current Positron energy and direction, return if energy too low
194 const G4DynamicParticle *aDynamicPositron = aTrack.GetDynamicParticle();
195 const G4double Mele = CLHEP::electron_mass_c2;
196 G4double Epos = aDynamicPositron->GetTotalEnergy();
197 G4double xs = CrossSectionPerVolume(Epos, aTrack.GetMaterial());
198
199 // test of cross section
200 if(xs > 0.0 && fCurrentSigma*G4UniformRand() > xs)
201 {
202 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
203 }
204
205 const G4ThreeVector PosiDirection = aDynamicPositron->GetMomentumDirection();
206 G4double xi = fLowEnergyLimit/Epos; // xi is always less than 1,
207 // goes to 0 at high Epos
208
209 // generate cost; probability function 1+cost**2 at high Epos
210 //
211 G4double cost;
212 do { cost = 2.*G4UniformRand()-1.; }
213 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
214 while (2.*G4UniformRand() > 1.+xi+cost*cost*(1.-xi) );
215 G4double sint = sqrt(1.-cost*cost);
216
217 // generate phi
218 //
219 G4double phi = 2.*CLHEP::pi*G4UniformRand();
220
221 G4double Ecm = std::sqrt(0.5*Mele*(Epos+Mele));
222 G4double Pcm = std::sqrt(Ecm*Ecm - fMass*fMass);
223 G4double beta = std::sqrt((Epos-Mele)/(Epos+Mele));
224 G4double gamma = Ecm/Mele;
225 G4double Pt = Pcm*sint;
226
227 // energy and momentum of the muons in the Lab
228 //
229 G4double EmuPlus = gamma*(Ecm + cost*beta*Pcm);
230 G4double EmuMinus = gamma*(Ecm - cost*beta*Pcm);
231 G4double PmuPlusZ = gamma*(beta*Ecm + cost*Pcm);
232 G4double PmuMinusZ = gamma*(beta*Ecm - cost*Pcm);
233 G4double PmuPlusX = Pt*std::cos(phi);
234 G4double PmuPlusY = Pt*std::sin(phi);
235 G4double PmuMinusX =-PmuPlusX;
236 G4double PmuMinusY =-PmuPlusY;
237 // absolute momenta
238 G4double PmuPlus = std::sqrt(Pt*Pt+PmuPlusZ *PmuPlusZ );
239 G4double PmuMinus = std::sqrt(Pt*Pt+PmuMinusZ*PmuMinusZ);
240
241 // mu+ mu- directions for Positron in z-direction
242 //
244 MuPlusDirection(PmuPlusX/PmuPlus, PmuPlusY/PmuPlus, PmuPlusZ/PmuPlus);
246 MuMinusDirection(PmuMinusX/PmuMinus,PmuMinusY/PmuMinus,PmuMinusZ/PmuMinus);
247
248 // rotate to actual Positron direction
249 //
250 MuPlusDirection.rotateUz(PosiDirection);
251 MuMinusDirection.rotateUz(PosiDirection);
252
254
255 // create G4DynamicParticle object for the particle1
256 G4DynamicParticle* aParticle1 =
257 new G4DynamicParticle(part1, MuPlusDirection, EmuPlus-fMass);
258 aParticleChange.AddSecondary(aParticle1);
259 // create G4DynamicParticle object for the particle2
260 G4DynamicParticle* aParticle2 =
261 new G4DynamicParticle(part2, MuMinusDirection, EmuMinus-fMass);
262 aParticleChange.AddSecondary(aParticle2);
263
264 // Kill the incident positron
265 //
268
269 return &aParticleChange;
270}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331

◆ PrintInfoDefinition()

void G4AnnihiToMuPair::PrintInfoDefinition ( )

Definition at line 274 of file G4AnnihiToMuPair.cc.

275{
276 G4String comments = fInfo + " annihilation, atomic e- at rest, SubType=";
277 G4cout << G4endl << GetProcessName() << ": " << comments
279 G4cout << " threshold at " << fLowEnergyLimit/CLHEP::GeV << " GeV"
280 << " good description up to "
281 << fHighEnergyLimit/CLHEP::TeV << " TeV for all Z." << G4endl;
282}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

Referenced by BuildPhysicsTable().

◆ SetCrossSecFactor()

void G4AnnihiToMuPair::SetCrossSecFactor ( G4double  fac)

Definition at line 112 of file G4AnnihiToMuPair.cc.

114{
115 fCrossSecFactor = fac;
116 //G4cout << "The cross section for AnnihiToMuPair is artificially "
117 // << "increased by the CrossSecFactor=" << fCrossSecFactor << G4endl;
118}

Referenced by G4EmExtraPhysics::ConstructProcess().


The documentation for this class was generated from the following files: