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

#include <G4PionDecayMakeSpin.hh>

+ Inheritance diagram for G4PionDecayMakeSpin:

Public Member Functions

 G4PionDecayMakeSpin (const G4String &processName="Decay")
 
virtual ~G4PionDecayMakeSpin ()
 
- Public Member Functions inherited from G4Decay
 G4Decay (const G4String &processName="Decay")
 
virtual ~G4Decay ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
void SetExtDecayer (G4VExtDecayer *)
 
const G4VExtDecayerGetExtDecayer () const
 
G4double GetRemainderLifeTime () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
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
 

Protected Member Functions

virtual void DaughterPolarization (const G4Track &aTrack, G4DecayProducts *products)
 
- Protected Member Functions inherited from G4Decay
virtual G4VParticleChangeDecayIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void DaughterPolarization (const G4Track &aTrack, G4DecayProducts *products)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4Decay
G4int verboseLevel
 
const G4double HighestValue
 
G4double fRemainderLifeTime
 
G4ParticleChangeForDecay fParticleChangeForDecay
 
G4VExtDecayerpExtDecayer
 
- 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 39 of file G4PionDecayMakeSpin.hh.

Constructor & Destructor Documentation

◆ G4PionDecayMakeSpin()

G4PionDecayMakeSpin::G4PionDecayMakeSpin ( const G4String processName = "Decay")

Definition at line 38 of file G4PionDecayMakeSpin.cc.

39 : G4Decay(processName)
40{
41 // set Process Sub Type
42 SetProcessSubType(static_cast<int>(DECAY_PionMakeSpin));
43
44}
@ DECAY_PionMakeSpin
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4PionDecayMakeSpin()

G4PionDecayMakeSpin::~G4PionDecayMakeSpin ( )
virtual

Definition at line 46 of file G4PionDecayMakeSpin.cc.

46{ }

Member Function Documentation

◆ DaughterPolarization()

void G4PionDecayMakeSpin::DaughterPolarization ( const G4Track aTrack,
G4DecayProducts products 
)
protectedvirtual

Reimplemented from G4Decay.

Definition at line 48 of file G4PionDecayMakeSpin.cc.

50{
51 // This routine deals only with particles that can decay into a muon
52 // pi+, pi-, K+, K- and K0_long
53
54 // get particle
55
56 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
57 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
58
59 G4ParticleDefinition* aMuonPlus =
61 G4ParticleDefinition* aMuonMinus =
63 G4ParticleDefinition* aPionPlus =
65 G4ParticleDefinition* aPionMinus =
67 G4ParticleDefinition* aKaonPlus =
69 G4ParticleDefinition* aKaonMinus =
71 G4ParticleDefinition* aKaon0Long =
73 G4ParticleDefinition* aNeutrinoMu =
75 G4ParticleDefinition* aAntiNeutrinoMu =
77
78 if( aParticleDef == aPionPlus ||
79 aParticleDef == aPionMinus ||
80 aParticleDef == aKaonPlus ||
81 aParticleDef == aKaonMinus ||
82 aParticleDef == aKaon0Long ) {
83 } else {
84 return;
85 }
86
87 G4DynamicParticle* aMuon = NULL;
88
89 G4double emu(0), eneutrino(0);
90 G4ThreeVector p_muon, p_neutrino;
91
92 G4int numberOfSecondaries = products->entries();
93
94 if (numberOfSecondaries > 0) {
95 for (G4int index=0; index < numberOfSecondaries; index++)
96 {
97 G4DynamicParticle* aSecondary = (*products)[index];
98 const G4ParticleDefinition* aSecondaryDef = aSecondary->GetDefinition();
99
100 if (aSecondaryDef == aMuonPlus ||
101 aSecondaryDef == aMuonMinus ) {
102 // Muon+ or Muon-
103 aMuon = aSecondary;
104 emu = aSecondary->GetTotalEnergy();
105 p_muon = aSecondary->GetMomentum();
106 } else if (aSecondaryDef == aNeutrinoMu ||
107 aSecondaryDef == aAntiNeutrinoMu ) {
108 // Muon-Neutrino / Muon-Anti-Neutrino
109 eneutrino = aSecondary->GetTotalEnergy();
110 p_neutrino = aSecondary->GetMomentum();
111 }
112 }
113 }
114
115 // This routine deals only with decays with a
116 // muon and mu-(anti)neutrinos in the final state
117
118 if (!aMuon) return;
119 if (eneutrino==0||emu==0) return;
120
121 G4ThreeVector spin(0,0,0);
122
123 const G4DynamicParticle* theParentParticle = products->GetParentParticle();
124
125 G4double amass = theParentParticle->GetMass();
126 G4double emmu = aMuonPlus->GetPDGMass();
127
128 if (numberOfSecondaries == 2 ) {
129
130 G4double scale = - (eneutrino - ( p_muon * p_neutrino )/(emu+emmu));
131
132 p_muon = scale * p_muon;
133 p_neutrino = emmu * p_neutrino;
134 spin = p_muon + p_neutrino;
135
136 scale = 2./(amass*amass-emmu*emmu);
137 spin = scale * spin;
138
139 if (aParticle->GetCharge() < 0.0) spin = -spin;
140
141 } else {
142
143 spin = G4RandomDirection();
144
145 }
146
147 spin = spin.unit();
148
149 aMuon->SetPolarization(spin.x(),spin.y(),spin.z());
150
151 return;
152}
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4int entries() const
const G4DynamicParticle * GetParentParticle() const
G4double GetMass() const
G4double GetCharge() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
const G4DynamicParticle * GetDynamicParticle() const

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