G4DecayWithSpin Class Reference

#include <G4DecayWithSpin.hh>

Inheritance diagram for G4DecayWithSpin:

G4Decay G4VRestDiscreteProcess G4VProcess

Public Member Functions

 G4DecayWithSpin (const G4String &processName="DecayWithSpin")
virtual ~G4DecayWithSpin ()

Protected Member Functions

virtual G4VParticleChangeDecayIt (const G4Track &aTrack, const G4Step &aStep)

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.

References DECAY_WithSpin, and G4VProcess::SetProcessSubType().

00052                                                            :G4Decay(processName)
00053 {
00054   // set Process Sub Type   
00055   SetProcessSubType(static_cast<int>(DECAY_WithSpin));
00056 
00057 }

G4DecayWithSpin::~G4DecayWithSpin (  )  [virtual]

Definition at line 59 of file G4DecayWithSpin.cc.

00059 {}


Member Function Documentation

G4VParticleChange * G4DecayWithSpin::DecayIt ( const G4Track aTrack,
const G4Step aStep 
) [protected, virtual]

Reimplemented from G4Decay.

Definition at line 61 of file G4DecayWithSpin.cc.

References G4Decay::DecayIt(), G4Decay::fRemainderLifeTime, G4UniformRand, G4PropagatorInField::GetCurrentFieldManager(), G4ParticleDefinition::GetDecayTable(), G4DynamicParticle::GetDefinition(), G4FieldManager::GetDetectorField(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetPolarization(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4TransportationManager::GetPropagatorInField(), G4Step::GetTrack(), G4TransportationManager::GetTransportationManager(), G4Track::GetVolume(), G4ParticleChangeForDecay::ProposePolarization(), G4DecayTable::SelectADecayChannel(), and G4MuonDecayChannelWithSpin::SetPolarization().

00062 {
00063 
00064 // get particle
00065   const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00066   const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
00067 
00068 // get parent_polarization
00069   G4ThreeVector parent_polarization = aParticle->GetPolarization();
00070 
00071   if(parent_polarization == G4ThreeVector(0,0,0))
00072   {
00073     // Generate random polarization direction
00074 
00075     G4double cost = 1. - 2.*G4UniformRand();
00076     G4double sint = std::sqrt((1.-cost)*(1.+cost));
00077 
00078     G4double phi = twopi*G4UniformRand();
00079     G4double sinp = std::sin(phi);
00080     G4double cosp = std::cos(phi);
00081 
00082     G4double px = sint*cosp;
00083     G4double py = sint*sinp;
00084     G4double pz = cost;
00085 
00086     parent_polarization.setX(px);
00087     parent_polarization.setY(py);
00088     parent_polarization.setZ(pz);
00089 
00090   }else{
00091 
00092     G4FieldManager* fieldMgr = aStep.GetTrack()->GetVolume()->
00093                                      GetLogicalVolume()->GetFieldManager();
00094 
00095     if (!fieldMgr) {
00096        G4TransportationManager *transportMgr =
00097                          G4TransportationManager::GetTransportationManager();
00098        G4PropagatorInField* fFieldPropagator = 
00099                                         transportMgr->GetPropagatorInField();
00100        if (fFieldPropagator) fieldMgr = 
00101                                   fFieldPropagator->GetCurrentFieldManager();
00102     }
00103 
00104     const G4Field* field = NULL;
00105     if(fieldMgr)field = fieldMgr->GetDetectorField();
00106 
00107     if (field) {
00108 
00109        G4double point[4];
00110        point[0] = (aStep.GetPostStepPoint()->GetPosition())[0];
00111        point[1] = (aStep.GetPostStepPoint()->GetPosition())[1];
00112        point[2] = (aStep.GetPostStepPoint()->GetPosition())[2];
00113        point[3] = aTrack.GetGlobalTime();
00114 
00115        G4double fieldValue[6];
00116        field -> GetFieldValue(point,fieldValue);
00117 
00118        G4ThreeVector B(fieldValue[0],fieldValue[1],fieldValue[2]);
00119 
00120        // Call the spin precession only for non-zero mag. field
00121        if (B.mag2() > 0.) parent_polarization = 
00122                                  Spin_Precession(aStep,B,fRemainderLifeTime);
00123 
00124     }
00125   }
00126 
00127 // decay table
00128   G4DecayTable *decaytable = aParticleDef->GetDecayTable();
00129 
00130   if (decaytable) {
00131      G4MuonDecayChannelWithSpin *decaychannel;
00132      decaychannel = (G4MuonDecayChannelWithSpin*)decaytable->SelectADecayChannel();
00133      if (decaychannel) decaychannel->SetPolarization(parent_polarization);
00134   }
00135 
00136   G4ParticleChangeForDecay* pParticleChangeForDecay;
00137 
00138   pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
00139 
00140   pParticleChangeForDecay->ProposePolarization(parent_polarization);
00141 
00142   return pParticleChangeForDecay;
00143 
00144 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:45 2013 for Geant4 by  doxygen 1.4.7