Geant4  10.02.p02
G4DecayWithSpin.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // ------------------------------------------------------------
27 // GEANT 4 class header file
28 //
29 // History:
30 // 17 August 2004 P. Gumplinger, T. MacPhail
31 // 11 April 2008 Kamil Sedlak (PSI), Toni Shiroka (PSI)
32 // ------------------------------------------------------------
33 //
34 #include "G4DecayWithSpin.hh"
35 
36 #include "G4Step.hh"
37 #include "G4Track.hh"
38 #include "G4DecayTable.hh"
40 
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4Vector3D.hh"
44 
46 #include "G4PropagatorInField.hh"
47 #include "G4FieldManager.hh"
48 #include "G4Field.hh"
49 
50 #include "G4Transform3D.hh"
51 
52 G4DecayWithSpin::G4DecayWithSpin(const G4String& processName):G4Decay(processName)
53 {
54  // set Process Sub Type
55  SetProcessSubType(static_cast<int>(DECAY_WithSpin));
56 
57 }
58 
60 
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 }
116 
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 =
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 }
199 
201  G4ThreeVector B, G4double deltatime )
202 {
203  G4double Bnorm = std::sqrt(sqr(B[0]) + sqr(B[1]) +sqr(B[2]) );
204 
205  G4double q = aStep.GetTrack()->GetDefinition()->GetPDGCharge();
206  G4double a = 1.165922e-3;
207  G4double s_omega = 8.5062e+7*rad/(s*kilogauss);
208 
209  G4double omega = -(q*s_omega)*(1.+a) * Bnorm;
210 
211  G4double rotationangle = deltatime * omega;
212 
213  G4Transform3D SpinRotation = G4Rotate3D(rotationangle,B.unit());
214 
215  G4Vector3D Spin = aStep.GetTrack() -> GetPolarization();
216 
217  G4Vector3D newSpin = SpinRotation * Spin;
218 
219 #ifdef G4VERBOSE
220  if (GetVerboseLevel()>2) {
221  G4double normspin = std::sqrt(Spin*Spin);
222  G4double normnewspin = std::sqrt(newSpin*newSpin);
223  //G4double cosalpha = Spin*newSpin/normspin/normnewspin;
224  //G4double alpha = std::acos(cosalpha);
225 
226  G4cout << "AT REST::: PARAMETERS " << G4endl;
227  G4cout << "Initial spin : " << Spin << G4endl;
228  G4cout << "Delta time : " << deltatime << G4endl;
229  G4cout << "Rotation angle: " << rotationangle/rad << G4endl;
230  G4cout << "New spin : " << newSpin << G4endl;
231  G4cout << "Checked norms : " << normspin <<" " << normnewspin << G4endl;
232  }
233 #endif
234 
235  return newSpin;
236 
237 }
G4ParticleDefinition * GetDefinition() const
virtual G4VParticleChange * AtRestDoIt(const G4Track &aTrack, const G4Step &aStep)
CLHEP::Hep3Vector G4ThreeVector
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus GetTrackStatus() const
G4VDecayChannel * GetDecayChannel(G4int index) const
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
G4ParticleDefinition * GetDefinition() const
G4double a
Definition: TRTMaterials.hh:39
double B(double temperature)
int G4int
Definition: G4Types.hh:78
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
static const double s
Definition: G4SIunits.hh:168
G4int entries() const
G4DecayTable * GetDecayTable() const
virtual void Initialize(const G4Track &)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4ThreeVector Spin_Precession(const G4Step &aStep, G4ThreeVector B, G4double deltatime)
const G4ThreeVector & GetPosition() const
static const double twopi
Definition: G4SIunits.hh:75
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
HepGeom::Transform3D G4Transform3D
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Step.hh:76
G4double GetGlobalTime() const
static const double rad
Definition: G4SIunits.hh:148
static G4TransportationManager * GetTransportationManager()
virtual ~G4DecayWithSpin()
HepGeom::Rotate3D G4Rotate3D
G4FieldManager * GetCurrentFieldManager()
G4StepPoint * GetPostStepPoint() const
const G4ThreeVector & GetPolarization() const
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
const G4Field * GetDetectorField() const
G4Track * GetTrack() const
G4double GetPDGCharge() const
G4PropagatorInField * GetPropagatorInField() const
static const double kilogauss
Definition: G4SIunits.hh:268
G4DecayWithSpin(const G4String &processName="DecayWithSpin")
void SetPolarization(const G4ThreeVector &)
G4ParticleChangeForDecay fParticleChangeForDecay
Definition: G4Decay.hh:178
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188
void ProposePolarization(G4double Px, G4double Py, G4double Pz)