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

#include <G4DecayWithSpin.hh>

+ Inheritance diagram for G4DecayWithSpin:

Public Member Functions

 G4DecayWithSpin (const G4String &processName="DecayWithSpin")
 
virtual ~G4DecayWithSpin ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
- Public Member Functions inherited from G4Decay
 G4Decay (const G4String &processName="Decay")
 
virtual ~G4Decay ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual void StartTracking (G4Track *) override
 
virtual void EndTracking () override
 
void SetExtDecayer (G4VExtDecayer *)
 
const G4VExtDecayerGetExtDecayer () const
 
G4double GetRemainderLifeTime () const
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
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 ()
 
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 &)
 

Protected Member Functions

virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
- 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) override
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition) override
 
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 prevStepSize)
 
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 = 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 43 of file G4DecayWithSpin.hh.

Constructor & Destructor Documentation

◆ G4DecayWithSpin()

G4DecayWithSpin::G4DecayWithSpin ( const G4String processName = "DecayWithSpin")

Definition at line 52 of file G4DecayWithSpin.cc.

52 :G4Decay(processName)
53{
54 // set Process Sub Type
55 SetProcessSubType(static_cast<int>(DECAY_WithSpin));
56
57}
@ DECAY_WithSpin
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4DecayWithSpin()

G4DecayWithSpin::~G4DecayWithSpin ( )
virtual

Definition at line 59 of file G4DecayWithSpin.cc.

59{}

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4DecayWithSpin::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overrideprotectedvirtual

Reimplemented from G4Decay.

Definition at line 109 of file G4DecayWithSpin.cc.

110{
111
112// get particle
113 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
114 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
115
116// get parent_polarization
117 G4ThreeVector parent_polarization = aParticle->GetPolarization();
118
119 if(parent_polarization == G4ThreeVector(0,0,0)) {
120 // Generate random polarization direction
121 G4double cost = 1. - 2.*G4UniformRand();
122 G4double sint = std::sqrt((1.-cost)*(1.+cost));
123
124 G4double phi = twopi*G4UniformRand();
125 G4double sinp = std::sin(phi);
126 G4double cosp = std::cos(phi);
127
128 G4double px = sint*cosp;
129 G4double py = sint*sinp;
130 G4double pz = cost;
131
132 parent_polarization.setX(px);
133 parent_polarization.setY(py);
134 parent_polarization.setZ(pz);
135
136 }else{
137
138 G4FieldManager* fieldMgr = aStep.GetTrack()->GetVolume()->
139 GetLogicalVolume()->GetFieldManager();
140 if (fieldMgr == nullptr) {
141 G4TransportationManager *transportMgr =
143 G4PropagatorInField* fFieldPropagator =
144 transportMgr->GetPropagatorInField();
145 if (fFieldPropagator) fieldMgr =
146 fFieldPropagator->GetCurrentFieldManager();
147 }
148
149 const G4Field* field = nullptr;
150 if (fieldMgr != nullptr) field = fieldMgr->GetDetectorField();
151
152 if ( field != nullptr ) {
153 G4double point[4];
154 point[0] = (aStep.GetPostStepPoint()->GetPosition())[0];
155 point[1] = (aStep.GetPostStepPoint()->GetPosition())[1];
156 point[2] = (aStep.GetPostStepPoint()->GetPosition())[2];
157 point[3] = aTrack.GetGlobalTime();
158
159 G4double fieldValue[6] ={ 0., 0., 0., 0., 0., 0.};
160 field -> GetFieldValue(point,fieldValue);
161 G4ThreeVector B(fieldValue[0],fieldValue[1],fieldValue[2]);
162
163 // Call the spin precession only for non-zero mag. field
164 if (B.mag2() > 0.) parent_polarization =
165 Spin_Precession(aStep,B,fRemainderLifeTime);
166
167 }
168 }
169
170// decay table
171 G4DecayTable *decaytable = aParticleDef->GetDecayTable();
172 if ( decaytable != nullptr) {
173 for (G4int ip=0; ip<decaytable->entries(); ip++){
174 decaytable->GetDecayChannel(ip)->SetPolarization(parent_polarization);
175 }
176 }
177
178 G4ParticleChangeForDecay* pParticleChangeForDecay;
179 pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
180 pParticleChangeForDecay->ProposePolarization(parent_polarization);
181
182 return pParticleChangeForDecay;
183
184}
G4double B(G4double temperature)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
void setY(double)
void setZ(double)
void setX(double)
G4VDecayChannel * GetDecayChannel(G4int index) const
G4int entries() const
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:180
G4double fRemainderLifeTime
Definition: G4Decay.hh:173
G4ParticleDefinition * GetDefinition() const
const G4ThreeVector & GetPolarization() const
const G4Field * GetDetectorField() const
void ProposePolarization(G4double Px, G4double Py, G4double Pz)
G4DecayTable * GetDecayTable() const
G4FieldManager * GetCurrentFieldManager()
const G4ThreeVector & GetPosition() const
G4Track * GetTrack() const
G4StepPoint * GetPostStepPoint() const
G4VPhysicalVolume * GetVolume() const
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
void SetPolarization(const G4ThreeVector &)

◆ PostStepDoIt()

G4VParticleChange * G4DecayWithSpin::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overrideprotectedvirtual

Reimplemented from G4Decay.

Definition at line 61 of file G4DecayWithSpin.cc.

62{
63 if ( (aTrack.GetTrackStatus() == fStopButAlive ) ||
64 (aTrack.GetTrackStatus() == fStopAndKill ) ){
67 }
68
69// get particle
70 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
71 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
72
73// get parent_polarization
74 G4ThreeVector parent_polarization = aParticle->GetPolarization();
75
76 if(parent_polarization == G4ThreeVector(0,0,0)){
77 // Generate random polarization direction
78 G4double cost = 1. - 2.*G4UniformRand();
79 G4double sint = std::sqrt((1.-cost)*(1.+cost));
80
81 G4double phi = twopi*G4UniformRand();
82 G4double sinp = std::sin(phi);
83 G4double cosp = std::cos(phi);
84
85 G4double px = sint*cosp;
86 G4double py = sint*sinp;
87 G4double pz = cost;
88
89 parent_polarization.setX(px);
90 parent_polarization.setY(py);
91 parent_polarization.setZ(pz);
92 }
93
94// decay table
95 G4DecayTable *decaytable = aParticleDef->GetDecayTable();
96 if (decaytable != nullptr) {
97 for (G4int ip=0; ip<decaytable->entries(); ip++){
98 decaytable->GetDecayChannel(ip)->SetPolarization(parent_polarization);
99 }
100 }
101
102 G4ParticleChangeForDecay* pParticleChangeForDecay;
103 pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
104 pParticleChangeForDecay->ProposePolarization(parent_polarization);
105
106 return pParticleChangeForDecay;
107}
@ fStopAndKill
@ fStopButAlive
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:176
void Initialize(const G4Track &) final
G4TrackStatus GetTrackStatus() const

◆ ProcessDescription()

void G4DecayWithSpin::ProcessDescription ( std::ostream &  outFile) const
overridevirtual

Reimplemented from G4Decay.

Definition at line 225 of file G4DecayWithSpin.cc.

226{
227 outFile << GetProcessName()
228 << ": Decay of particles considering parent polarization \n"
229 << "kinematics of daughters are dertermined by DecayChannels \n";
230}
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

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