Geant4 11.2.2
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 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
 
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 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)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- 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 59 of file G4AnnihiToMuPair.cc.

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

◆ ~G4AnnihiToMuPair()

G4AnnihiToMuPair::~G4AnnihiToMuPair ( )
override

Definition at line 87 of file G4AnnihiToMuPair.cc.

88{
89 fManager->DeRegister(this);
90}
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 101 of file G4AnnihiToMuPair.cc.

102{
104}

◆ ComputeCrossSectionPerAtom()

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

Definition at line 150 of file G4AnnihiToMuPair.cc.

152{
153 return ComputeCrossSectionPerElectron(energy)*Z;
154}
G4double ComputeCrossSectionPerElectron(const G4double energy)

◆ ComputeCrossSectionPerElectron()

G4double G4AnnihiToMuPair::ComputeCrossSectionPerElectron ( const G4double energy)

Definition at line 118 of file G4AnnihiToMuPair.cc.

121{
122 G4double rmuon = CLHEP::elm_coupling/fMass; //classical particle radius
123 G4double sig0 = CLHEP::pi*rmuon*rmuon/3.; //constant in crossSection
124 const G4double pial = CLHEP::pi*CLHEP::fine_structure_const; // pi * alphaQED
125
126 if (e <= fLowEnergyLimit) return 0.0;
127
128 const G4double xi = fLowEnergyLimit/e;
129 const G4double piaxi = pial * std::sqrt(xi);
130 G4double sigma = sig0 * xi * (1. + xi*0.5);
131 //G4cout << "### xi= " << xi << " piaxi=" << piaxi << G4endl;
132
133 // argument of the exponent below 0.1 or above 10
134 // Sigma per electron * number of electrons per atom
135 if(xi <= 1.0 - 100*piaxi*piaxi) {
136 sigma *= std::sqrt(1.0 - xi);
137 }
138 else if (xi >= 1.0 - 0.01 * piaxi * piaxi) {
139 sigma *= piaxi;
140 }
141 else {
142 sigma *= piaxi / (1. - G4Exp(-piaxi / std::sqrt(1 - xi)));
143 }
144 // G4cout << "### sigma= " << sigma << G4endl;
145 return sigma;
146}
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 158 of file G4AnnihiToMuPair.cc.

160{
162}
G4double GetTotNbOfElectPerVolume() const

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 166 of file G4AnnihiToMuPair.cc.

169{
170 const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
171 G4double energy = aDynamicPositron->GetTotalEnergy();
172 const G4Material* aMaterial = aTrack.GetMaterial();
173
174 // cross section before step
175 fCurrentSigma = CrossSectionPerVolume(energy, aMaterial);
176
177 // increase the CrossSection by CrossSecFactor (default 1)
178 return (fCurrentSigma > 0.0) ? 1.0/(fCurrentSigma*fCrossSecFactor) : DBL_MAX;
179}
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 94 of file G4AnnihiToMuPair.cc.

95{
96 return ( &particle == G4Positron::Positron() );
97}
static G4Positron * Positron()
Definition G4Positron.cc:90

◆ 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 183 of file G4AnnihiToMuPair.cc.

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

◆ PrintInfoDefinition()

void G4AnnihiToMuPair::PrintInfoDefinition ( )

Definition at line 267 of file G4AnnihiToMuPair.cc.

268{
269 G4String comments = fInfo + " annihilation, atomic e- at rest, SubType=";
270 G4cout << G4endl << GetProcessName() << ": " << comments << GetProcessSubType() << G4endl;
271 G4cout << " threshold at " << fLowEnergyLimit / CLHEP::GeV << " GeV"
272 << " good description up to " << fHighEnergyLimit / CLHEP::TeV << " TeV for all Z."
273 << G4endl;
274}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4int GetProcessSubType() const
const G4String & GetProcessName() const

Referenced by BuildPhysicsTable().

◆ SetCrossSecFactor()

void G4AnnihiToMuPair::SetCrossSecFactor ( G4double fac)

Definition at line 108 of file G4AnnihiToMuPair.cc.

110{
111 fCrossSecFactor = fac;
112 //G4cout << "The cross section for AnnihiToMuPair is artificially "
113 // << "increased by the CrossSecFactor=" << fCrossSecFactor << G4endl;
114}

Referenced by G4EmExtraPhysics::ConstructProcess().


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