Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4DecayWithSpin Class Reference

#include <G4DecayWithSpin.hh>

Inheritance diagram for G4DecayWithSpin:
Collaboration diagram for G4DecayWithSpin:

Public Member Functions

 G4DecayWithSpin (const G4String &processName="DecayWithSpin")
 
virtual ~G4DecayWithSpin ()
 
- Public Member Functions inherited from G4Decay
 G4Decay (const G4String &processName="Decay")
 
virtual ~G4Decay ()
 
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 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
 
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 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
 
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)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
- 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)
 
- 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 43 of file G4DecayWithSpin.hh.

Constructor & Destructor Documentation

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 }
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4Decay(const G4String &processName="Decay")
Definition: G4Decay.cc:63

Here is the call graph for this function:

G4DecayWithSpin::~G4DecayWithSpin ( )
virtual

Definition at line 59 of file G4DecayWithSpin.cc.

59 {}

Member Function Documentation

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

Reimplemented from G4Decay.

Definition at line 117 of file G4DecayWithSpin.cc.

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

Here is the call graph for this function:

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

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  {
78  // Generate random polarization direction
79 
80  G4double cost = 1. - 2.*G4UniformRand();
81  G4double sint = std::sqrt((1.-cost)*(1.+cost));
82 
84  G4double sinp = std::sin(phi);
85  G4double cosp = std::cos(phi);
86 
87  G4double px = sint*cosp;
88  G4double py = sint*sinp;
89  G4double pz = cost;
90 
91  parent_polarization.setX(px);
92  parent_polarization.setY(py);
93  parent_polarization.setZ(pz);
94 
95  }
96 
97 // decay table
98  G4DecayTable *decaytable = aParticleDef->GetDecayTable();
99  if (decaytable) {
100  for (G4int ip=0; ip<decaytable->entries(); ip++){
101  decaytable->GetDecayChannel(ip)->SetPolarization(parent_polarization);
102  }
103  }
104 
105  G4ParticleChangeForDecay* pParticleChangeForDecay;
106  pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
107 
108  pParticleChangeForDecay->ProposePolarization(parent_polarization);
109 
110  //G4cout << parent_polarization.x() << ", "
111  // << parent_polarization.y() << ", "
112  // << parent_polarization.z() << G4endl;
113 
114  return pParticleChangeForDecay;
115 }
CLHEP::Hep3Vector G4ThreeVector
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus GetTrackStatus() const
G4VDecayChannel * GetDecayChannel(G4int index) const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void setY(double)
void setZ(double)
void setX(double)
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181
static constexpr double twopi
Definition: G4SIunits.hh:76
G4int entries() const
G4DecayTable * GetDecayTable() const
virtual void Initialize(const G4Track &)
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetPolarization() const
double G4double
Definition: G4Types.hh:76
void SetPolarization(const G4ThreeVector &)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
void ProposePolarization(G4double Px, G4double Py, G4double Pz)

Here is the call graph for this function:


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