Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
64 // get particle
65  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
66  const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
67 
68 // get parent_polarization
69  G4ThreeVector parent_polarization = aParticle->GetPolarization();
70 
71  if(parent_polarization == G4ThreeVector(0,0,0))
72  {
73  // Generate random polarization direction
74 
75  G4double cost = 1. - 2.*G4UniformRand();
76  G4double sint = std::sqrt((1.-cost)*(1.+cost));
77 
79  G4double sinp = std::sin(phi);
80  G4double cosp = std::cos(phi);
81 
82  G4double px = sint*cosp;
83  G4double py = sint*sinp;
84  G4double pz = cost;
85 
86  parent_polarization.setX(px);
87  parent_polarization.setY(py);
88  parent_polarization.setZ(pz);
89 
90  }else{
91 
92  G4FieldManager* fieldMgr = aStep.GetTrack()->GetVolume()->
93  GetLogicalVolume()->GetFieldManager();
94 
95  if (!fieldMgr) {
96  G4TransportationManager *transportMgr =
98  G4PropagatorInField* fFieldPropagator =
99  transportMgr->GetPropagatorInField();
100  if (fFieldPropagator) fieldMgr =
101  fFieldPropagator->GetCurrentFieldManager();
102  }
103 
104  const G4Field* field = NULL;
105  if(fieldMgr)field = fieldMgr->GetDetectorField();
106 
107  if (field) {
108 
109  G4double point[4];
110  point[0] = (aStep.GetPostStepPoint()->GetPosition())[0];
111  point[1] = (aStep.GetPostStepPoint()->GetPosition())[1];
112  point[2] = (aStep.GetPostStepPoint()->GetPosition())[2];
113  point[3] = aTrack.GetGlobalTime();
114 
115  G4double fieldValue[6];
116  field -> GetFieldValue(point,fieldValue);
117 
118  G4ThreeVector B(fieldValue[0],fieldValue[1],fieldValue[2]);
119 
120  // Call the spin precession only for non-zero mag. field
121  if (B.mag2() > 0.) parent_polarization =
122  Spin_Precession(aStep,B,fRemainderLifeTime);
123 
124  }
125  }
126 
127 // decay table
128  G4DecayTable *decaytable = aParticleDef->GetDecayTable();
129 
130  if (decaytable) {
131  G4MuonDecayChannelWithSpin *decaychannel;
132  decaychannel = (G4MuonDecayChannelWithSpin*)decaytable->SelectADecayChannel();
133  if (decaychannel) decaychannel->SetPolarization(parent_polarization);
134  }
135 
136  G4ParticleChangeForDecay* pParticleChangeForDecay;
137 
138  pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
139 
140  pParticleChangeForDecay->ProposePolarization(parent_polarization);
141 
142  return pParticleChangeForDecay;
143 
144 }
145 
146 G4ThreeVector G4DecayWithSpin::Spin_Precession( const G4Step& aStep,
147  G4ThreeVector B, G4double deltatime )
148 {
149  G4double Bnorm = std::sqrt(sqr(B[0]) + sqr(B[1]) +sqr(B[2]) );
150 
151  G4double q = aStep.GetTrack()->GetDefinition()->GetPDGCharge();
152  G4double a = 1.165922e-3;
153  G4double s_omega = 8.5062e+7*rad/(s*kilogauss);
154 
155  G4double omega = -(q*s_omega)*(1.+a) * Bnorm;
156 
157  G4double rotationangle = deltatime * omega;
158 
159  G4Transform3D SpinRotation = G4Rotate3D(rotationangle,B.unit());
160 
161  G4Vector3D Spin = aStep.GetTrack() -> GetPolarization();
162 
163  G4Vector3D newSpin = SpinRotation * Spin;
164 
165 #ifdef G4VERBOSE
166  if (GetVerboseLevel()>2) {
167  G4double normspin = std::sqrt(Spin*Spin);
168  G4double normnewspin = std::sqrt(newSpin*newSpin);
169  //G4double cosalpha = Spin*newSpin/normspin/normnewspin;
170  //G4double alpha = std::acos(cosalpha);
171 
172  G4cout << "AT REST::: PARAMETERS " << G4endl;
173  G4cout << "Initial spin : " << Spin << G4endl;
174  G4cout << "Delta time : " << deltatime << G4endl;
175  G4cout << "Rotation angle: " << rotationangle/rad << G4endl;
176  G4cout << "New spin : " << newSpin << G4endl;
177  G4cout << "Checked norms : " << normspin <<" " << normnewspin << G4endl;
178  }
179 #endif
180 
181  return newSpin;
182 
183 }
G4ParticleDefinition * GetDefinition() const
CLHEP::Hep3Vector G4ThreeVector
const G4DynamicParticle * GetDynamicParticle() const
const XML_Char * s
G4ParticleDefinition * GetDefinition() const
void setY(double)
void setZ(double)
void setX(double)
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Decay.cc:181
G4double fRemainderLifeTime
Definition: G4Decay.hh:175
G4DecayTable * GetDecayTable() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetPosition() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
G4double GetGlobalTime() const
G4VDecayChannel * SelectADecayChannel()
Definition: G4DecayTable.cc:81
static G4TransportationManager * GetTransportationManager()
virtual ~G4DecayWithSpin()
HepGeom::Rotate3D G4Rotate3D
virtual G4VParticleChange * DecayIt(const G4Track &aTrack, const G4Step &aStep)
G4FieldManager * GetCurrentFieldManager()
Hep3Vector unit() const
G4StepPoint * GetPostStepPoint() const
const G4ThreeVector & GetPolarization() const
double mag2() 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
G4DecayWithSpin(const G4String &processName="DecayWithSpin")
G4int GetVerboseLevel() const
Definition: G4Decay.hh:188
void ProposePolarization(G4double Px, G4double Py, G4double Pz)