Geant4 10.7.0
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 ()
 
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 ComputeCrossSectionPerAtom (G4double PositronEnergy, G4double AtomicZ)
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
- 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 &)
 
- 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 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 56 of file G4AnnihiToMuPair.hh.

Constructor & Destructor Documentation

◆ G4AnnihiToMuPair()

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 const G4double Mu_massc2 = G4MuonPlus::MuonPlus()->GetPDGMass();
64 LowestEnergyLimit = 2.*Mu_massc2*Mu_massc2/electron_mass_c2 - electron_mass_c2;
65
66 //modele ok up to 1000 TeV due to neglected Z-interference
67 HighestEnergyLimit = 1000.*TeV;
68
69 CurrentSigma = 0.0;
70 CrossSecFactor = 1.;
73}
double G4double
Definition: G4Types.hh:83
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4AnnihiToMuPair()

G4AnnihiToMuPair::~G4AnnihiToMuPair ( )

Definition at line 77 of file G4AnnihiToMuPair.cc.

78{
80}
void DeRegister(G4VEnergyLossProcess *p)

Member Function Documentation

◆ BuildPhysicsTable()

void G4AnnihiToMuPair::BuildPhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4AnnihiToMuPair.cc.

94{
95 CurrentSigma = 0.0;
97}

◆ ComputeCrossSectionPerAtom()

G4double G4AnnihiToMuPair::ComputeCrossSectionPerAtom ( G4double  PositronEnergy,
G4double  AtomicZ 
)

Definition at line 111 of file G4AnnihiToMuPair.cc.

114{
115 static const G4double Mmuon = G4MuonPlus::MuonPlus()->GetPDGMass();
116 static const G4double Rmuon = CLHEP::elm_coupling/Mmuon; //classical particle radius
117 static const G4double Sig0 = CLHEP::pi*Rmuon*Rmuon/3.; //constant in crossSection
118 static const G4double pia = CLHEP::pi * CLHEP::fine_structure_const; // pi * alphaQED
119
120 G4double CrossSection = 0.;
121 if (Epos < LowestEnergyLimit) return CrossSection;
122
123 G4double xi = LowestEnergyLimit/Epos;
124 G4double piaxi = pia * sqrt(xi);
125 G4double SigmaEl = Sig0 * xi * (1.+xi/2.) * piaxi;
126 if( Epos>LowestEnergyLimit+1.e-5 ) SigmaEl /= (1.-std::exp( -piaxi/std::sqrt(1-xi) ));
127 CrossSection = SigmaEl*Z; // SigmaEl per electron * number of electrons per atom
128 return CrossSection;
129}

Referenced by CrossSectionPerVolume().

◆ CrossSectionPerVolume()

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

Definition at line 133 of file G4AnnihiToMuPair.cc.

135{
136 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
137 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
138
139 G4double SIGMA = 0.0;
140
141 for ( size_t i=0 ; i < aMaterial->GetNumberOfElements() ; ++i )
142 {
143 G4double AtomicZ = (*theElementVector)[i]->GetZ();
144 SIGMA += NbOfAtomsPerVolume[i] *
145 ComputeCrossSectionPerAtom(PositronEnergy,AtomicZ);
146 }
147 return SIGMA;
148}
std::vector< G4Element * > G4ElementVector
G4double ComputeCrossSectionPerAtom(G4double PositronEnergy, G4double AtomicZ)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204

Referenced by GetMeanFreePath(), and PostStepDoIt().

◆ GetCrossSecFactor()

G4double G4AnnihiToMuPair::GetCrossSecFactor ( )
inline

Definition at line 80 of file G4AnnihiToMuPair.hh.

80{return CrossSecFactor;};

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 152 of file G4AnnihiToMuPair.cc.

157{
158 const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
159 G4double PositronEnergy = aDynamicPositron->GetKineticEnergy()
160 +electron_mass_c2;
161 G4Material* aMaterial = aTrack.GetMaterial();
162 CurrentSigma = CrossSectionPerVolume(PositronEnergy, aMaterial);
163
164 // increase the CrossSection by CrossSecFactor (default 1)
165 G4double mfp = DBL_MAX;
166 if(CurrentSigma > DBL_MIN) mfp = 1.0/(CurrentSigma*CrossSecFactor);
167
168 return mfp;
169}
G4double CrossSectionPerVolume(G4double PositronEnergy, const G4Material *)
G4double GetKineticEnergy() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62

◆ IsApplicable()

G4bool G4AnnihiToMuPair::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 84 of file G4AnnihiToMuPair.cc.

85{
86 return ( &particle == G4Positron::Positron() );
87}
static G4Positron * Positron()
Definition: G4Positron.cc:93

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 173 of file G4AnnihiToMuPair.cc.

178{
179
181 static const G4double Mele=electron_mass_c2;
182 static const G4double Mmuon=G4MuonPlus::MuonPlus()->GetPDGMass();
183
184 // current Positron energy and direction, return if energy too low
185 const G4DynamicParticle *aDynamicPositron = aTrack.GetDynamicParticle();
186 G4double Epos = aDynamicPositron->GetKineticEnergy() + Mele;
187
188 // test of cross section
189 if(CurrentSigma*G4UniformRand() >
190 CrossSectionPerVolume(Epos, aTrack.GetMaterial()))
191 {
192 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
193 }
194
195 if (Epos < LowestEnergyLimit) {
196 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
197 }
198
199 G4ParticleMomentum PositronDirection =
200 aDynamicPositron->GetMomentumDirection();
201 G4double xi = LowestEnergyLimit/Epos; // xi is always less than 1,
202 // goes to 0 at high Epos
203
204 // generate cost
205 //
206 G4double cost;
207 do { cost = 2.*G4UniformRand()-1.; }
208 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
209 while (2.*G4UniformRand() > 1.+xi+cost*cost*(1.-xi) );
210 //1+cost**2 at high Epos
211 G4double sint = sqrt(1.-cost*cost);
212
213 // generate phi
214 //
215 G4double phi=2.*pi*G4UniformRand();
216
217 G4double Ecm = sqrt(0.5*Mele*(Epos+Mele));
218 G4double Pcm = sqrt(Ecm*Ecm-Mmuon*Mmuon);
219 G4double beta = sqrt((Epos-Mele)/(Epos+Mele));
220 G4double gamma = Ecm/Mele; // =sqrt((Epos+Mele)/(2.*Mele));
221 G4double Pt = Pcm*sint;
222
223 // energy and momentum of the muons in the Lab
224 //
225 G4double EmuPlus = gamma*( Ecm+cost*beta*Pcm);
226 G4double EmuMinus = gamma*( Ecm-cost*beta*Pcm);
227 G4double PmuPlusZ = gamma*(beta*Ecm+cost* Pcm);
228 G4double PmuMinusZ = gamma*(beta*Ecm-cost* Pcm);
229 G4double PmuPlusX = Pt*cos(phi);
230 G4double PmuPlusY = Pt*sin(phi);
231 G4double PmuMinusX =-Pt*cos(phi);
232 G4double PmuMinusY =-Pt*sin(phi);
233 // absolute momenta
234 G4double PmuPlus = sqrt(Pt*Pt+PmuPlusZ *PmuPlusZ );
235 G4double PmuMinus = sqrt(Pt*Pt+PmuMinusZ*PmuMinusZ);
236
237 // mu+ mu- directions for Positron in z-direction
238 //
240 MuPlusDirection ( PmuPlusX/PmuPlus, PmuPlusY/PmuPlus, PmuPlusZ/PmuPlus );
242 MuMinusDirection(PmuMinusX/PmuMinus,PmuMinusY/PmuMinus,PmuMinusZ/PmuMinus);
243
244 // rotate to actual Positron direction
245 //
246 MuPlusDirection.rotateUz(PositronDirection);
247 MuMinusDirection.rotateUz(PositronDirection);
248
250 // create G4DynamicParticle object for the particle1
251 G4DynamicParticle* aParticle1= new G4DynamicParticle(
252 G4MuonPlus::MuonPlus(),MuPlusDirection,EmuPlus-Mmuon);
253 aParticleChange.AddSecondary(aParticle1);
254 // create G4DynamicParticle object for the particle2
255 G4DynamicParticle* aParticle2= new G4DynamicParticle(
256 G4MuonMinus::MuonMinus(),MuMinusDirection,EmuMinus-Mmuon);
257 aParticleChange.AddSecondary(aParticle2);
258
259 // Kill the incident positron
260 //
263
264 return &aParticleChange;
265}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
const G4ThreeVector & GetMomentumDirection() const
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
const G4double pi

◆ PrintInfoDefinition()

void G4AnnihiToMuPair::PrintInfoDefinition ( )

Definition at line 269 of file G4AnnihiToMuPair.cc.

270{
271 G4String comments ="e+e->mu+mu- annihilation, atomic e- at rest, SubType=.";
272 G4cout << G4endl << GetProcessName() << ": " << comments
274 G4cout << " threshold at " << LowestEnergyLimit/GeV << " GeV"
275 << " good description up to "
276 << HighestEnergyLimit/TeV << " TeV for all Z." << G4endl;
277}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

Referenced by BuildPhysicsTable().

◆ SetCrossSecFactor()

void G4AnnihiToMuPair::SetCrossSecFactor ( G4double  fac)

Definition at line 101 of file G4AnnihiToMuPair.cc.

103{
104 CrossSecFactor = fac;
105 G4cout << "The cross section for AnnihiToMuPair is artificially "
106 << "increased by the CrossSecFactor=" << CrossSecFactor << G4endl;
107}

Referenced by G4EmExtraPhysics::ConstructProcess().


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