source: trunk/source/processes/decay/src/G4DecayWithSpin.cc @ 961

Last change on this file since 961 was 961, checked in by garnier, 15 years ago

update processes

File size: 6.3 KB
Line 
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"
39#include "G4MuonDecayChannelWithSpin.hh"
40
41#include "G4Vector3D.hh"
42
43#include "G4TransportationManager.hh"
44#include "G4PropagatorInField.hh"
45#include "G4FieldManager.hh"
46#include "G4Field.hh"
47
48#include "G4Transform3D.hh"
49
50G4DecayWithSpin::G4DecayWithSpin(const G4String& processName):G4Decay(processName)
51{
52  // set Process Sub Type   
53  SetProcessSubType(static_cast<int>(DECAY_WithSpin));
54
55}
56
57G4DecayWithSpin::~G4DecayWithSpin(){}
58
59G4VParticleChange* G4DecayWithSpin::DecayIt(const G4Track& aTrack, const G4Step& aStep)
60{
61
62// get particle
63  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
64  G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
65
66// get parent_polarization
67  G4ThreeVector parent_polarization = aParticle->GetPolarization();
68
69  if(parent_polarization == G4ThreeVector(0,0,0))
70  {
71    // Generate random polarization direction
72
73    G4double cost = 1. - 2.*G4UniformRand();
74    G4double sint = std::sqrt((1.-cost)*(1.+cost));
75
76    G4double phi = twopi*G4UniformRand();
77    G4double sinp = std::sin(phi);
78    G4double cosp = std::cos(phi);
79
80    G4double px = sint*cosp;
81    G4double py = sint*sinp;
82    G4double pz = cost;
83
84    parent_polarization.setX(px);
85    parent_polarization.setY(py);
86    parent_polarization.setZ(pz);
87
88  }else{
89
90    G4FieldManager* fieldMgr = aStep.GetTrack()->GetVolume()->
91                                     GetLogicalVolume()->GetFieldManager();
92
93    if (!fieldMgr) {
94       G4TransportationManager *transportMgr =
95                         G4TransportationManager::GetTransportationManager();
96       G4PropagatorInField* fFieldPropagator = 
97                                        transportMgr->GetPropagatorInField();
98       if (fFieldPropagator) fieldMgr = 
99                                  fFieldPropagator->GetCurrentFieldManager();
100    }
101
102    const G4Field* field = NULL;
103    if(fieldMgr)field = fieldMgr->GetDetectorField();
104
105    if (field) {
106
107       G4double point[4];
108       point[0] = (aStep.GetPostStepPoint()->GetPosition())[0];
109       point[1] = (aStep.GetPostStepPoint()->GetPosition())[1];
110       point[2] = (aStep.GetPostStepPoint()->GetPosition())[2];
111       point[3] = aTrack.GetGlobalTime();
112
113       G4double fieldValue[6];
114       field -> GetFieldValue(point,fieldValue);
115
116       G4ThreeVector B(fieldValue[0],fieldValue[1],fieldValue[2]);
117
118       // Call the spin precession only for non-zero mag. field
119       if (B.mag2() > 0.) parent_polarization = 
120                                 Spin_Precession(aStep,B,fRemainderLifeTime);
121
122    }
123  }
124
125// decay table
126  G4DecayTable *decaytable = aParticleDef->GetDecayTable();
127
128  if (decaytable) {
129     G4MuonDecayChannelWithSpin *decaychannel;
130     decaychannel = (G4MuonDecayChannelWithSpin*)decaytable->SelectADecayChannel();
131     if (decaychannel) decaychannel->SetPolarization(parent_polarization);
132  }
133
134  G4ParticleChangeForDecay* pParticleChangeForDecay;
135
136  pParticleChangeForDecay = (G4ParticleChangeForDecay*)G4Decay::DecayIt(aTrack,aStep);
137
138  pParticleChangeForDecay->ProposePolarization(parent_polarization);
139
140  return pParticleChangeForDecay;
141
142}
143
144G4ThreeVector G4DecayWithSpin::Spin_Precession( const G4Step& aStep,
145                                        G4ThreeVector B, G4double deltatime )
146{
147  G4double Bnorm = std::sqrt(sqr(B[0])  + sqr(B[1]) +sqr(B[2]) );
148
149  G4double q = aStep.GetTrack()->GetDefinition()->GetPDGCharge();
150  G4double a = 1.165922e-3;
151  G4double s_omega = 8.5062e+7*rad/(s*kilogauss);
152
153  G4double omega = -(q*s_omega)*(1.+a) * Bnorm;
154
155  G4double rotationangle = deltatime * omega;
156
157  G4Transform3D SpinRotation = G4Rotate3D(rotationangle,B.unit());
158
159  G4Vector3D Spin = aStep.GetTrack() -> GetPolarization();
160
161  G4Vector3D newSpin = SpinRotation * Spin;
162
163#ifdef G4VERBOSE
164  if (GetVerboseLevel()>2) {
165    G4double normspin = std::sqrt(Spin*Spin);
166    G4double normnewspin = std::sqrt(newSpin*newSpin);
167    //G4double cosalpha = Spin*newSpin/normspin/normnewspin;
168    //G4double alpha = std::acos(cosalpha);
169
170    G4cout << "AT REST::: PARAMETERS " << G4endl;
171    G4cout << "Initial spin  : " << Spin  << G4endl;
172    G4cout << "Delta time    : " << deltatime  << G4endl;
173    G4cout << "Rotation angle: " << rotationangle/rad  << G4endl;
174    G4cout << "New spin      : " << newSpin  << G4endl;
175    G4cout << "Checked norms : " << normspin <<" " << normnewspin << G4endl;
176  }
177#endif
178
179  return newSpin;
180
181}
Note: See TracBrowser for help on using the repository browser.