Geant4 9.6.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 &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
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 *)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 57 of file G4AnnihiToMuPair.hh.

Constructor & Destructor Documentation

◆ G4AnnihiToMuPair()

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

Definition at line 58 of file G4AnnihiToMuPair.cc.

59 :G4VDiscreteProcess (processName, type)
60{
61 //e+ Energy threshold
62 const G4double Mu_massc2 = G4MuonPlus::MuonPlus()->GetPDGMass();
63 LowestEnergyLimit = 2*Mu_massc2*Mu_massc2/electron_mass_c2 - electron_mass_c2;
64
65 //modele ok up to 1000 TeV due to neglected Z-interference
66 HighestEnergyLimit = 1000*TeV;
67
68 CurrentSigma = 0.0;
69 CrossSecFactor = 1.;
71
72}
double G4double
Definition: G4Types.hh:64
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4AnnihiToMuPair()

G4AnnihiToMuPair::~G4AnnihiToMuPair ( )

Definition at line 76 of file G4AnnihiToMuPair.cc.

77{ }

Member Function Documentation

◆ BuildPhysicsTable()

void G4AnnihiToMuPair::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 88 of file G4AnnihiToMuPair.cc.

91{
92 CurrentSigma = 0.0;
94}

◆ ComputeCrossSectionPerAtom()

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

Definition at line 108 of file G4AnnihiToMuPair.cc.

111{
112 static const G4double Mmuon = G4MuonPlus::MuonPlus()->GetPDGMass();
113 static const G4double Rmuon = elm_coupling/Mmuon; //classical particle radius
114 static const G4double Sig0 = pi*Rmuon*Rmuon/3.; //constant in crossSection
115
116 G4double CrossSection = 0.;
117 if (Epos < LowestEnergyLimit) return CrossSection;
118
119 G4double xi = LowestEnergyLimit/Epos;
120 G4double SigmaEl = Sig0*xi*(1.+xi/2.)*sqrt(1.-xi); // per electron
121 CrossSection = SigmaEl*Z; // number of electrons per atom
122 return CrossSection;
123}
const G4double pi

Referenced by CrossSectionPerVolume().

◆ CrossSectionPerVolume()

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

Definition at line 127 of file G4AnnihiToMuPair.cc.

129{
130 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
131 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
132
133 G4double SIGMA = 0.0;
134
135 for ( size_t i=0 ; i < aMaterial->GetNumberOfElements() ; ++i )
136 {
137 G4double AtomicZ = (*theElementVector)[i]->GetZ();
138 SIGMA += NbOfAtomsPerVolume[i] *
139 ComputeCrossSectionPerAtom(PositronEnergy,AtomicZ);
140 }
141 return SIGMA;
142}
std::vector< G4Element * > G4ElementVector
G4double ComputeCrossSectionPerAtom(G4double PositronEnergy, G4double AtomicZ)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205

Referenced by GetMeanFreePath(), and PostStepDoIt().

◆ GetCrossSecFactor()

G4double G4AnnihiToMuPair::GetCrossSecFactor ( )
inline

Definition at line 81 of file G4AnnihiToMuPair.hh.

81{return CrossSecFactor;};

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 146 of file G4AnnihiToMuPair.cc.

151{
152 const G4DynamicParticle* aDynamicPositron = aTrack.GetDynamicParticle();
153 G4double PositronEnergy = aDynamicPositron->GetKineticEnergy()
154 +electron_mass_c2;
155 G4Material* aMaterial = aTrack.GetMaterial();
156 CurrentSigma = CrossSectionPerVolume(PositronEnergy, aMaterial);
157
158 // increase the CrossSection by CrossSecFactor (default 1)
159 G4double mfp = DBL_MAX;
160 if(CurrentSigma > DBL_MIN) mfp = 1.0/(CurrentSigma*CrossSecFactor);
161
162 return mfp;
163}
G4double CrossSectionPerVolume(G4double PositronEnergy, const G4Material *)
G4double GetKineticEnergy() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MIN
Definition: templates.hh:75
#define DBL_MAX
Definition: templates.hh:83

◆ IsApplicable()

G4bool G4AnnihiToMuPair::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 81 of file G4AnnihiToMuPair.cc.

82{
83 return ( &particle == G4Positron::Positron() );
84}
static G4Positron * Positron()
Definition: G4Positron.cc:94

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 167 of file G4AnnihiToMuPair.cc.

172{
173
175 static const G4double Mele=electron_mass_c2;
176 static const G4double Mmuon=G4MuonPlus::MuonPlus()->GetPDGMass();
177
178 // current Positron energy and direction, return if energy too low
179 const G4DynamicParticle *aDynamicPositron = aTrack.GetDynamicParticle();
180 G4double Epos = aDynamicPositron->GetKineticEnergy() + Mele;
181
182 // test of cross section
183 if(CurrentSigma*G4UniformRand() >
184 CrossSectionPerVolume(Epos, aTrack.GetMaterial()))
185 {
186 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
187 }
188
189 if (Epos < LowestEnergyLimit) {
190 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
191 }
192
193 G4ParticleMomentum PositronDirection =
194 aDynamicPositron->GetMomentumDirection();
195 G4double xi = LowestEnergyLimit/Epos; // xi is always less than 1,
196 // goes to 0 at high Epos
197
198 // generate cost
199 //
200 G4double cost;
201 do cost = 2.*G4UniformRand()-1.;
202 while (2.*G4UniformRand() > 1.+xi+cost*cost*(1.-xi) );
203 //1+cost**2 at high Epos
204 G4double sint = sqrt(1.-cost*cost);
205
206 // generate phi
207 //
208 G4double phi=2.*pi*G4UniformRand();
209
210 G4double Ecm = sqrt(0.5*Mele*(Epos+Mele));
211 G4double Pcm = sqrt(Ecm*Ecm-Mmuon*Mmuon);
212 G4double beta = sqrt((Epos-Mele)/(Epos+Mele));
213 G4double gamma = Ecm/Mele; // =sqrt((Epos+Mele)/(2.*Mele));
214 G4double Pt = Pcm*sint;
215
216 // energy and momentum of the muons in the Lab
217 //
218 G4double EmuPlus = gamma*( Ecm+cost*beta*Pcm);
219 G4double EmuMinus = gamma*( Ecm-cost*beta*Pcm);
220 G4double PmuPlusZ = gamma*(beta*Ecm+cost* Pcm);
221 G4double PmuMinusZ = gamma*(beta*Ecm-cost* Pcm);
222 G4double PmuPlusX = Pt*cos(phi);
223 G4double PmuPlusY = Pt*sin(phi);
224 G4double PmuMinusX =-Pt*cos(phi);
225 G4double PmuMinusY =-Pt*sin(phi);
226 // absolute momenta
227 G4double PmuPlus = sqrt(Pt*Pt+PmuPlusZ *PmuPlusZ );
228 G4double PmuMinus = sqrt(Pt*Pt+PmuMinusZ*PmuMinusZ);
229
230 // mu+ mu- directions for Positron in z-direction
231 //
233 MuPlusDirection ( PmuPlusX/PmuPlus, PmuPlusY/PmuPlus, PmuPlusZ/PmuPlus );
235 MuMinusDirection(PmuMinusX/PmuMinus,PmuMinusY/PmuMinus,PmuMinusZ/PmuMinus);
236
237 // rotate to actual Positron direction
238 //
239 MuPlusDirection.rotateUz(PositronDirection);
240 MuMinusDirection.rotateUz(PositronDirection);
241
243 // create G4DynamicParticle object for the particle1
244 G4DynamicParticle* aParticle1= new G4DynamicParticle(
245 G4MuonPlus::MuonPlus(),MuPlusDirection,EmuPlus-Mmuon);
246 aParticleChange.AddSecondary(aParticle1);
247 // create G4DynamicParticle object for the particle2
248 G4DynamicParticle* aParticle2= new G4DynamicParticle(
249 G4MuonMinus::MuonMinus(),MuMinusDirection,EmuMinus-Mmuon);
250 aParticleChange.AddSecondary(aParticle2);
251
252 // Kill the incident positron
253 //
256
257 return &aParticleChange;
258}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
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:289

◆ PrintInfoDefinition()

void G4AnnihiToMuPair::PrintInfoDefinition ( )

Definition at line 262 of file G4AnnihiToMuPair.cc.

263{
264 G4String comments ="e+e->mu+mu- annihilation, atomic e- at rest, SubType=.";
265 G4cout << G4endl << GetProcessName() << ": " << comments
267 G4cout << " threshold at " << LowestEnergyLimit/GeV << " GeV"
268 << " good description up to "
269 << HighestEnergyLimit/TeV << " TeV for all Z." << G4endl;
270}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

Referenced by BuildPhysicsTable().

◆ SetCrossSecFactor()

void G4AnnihiToMuPair::SetCrossSecFactor ( G4double  fac)

Definition at line 98 of file G4AnnihiToMuPair.cc.

100{
101 CrossSecFactor = fac;
102 G4cout << "The cross section for AnnihiToMuPair is artificially "
103 << "increased by the CrossSecFactor=" << CrossSecFactor << G4endl;
104}

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